Actual source code: test15.c
slepc-3.6.1 2015-09-03
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Test DSPEP.\n\n";
24: #include <slepcds.h>
28: int main(int argc,char **argv)
29: {
31: DS ds;
32: SlepcSC sc;
33: PetscScalar *K,*C,*M,*wr,*wi,z=1.0;
34: PetscReal re,im;
35: PetscInt i,n=10,d=2,ld;
36: PetscViewer viewer;
37: PetscBool verbose;
39: SlepcInitialize(&argc,&argv,(char*)0,help);
40: PetscOptionsGetInt(NULL,"-n",&n,NULL);
41: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type PEP - n=%D.\n",n);
42: PetscOptionsHasName(NULL,"-verbose",&verbose);
44: /* Create DS object */
45: DSCreate(PETSC_COMM_WORLD,&ds);
46: DSSetType(ds,DSPEP);
47: DSSetFromOptions(ds);
48: DSPEPSetDegree(ds,d);
50: /* Set dimensions */
51: ld = n+2; /* test leading dimension larger than n */
52: DSAllocate(ds,ld);
53: DSSetDimensions(ds,n,0,0,0);
55: /* Set up viewer */
56: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
57: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
58: DSView(ds,viewer);
59: PetscViewerPopFormat(viewer);
60: if (verbose) {
61: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
62: }
64: /* Fill matrices */
65: DSGetArray(ds,DS_MAT_E0,&K);
66: for (i=0;i<n-1;i++) K[i+i*ld] = 2.0*n;
67: K[n-1+(n-1)*ld] = 1.0*n;
68: for (i=1;i<n;i++) {
69: K[i+(i-1)*ld] = -1.0*n;
70: K[(i-1)+i*ld] = -1.0*n;
71: }
72: DSRestoreArray(ds,DS_MAT_E0,&K);
73: DSGetArray(ds,DS_MAT_E1,&C);
74: C[n-1+(n-1)*ld] = 2.0*PETSC_PI/z;
75: DSRestoreArray(ds,DS_MAT_E1,&C);
76: DSGetArray(ds,DS_MAT_E2,&M);
77: for (i=0;i<n-1;i++) M[i+i*ld] = -4.0*PETSC_PI*PETSC_PI/n;
78: M[i+i*ld] = -2.0*PETSC_PI*PETSC_PI/n;
79: DSRestoreArray(ds,DS_MAT_E2,&M);
81: if (verbose) {
82: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
83: DSView(ds,viewer);
84: }
86: /* Solve */
87: PetscMalloc2(d*n,&wr,d*n,&wi);
88: DSGetSlepcSC(ds,&sc);
89: sc->comparison = SlepcCompareLargestReal;
90: sc->comparisonctx = NULL;
91: sc->map = NULL;
92: sc->mapobj = NULL;
93: DSSolve(ds,wr,wi);
94: DSSort(ds,wr,wi,NULL,NULL,NULL);
95: if (verbose) {
96: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
97: DSView(ds,viewer);
98: }
100: /* Print first eigenvalue */
101: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
102: for (i=0;i<d*n;i++) {
103: #if defined(PETSC_USE_COMPLEX)
104: re = PetscRealPart(wr[i]);
105: im = PetscImaginaryPart(wr[i]);
106: #else
107: re = wr[i];
108: im = wi[i];
109: #endif
110: if (PetscAbs(im)<1e-10) {
111: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
112: } else {
113: PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im);
114: }
115: }
117: PetscFree2(wr,wi);
118: DSDestroy(&ds);
119: SlepcFinalize();
120: return 0;
121: }