Actual source code: ex1.c
2: static char help[] = "Standard symmetric eigenproblem corresponding to the Laplacian operator in 1 dimension.\n\n"
3: "The command line options are:\n"
4: " -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n\n";
6: #include slepceps.h
10: int main( int argc, char **argv )
11: {
12: Mat A; /* operator matrix */
13: EPS eps; /* eigenproblem solver context */
14: EPSType type;
15: PetscReal error, tol,re, im;
16: PetscScalar kr, ki;
17: Vec xr, xi;
19: PetscInt n=30, i, Istart, Iend, col[3];
20: int nev, maxit,its, nconv;
21: PetscTruth FirstBlock=PETSC_FALSE, LastBlock=PETSC_FALSE;
22: PetscScalar value[3];
24: SlepcInitialize(&argc,&argv,(char*)0,help);
26: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
27: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%d\n\n",n);
29: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30: Compute the operator matrix that defines the eigensystem, Ax=kx
31: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
33: MatCreate(PETSC_COMM_WORLD,&A);
34: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
35: MatSetFromOptions(A);
36:
37: MatGetOwnershipRange(A,&Istart,&Iend);
38: if (Istart==0) FirstBlock=PETSC_TRUE;
39: if (Iend==n) LastBlock=PETSC_TRUE;
40: value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
41: for( i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++ ) {
42: col[0]=i-1; col[1]=i; col[2]=i+1;
43: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
44: }
45: if (LastBlock) {
46: i=n-1; col[0]=n-2; col[1]=n-1;
47: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
48: }
49: if (FirstBlock) {
50: i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
51: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
52: }
54: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
57: MatGetVecs(A,PETSC_NULL,&xr);
58: MatGetVecs(A,PETSC_NULL,&xi);
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Create the eigensolver and set various options
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: /*
64: Create eigensolver context
65: */
66: EPSCreate(PETSC_COMM_WORLD,&eps);
68: /*
69: Set operators. In this case, it is a standard eigenvalue problem
70: */
71: EPSSetOperators(eps,A,PETSC_NULL);
72: EPSSetProblemType(eps,EPS_HEP);
74: /*
75: Set solver parameters at runtime
76: */
77: EPSSetFromOptions(eps);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Solve the eigensystem
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: EPSSolve(eps);
84: EPSGetIterationNumber(eps, &its);
85: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
86: /*
87: Optional: Get some information from the solver and display it
88: */
89: EPSGetType(eps,&type);
90: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
91: EPSGetDimensions(eps,&nev,PETSC_NULL);
92: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
93: EPSGetTolerances(eps,&tol,&maxit);
94: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Display solution and clean up
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: /*
100: Get number of converged approximate eigenpairs
101: */
102: EPSGetConverged(eps,&nconv);
103: PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %d\n\n",nconv);
105: if (nconv>0) {
106: /*
107: Display eigenvalues and relative errors
108: */
109: PetscPrintf(PETSC_COMM_WORLD,
110: " k ||Ax-kx||/||kx||\n"
111: " ----------------- ------------------\n" );
113: for( i=0; i<nconv; i++ ) {
114: /*
115: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
116: ki (imaginary part)
117: */
118: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
119: /*
120: Compute the relative error associated to each eigenpair
121: */
122: EPSComputeRelativeError(eps,i,&error);
124: #ifdef PETSC_USE_COMPLEX
125: re = PetscRealPart(kr);
126: im = PetscImaginaryPart(kr);
127: #else
128: re = kr;
129: im = ki;
130: #endif
131: if (im!=0.0) {
132: PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",re,im,error);
133: } else {
134: PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",re,error);
135: }
136: }
137: PetscPrintf(PETSC_COMM_WORLD,"\n" );
138: }
139:
140: /*
141: Free work space
142: */
143: EPSDestroy(eps);
144: MatDestroy(A);
145: VecDestroy(xr);
146: VecDestroy(xi);
147: SlepcFinalize();
148: return 0;
149: }