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);
83: EPSSetProblemType(eps,EPS_HEP);
85: /*
86: Set solver parameters at runtime
87: */
88: EPSSetFromOptions(eps);
90: /*
91: Initialize shell spectral transformation if selected by user
92: */
93: EPSGetST(eps,&st);
94: PetscTypeCompare((PetscObject)st,STSHELL,&isShell);
95: if (isShell) {
96: /* (Optional) Create a context for the user-defined spectral tranform;
97: this context can be defined to contain any application-specific data. */
98: SampleShellSTCreate(&shell);
100: /* (Required) Set the user-defined routine for applying the operator */
101: STShellSetApply(st,SampleShellSTApply,(void*)shell);
103: /* (Optional) Set the user-defined routine for back-transformation */
104: STShellSetBackTransform(st,SampleShellSTBackTransform);
106: /* (Optional) Set a name for the transformation, used for STView() */
107: STShellSetName(st,"MyTransformation");
109: /* (Optional) Do any setup required for the new transformation */
110: SampleShellSTSetUp(shell,st);
111: }
113: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114: Solve the eigensystem
115: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117: EPSSolve(eps);
118: EPSGetIterationNumber(eps, &its);
119: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
120:
122: /*
123: Optional: Get some information from the solver and display it
124: */
125: EPSGetType(eps,&type);
126: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
127: EPSGetDimensions(eps,&nev,PETSC_NULL);
128: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
129:
130: EPSGetTolerances(eps,&tol,&maxit);
131: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
132:
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Display solution and clean up
136: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138: /*
139: Get number of converged approximate eigenpairs
140: */
141: EPSGetConverged(eps,&nconv);
142: PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %d\n\n",nconv);
143:
145: if (nconv>0) {
146: /*
147: Display eigenvalues and relative errors
148: */
149: PetscPrintf(PETSC_COMM_WORLD,
150: " k ||Ax-kx||/||kx||\n"
151: " ----------------- ------------------\n" );
153: for( i=0; i<nconv; i++ ) {
154: /*
155: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
156: ki (imaginary part)
157: */
158: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
159: /*
160: Compute the relative error associated to each eigenpair
161: */
162: EPSComputeRelativeError(eps,i,&error);
164: #ifdef PETSC_USE_COMPLEX
165: re = PetscRealPart(kr);
166: im = PetscImaginaryPart(kr);
167: #else
168: re = kr;
169: im = ki;
170: #endif
171: if (im!=0.0) {
172: PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12f\n",re,im,error);
173: } else {
174: PetscPrintf(PETSC_COMM_WORLD," %12f %12f\n",re,error);
175: }
176: }
177: PetscPrintf(PETSC_COMM_WORLD,"\n" );
178: }
179:
180: /*
181: Free work space
182: */
183: if (isShell) {
184: SampleShellSTDestroy(shell);
185: }
186: EPSDestroy(eps);
187: MatDestroy(A);
188: SlepcFinalize();
189: return 0;
190: }
192: /***********************************************************************/
193: /* Routines for a user-defined shell transformation */
194: /***********************************************************************/
198: /*
199: SampleShellSTCreate - This routine creates a user-defined
200: spectral transformation context.
202: Output Parameter:
203: . shell - user-defined spectral transformation context
204: */
205: int SampleShellSTCreate(SampleShellST **shell)
206: {
207: SampleShellST *newctx;
208: int ierr;
210: PetscNew(SampleShellST,&newctx);
211: KSPCreate(PETSC_COMM_WORLD,&newctx->ksp);
212: KSPAppendOptionsPrefix(newctx->ksp,"st_");
213: *shell = newctx;
214: return 0;
215: }
216: /* ------------------------------------------------------------------- */
219: /*
220: SampleShellSTSetUp - This routine sets up a user-defined
221: spectral transformation context.
223: Input Parameters:
224: . shell - user-defined spectral transformation context
225: . st - spectral transformation context containing the operator matrices
227: Output Parameter:
228: . shell - fully set up user-defined transformation context
230: Notes:
231: In this example, the user-defined transformation is simply OP=A^-1.
232: Therefore, the eigenpairs converge in reversed order. The KSP object
233: used for the solution of linear systems with A is handled via the
234: user-defined context SampleShellST.
235: */
236: int SampleShellSTSetUp(SampleShellST *shell,ST st)
237: {
238: Mat A,B;
239: int ierr;
241: STGetOperators( st, &A, &B );
242: if (B) { SETERRQ(0,"Warning: This transformation is not intended for generalized problems"); }
243: KSPSetOperators(shell->ksp,A,A,DIFFERENT_NONZERO_PATTERN);
244: KSPSetFromOptions(shell->ksp);
246: return 0;
247: }
248: /* ------------------------------------------------------------------- */
251: /*
252: SampleShellSTApply - This routine demonstrates the use of a
253: user-provided spectral transformation.
255: Input Parameters:
256: . ctx - optional user-defined context, as set by STShellSetApply()
257: . x - input vector
259: Output Parameter:
260: . y - output vector
262: Notes:
263: The transformation implemented in this code is just OP=A^-1 and
264: therefore it is of little use, merely as an example of working with
265: a STSHELL.
266: */
267: int SampleShellSTApply(void *ctx,Vec x,Vec y)
268: {
269: SampleShellST *shell = (SampleShellST*)ctx;
270: int ierr;
272: KSPSolve(shell->ksp,x,y);
274: return 0;
275: }
276: /* ------------------------------------------------------------------- */
279: /*
280: SampleShellSTBackTransform - This routine demonstrates the use of a
281: user-provided spectral transformation.
283: Input Parameters:
284: . ctx - optional user-defined context, as set by STShellSetApply()
285: . eigr - pointer to real part of eigenvalues
286: . eigi - pointer to imaginary part of eigenvalues
288: Output Parameters:
289: . eigr - modified real part of eigenvalues
290: . eigi - modified imaginary part of eigenvalues
292: Notes:
293: This code implements the back transformation of eigenvalues in
294: order to retrieve the eigenvalues of the original problem. In this
295: example, simply set k_i = 1/k_i.
296: */
297: int SampleShellSTBackTransform(void *ctx,PetscScalar *eigr,PetscScalar *eigi)
298: {
299: *eigr = 1.0 / *eigr;
301: return 0;
302: }
303: /* ------------------------------------------------------------------- */
306: /*
307: SampleShellSTDestroy - This routine destroys a user-defined
308: spectral transformation context.
310: Input Parameter:
311: . shell - user-defined spectral transformation context
312: */
313: int SampleShellSTDestroy(SampleShellST *shell)
314: {
317: KSPDestroy(shell->ksp);
318: PetscFree(shell);
320: return 0;
321: }