Actual source code: ex10.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
6: This file is part of SLEPc. See the README file for conditions of use
7: and additional information.
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Illustrates the use of shell spectral transformations. "
12: "The problem to be solved is the same as ex1.c and"
13: "corresponds to the Laplacian operator in 1 dimension.\n\n"
14: "The command line options are:\n"
15: " -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n\n";
17: #include slepceps.h
19: /* Define context for user-provided spectral transformation */
20: typedef struct {
21: KSP ksp;
22: } SampleShellST;
24: /* Declare routines for user-provided spectral transformation */
25: PetscErrorCode SampleShellSTCreate(SampleShellST**);
26: PetscErrorCode SampleShellSTSetUp(SampleShellST*,ST);
27: PetscErrorCode SampleShellSTApply(void*,Vec,Vec);
28: PetscErrorCode SampleShellSTBackTransform(void*,PetscScalar*,PetscScalar*);
29: PetscErrorCode SampleShellSTDestroy(SampleShellST*);
33: int main( int argc, char **argv )
34: {
35: Mat A; /* operator matrix */
36: EPS eps; /* eigenproblem solver context */
37: ST st; /* spectral transformation context */
38: SampleShellST *shell; /* user-defined spectral transform context */
39: EPSType type;
40: PetscReal error, tol, re, im;
41: PetscScalar kr, ki;
43: int nev, maxit, its, nconv;
44: PetscInt n=30, i, col[3], Istart, Iend, FirstBlock=0, LastBlock=0;
45: PetscScalar value[3];
46: PetscTruth isShell;
48: SlepcInitialize(&argc,&argv,(char*)0,help);
50: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
51: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem (shell-enabled), n=%d\n\n",n);
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Compute the operator matrix that defines the eigensystem, Ax=kx
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: MatCreate(PETSC_COMM_WORLD,&A);
58: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
59: MatSetFromOptions(A);
60:
61: MatGetOwnershipRange(A,&Istart,&Iend);
62: if (Istart==0) FirstBlock=PETSC_TRUE;
63: if (Iend==n) LastBlock=PETSC_TRUE;
64: value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
65: for( i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++ ) {
66: col[0]=i-1; col[1]=i; col[2]=i+1;
67: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
68: }
69: if (LastBlock) {
70: i=n-1; col[0]=n-2; col[1]=n-1;
71: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
72: }
73: if (FirstBlock) {
74: i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
75: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
76: }
78: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
79: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
81: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: Create the eigensolver and set various options
83: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: /*
86: Create eigensolver context
87: */
88: EPSCreate(PETSC_COMM_WORLD,&eps);
90: /*
91: Set operators. In this case, it is a standard eigenvalue problem
92: */
93: EPSSetOperators(eps,A,PETSC_NULL);
94: EPSSetProblemType(eps,EPS_HEP);
96: /*
97: Set solver parameters at runtime
98: */
99: EPSSetFromOptions(eps);
101: /*
102: Initialize shell spectral transformation if selected by user
103: */
104: EPSGetST(eps,&st);
105: PetscTypeCompare((PetscObject)st,STSHELL,&isShell);
106: if (isShell) {
107: /* (Optional) Create a context for the user-defined spectral tranform;
108: this context can be defined to contain any application-specific data. */
109: SampleShellSTCreate(&shell);
111: /* (Required) Set the user-defined routine for applying the operator */
112: STShellSetApply(st,SampleShellSTApply);
113: STShellSetContext(st,shell);
115: /* (Optional) Set the user-defined routine for back-transformation */
116: STShellSetBackTransform(st,SampleShellSTBackTransform);
118: /* (Optional) Set a name for the transformation, used for STView() */
119: STShellSetName(st,"MyTransformation");
121: /* (Optional) Do any setup required for the new transformation */
122: SampleShellSTSetUp(shell,st);
123: }
125: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126: Solve the eigensystem
127: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129: EPSSolve(eps);
130: EPSGetIterationNumber(eps, &its);
131: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
133: /*
134: Optional: Get some information from the solver and display it
135: */
136: EPSGetType(eps,&type);
137: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
138: EPSGetDimensions(eps,&nev,PETSC_NULL);
139: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
140: EPSGetTolerances(eps,&tol,&maxit);
141: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
143: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144: Display solution and clean up
145: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147: /*
148: Get number of converged approximate eigenpairs
149: */
150: EPSGetConverged(eps,&nconv);
151: PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %d\n\n",nconv);
153: if (nconv>0) {
154: /*
155: Display eigenvalues and relative errors
156: */
157: PetscPrintf(PETSC_COMM_WORLD,
158: " k ||Ax-kx||/||kx||\n"
159: " ----------------- ------------------\n" );
161: for( i=0; i<nconv; i++ ) {
162: /*
163: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
164: ki (imaginary part)
165: */
166: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
167: /*
168: Compute the relative error associated to each eigenpair
169: */
170: EPSComputeRelativeError(eps,i,&error);
172: #ifdef PETSC_USE_COMPLEX
173: re = PetscRealPart(kr);
174: im = PetscImaginaryPart(kr);
175: #else
176: re = kr;
177: im = ki;
178: #endif
179: if (im!=0.0) {
180: PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",re,im,error);
181: } else {
182: PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",re,error);
183: }
184: }
185: PetscPrintf(PETSC_COMM_WORLD,"\n" );
186: }
187:
188: /*
189: Free work space
190: */
191: if (isShell) {
192: SampleShellSTDestroy(shell);
193: }
194: EPSDestroy(eps);
195: MatDestroy(A);
196: SlepcFinalize();
197: return 0;
198: }
200: /***********************************************************************/
201: /* Routines for a user-defined shell spectral transformation */
202: /***********************************************************************/
206: /*
207: SampleShellSTCreate - This routine creates a user-defined
208: spectral transformation context.
210: Output Parameter:
211: . shell - user-defined spectral transformation context
212: */
213: PetscErrorCode SampleShellSTCreate(SampleShellST **shell)
214: {
215: SampleShellST *newctx;
219: PetscNew(SampleShellST,&newctx);
220: KSPCreate(PETSC_COMM_WORLD,&newctx->ksp);
221: KSPAppendOptionsPrefix(newctx->ksp,"st_");
222: *shell = newctx;
223: return(0);
224: }
225: /* ------------------------------------------------------------------- */
228: /*
229: SampleShellSTSetUp - This routine sets up a user-defined
230: spectral transformation context.
232: Input Parameters:
233: . shell - user-defined spectral transformation context
234: . st - spectral transformation context containing the operator matrices
236: Output Parameter:
237: . shell - fully set up user-defined transformation context
239: Notes:
240: In this example, the user-defined transformation is simply OP=A^-1.
241: Therefore, the eigenpairs converge in reversed order. The KSP object
242: used for the solution of linear systems with A is handled via the
243: user-defined context SampleShellST.
244: */
245: PetscErrorCode SampleShellSTSetUp(SampleShellST *shell,ST st)
246: {
247: Mat A,B;
251: STGetOperators(st,&A,&B);
252: if (B) { SETERRQ(0,"Warning: This transformation is not intended for generalized problems"); }
253: KSPSetOperators(shell->ksp,A,A,DIFFERENT_NONZERO_PATTERN);
254: KSPSetFromOptions(shell->ksp);
255: return(0);
256: }
257: /* ------------------------------------------------------------------- */
260: /*
261: SampleShellSTApply - This routine demonstrates the use of a
262: user-provided spectral transformation.
264: Input Parameters:
265: . ctx - optional user-defined context, as set by STShellSetContext()
266: . x - input vector
268: Output Parameter:
269: . y - output vector
271: Notes:
272: The transformation implemented in this code is just OP=A^-1 and
273: therefore it is of little use, merely as an example of working with
274: a STSHELL.
275: */
276: PetscErrorCode SampleShellSTApply(void *ctx,Vec x,Vec y)
277: {
278: SampleShellST *shell = (SampleShellST*)ctx;
282: KSPSolve(shell->ksp,x,y);
283: return(0);
284: }
285: /* ------------------------------------------------------------------- */
288: /*
289: SampleShellSTBackTransform - This routine demonstrates the use of a
290: user-provided spectral transformation.
292: Input Parameters:
293: . ctx - optional user-defined context, as set by STShellSetContext()
294: . eigr - pointer to real part of eigenvalues
295: . eigi - pointer to imaginary part of eigenvalues
297: Output Parameters:
298: . eigr - modified real part of eigenvalues
299: . eigi - modified imaginary part of eigenvalues
301: Notes:
302: This code implements the back transformation of eigenvalues in
303: order to retrieve the eigenvalues of the original problem. In this
304: example, simply set k_i = 1/k_i.
305: */
306: PetscErrorCode SampleShellSTBackTransform(void *ctx,PetscScalar *eigr,PetscScalar *eigi)
307: {
309: *eigr = 1.0 / *eigr;
310: return(0);
311: }
312: /* ------------------------------------------------------------------- */
315: /*
316: SampleShellSTDestroy - This routine destroys a user-defined
317: spectral transformation context.
319: Input Parameter:
320: . shell - user-defined spectral transformation context
321: */
322: PetscErrorCode SampleShellSTDestroy(SampleShellST *shell)
323: {
327: KSPDestroy(shell->ksp);
328: PetscFree(shell);
329: return(0);
330: }