Actual source code: ex10.c
1: static char help[] = "Illustrates the use of shell spectral transformations. "
2: "The problem to be solved is the same as ex1.c and"
3: "corresponds to the Laplacian operator in 1 dimension.\n\n"
4: "The command line options are:\n\n"
5: " -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n\n";
7: #include slepceps.h
9: /* Define context for user-provided spectral transformation */
10: typedef struct {
11: KSP ksp;
12: } SampleShellST;
14: /* Declare routines for user-provided spectral transformation */
15: extern int SampleShellSTCreate(SampleShellST**);
16: extern int SampleShellSTSetUp(SampleShellST*,ST);
17: extern int SampleShellSTApply(void*,Vec,Vec);
18: extern int SampleShellSTBackTransform(void*,PetscScalar*,PetscScalar*);
19: extern int SampleShellSTDestroy(SampleShellST*);
23: int main( int argc, char **argv )
24: {
25: Mat A; /* operator matrix */
26: EPS eps; /* eigenproblem solver context */
27: ST st; /* spectral transformation context */
28: SampleShellST *shell; /* user-defined spectral transform context */
29: EPSType type;
30: PetscReal error, tol, re, im;
31: PetscScalar kr, ki;
32: int n=30, nev, ierr, maxit, i, its, nconv,
33: col[3], Istart, Iend, FirstBlock=0, LastBlock=0;
34: PetscScalar value[3];
35: PetscTruth isShell;
37: SlepcInitialize(&argc,&argv,(char*)0,help);
39: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
40: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem (shell-enabled), n=%d\n\n",n);
41:
43: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44: Compute the operator matrix that defines the eigensystem, Ax=kx
45: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&A);
48: MatSetFromOptions(A);
49:
50: MatGetOwnershipRange(A,&Istart,&Iend);
51: if (Istart==0) FirstBlock=PETSC_TRUE;
52: if (Iend==n) LastBlock=PETSC_TRUE;
53: value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
54: for( i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++ ) {
55: col[0]=i-1; col[1]=i; col[2]=i+1;
56: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
57: }
58: if (LastBlock) {
59: i=n-1; col[0]=n-2; col[1]=n-1;
60: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
61: }
62: if (FirstBlock) {
63: i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
64: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
65: }
67: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
68: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Create the eigensolver and set various options
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: /*
75: Create eigensolver context
76: */
77: EPSCreate(PETSC_COMM_WORLD,&eps);
79: /*
80: Set operators. In this case, it is a standard eigenvalue problem
81: */
82: EPSSetOperators(eps,A,PETSC_NULL);
84: /*
85: Set solver parameters at runtime
86: */
87: EPSSetFromOptions(eps);
89: /*
90: Initialize shell spectral transformation if selected by user
91: */
92: EPSGetST(eps,&st);
93: PetscTypeCompare((PetscObject)st,STSHELL,&isShell);
94: if (isShell) {
95: /* (Optional) Create a context for the user-defined spectral tranform;
96: this context can be defined to contain any application-specific data. */
97: SampleShellSTCreate(&shell);
99: /* (Required) Set the user-defined routine for applying the operator */
100: STShellSetApply(st,SampleShellSTApply,(void*)shell);
102: /* (Optional) Set the user-defined routine for back-transformation */
103: STShellSetBackTransform(st,SampleShellSTBackTransform);
105: /* (Optional) Set a name for the transformation, used for STView() */
106: STShellSetName(st,"MyTransformation");
108: /* (Optional) Do any setup required for the new transformation */
109: SampleShellSTSetUp(shell,st);
110: }
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Solve the eigensystem
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: EPSSolve(eps);
117: EPSGetIterationNumber(eps, &its);
118: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
119:
121: /*
122: Optional: Get some information from the solver and display it
123: */
124: EPSGetType(eps,&type);
125: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
126: EPSGetDimensions(eps,&nev,PETSC_NULL);
127: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
128:
129: EPSGetTolerances(eps,&tol,&maxit);
130: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
131:
133: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134: Display solution and clean up
135: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137: /*
138: Get number of converged approximate eigenpairs
139: */
140: EPSGetConverged(eps,&nconv);
141: PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %d\n\n",nconv);
142:
144: if (nconv>0) {
145: /*
146: Display eigenvalues and relative errors
147: */
148: PetscPrintf(PETSC_COMM_WORLD,
149: " k ||Ax-kx||/||kx||\n"
150: " ----------------- ------------------\n" );
152: for( i=0; i<nconv; i++ ) {
153: /*
154: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
155: ki (imaginary part)
156: */
157: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
158: /*
159: Compute the relative error associated to each eigenpair
160: */
161: EPSComputeRelativeError(eps,i,&error);
163: #ifdef PETSC_USE_COMPLEX
164: re = PetscRealPart(kr);
165: im = PetscImaginaryPart(kr);
166: #else
167: re = kr;
168: im = ki;
169: #endif
170: if (im!=0.0) {
171: PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12f\n",re,im,error);
172: } else {
173: PetscPrintf(PETSC_COMM_WORLD," %12f %12f\n",re,error);
174: }
175: }
176: PetscPrintf(PETSC_COMM_WORLD,"\n" );
177: }
178:
179: /*
180: Free work space
181: */
182: if (isShell) {
183: SampleShellSTDestroy(shell);
184: }
185: EPSDestroy(eps);
186: MatDestroy(A);
187: SlepcFinalize();
188: return 0;
189: }
191: /***********************************************************************/
192: /* Routines for a user-defined shell transformation */
193: /***********************************************************************/
197: /*
198: SampleShellSTCreate - This routine creates a user-defined
199: spectral transformation context.
201: Output Parameter:
202: . shell - user-defined spectral transformation context
203: */
204: int SampleShellSTCreate(SampleShellST **shell)
205: {
206: SampleShellST *newctx;
207: int ierr;
209: PetscNew(SampleShellST,&newctx);
210: KSPCreate(PETSC_COMM_WORLD,&newctx->ksp);
211: KSPAppendOptionsPrefix(newctx->ksp,"st_");
212: *shell = newctx;
213: return 0;
214: }
215: /* ------------------------------------------------------------------- */
218: /*
219: SampleShellSTSetUp - This routine sets up a user-defined
220: spectral transformation context.
222: Input Parameters:
223: . shell - user-defined spectral transformation context
224: . st - spectral transformation context containing the operator matrices
226: Output Parameter:
227: . shell - fully set up user-defined transformation context
229: Notes:
230: In this example, the user-defined transformation is simply OP=A^-1.
231: Therefore, the eigenpairs converge in reversed order. The KSP object
232: used for the solution of linear systems with A is handled via the
233: user-defined context SampleShellST.
234: */
235: int SampleShellSTSetUp(SampleShellST *shell,ST st)
236: {
237: Mat A,B;
238: int ierr;
240: STGetOperators( st, &A, &B );
241: if (B) { SETERRQ(0,"Warning: This transformation is not intended for generalized problems"); }
242: KSPSetOperators(shell->ksp,A,A,DIFFERENT_NONZERO_PATTERN);
243: KSPSetFromOptions(shell->ksp);
245: return 0;
246: }
247: /* ------------------------------------------------------------------- */
250: /*
251: SampleShellSTApply - This routine demonstrates the use of a
252: user-provided spectral transformation.
254: Input Parameters:
255: . ctx - optional user-defined context, as set by STShellSetApply()
256: . x - input vector
258: Output Parameter:
259: . y - output vector
261: Notes:
262: The transformation implemented in this code is just OP=A^-1 and
263: therefore it is of little use, merely as an example of working with
264: a STSHELL.
265: */
266: int SampleShellSTApply(void *ctx,Vec x,Vec y)
267: {
268: SampleShellST *shell = (SampleShellST*)ctx;
269: int ierr;
271: KSPSetRhs(shell->ksp,x);
272: KSPSetSolution(shell->ksp,y);
273: KSPSolve(shell->ksp);
275: return 0;
276: }
277: /* ------------------------------------------------------------------- */
280: /*
281: SampleShellSTBackTransform - This routine demonstrates the use of a
282: user-provided spectral transformation.
284: Input Parameters:
285: . ctx - optional user-defined context, as set by STShellSetApply()
286: . eigr - pointer to real part of eigenvalues
287: . eigi - pointer to imaginary part of eigenvalues
289: Output Parameters:
290: . eigr - modified real part of eigenvalues
291: . eigi - modified imaginary part of eigenvalues
293: Notes:
294: This code implements the back transformation of eigenvalues in
295: order to retrieve the eigenvalues of the original problem. In this
296: example, simply set k_i = 1/k_i.
297: */
298: int SampleShellSTBackTransform(void *ctx,PetscScalar *eigr,PetscScalar *eigi)
299: {
300: *eigr = 1.0 / *eigr;
302: return 0;
303: }
304: /* ------------------------------------------------------------------- */
307: /*
308: SampleShellSTDestroy - This routine destroys a user-defined
309: spectral transformation context.
311: Input Parameter:
312: . shell - user-defined spectral transformation context
313: */
314: int SampleShellSTDestroy(SampleShellST *shell)
315: {
318: KSPDestroy(shell->ksp);
319: PetscFree(shell);
321: return 0;
322: }