Actual source code: ex10.c
slepc-3.5.2 2014-10-10
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2014, 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[] = "Illustrates the use of shell spectral transformations. "
23: "The problem to be solved is the same as ex1.c and"
24: "corresponds to the Laplacian operator in 1 dimension.\n\n"
25: "The command line options are:\n"
26: " -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n\n";
28: #include <slepceps.h>
30: /* Define context for user-provided spectral transformation */
31: typedef struct {
32: KSP ksp;
33: } SampleShellST;
35: /* Declare routines for user-provided spectral transformation */
36: PetscErrorCode STCreate_User(SampleShellST**);
37: PetscErrorCode STSetUp_User(SampleShellST*,ST);
38: PetscErrorCode STApply_User(ST,Vec,Vec);
39: PetscErrorCode STBackTransform_User(ST,PetscInt,PetscScalar*,PetscScalar*);
40: PetscErrorCode STDestroy_User(SampleShellST*);
44: int main (int argc,char **argv)
45: {
46: Mat A; /* operator matrix */
47: EPS eps; /* eigenproblem solver context */
48: ST st; /* spectral transformation context */
49: SampleShellST *shell; /* user-defined spectral transform context */
50: EPSType type;
51: PetscScalar value[3];
52: PetscInt n=30,i,col[3],Istart,Iend,FirstBlock=0,LastBlock=0,nev;
53: PetscBool isShell;
56: SlepcInitialize(&argc,&argv,(char*)0,help);
58: PetscOptionsGetInt(NULL,"-n",&n,NULL);
59: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem (shell-enabled), n=%D\n\n",n);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Compute the operator matrix that defines the eigensystem, Ax=kx
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: MatCreate(PETSC_COMM_WORLD,&A);
66: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
67: MatSetFromOptions(A);
68: MatSetUp(A);
70: MatGetOwnershipRange(A,&Istart,&Iend);
71: if (Istart==0) FirstBlock=PETSC_TRUE;
72: if (Iend==n) LastBlock=PETSC_TRUE;
73: value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
74: for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
75: col[0]=i-1; col[1]=i; col[2]=i+1;
76: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
77: }
78: if (LastBlock) {
79: i=n-1; col[0]=n-2; col[1]=n-1;
80: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
81: }
82: if (FirstBlock) {
83: i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
84: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
85: }
87: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
88: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Create the eigensolver and set various options
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: /*
95: Create eigensolver context
96: */
97: EPSCreate(PETSC_COMM_WORLD,&eps);
99: /*
100: Set operators. In this case, it is a standard eigenvalue problem
101: */
102: EPSSetOperators(eps,A,NULL);
103: EPSSetProblemType(eps,EPS_HEP);
105: /*
106: Set solver parameters at runtime
107: */
108: EPSSetFromOptions(eps);
110: /*
111: Initialize shell spectral transformation if selected by user
112: */
113: EPSGetST(eps,&st);
114: PetscObjectTypeCompare((PetscObject)st,STSHELL,&isShell);
115: if (isShell) {
116: /* Change sorting criterion since this ST example computes values
117: closest to 0 */
118: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
120: /* (Optional) Create a context for the user-defined spectral tranform;
121: this context can be defined to contain any application-specific data. */
122: STCreate_User(&shell);
124: /* (Required) Set the user-defined routine for applying the operator */
125: STShellSetApply(st,STApply_User);
126: STShellSetContext(st,shell);
128: /* (Optional) Set the user-defined routine for back-transformation */
129: STShellSetBackTransform(st,STBackTransform_User);
131: /* (Optional) Set a name for the transformation, used for STView() */
132: PetscObjectSetName((PetscObject)st,"MyTransformation");
134: /* (Optional) Do any setup required for the new transformation */
135: STSetUp_User(shell,st);
136: }
138: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139: Solve the eigensystem
140: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142: EPSSolve(eps);
144: /*
145: Optional: Get some information from the solver and display it
146: */
147: EPSGetType(eps,&type);
148: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
149: EPSGetDimensions(eps,&nev,NULL,NULL);
150: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
152: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153: Display solution and clean up
154: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156: EPSPrintSolution(eps,NULL);
157: if (isShell) {
158: STDestroy_User(shell);
159: }
160: EPSDestroy(&eps);
161: MatDestroy(&A);
162: SlepcFinalize();
163: return 0;
164: }
166: /***********************************************************************/
167: /* Routines for a user-defined shell spectral transformation */
168: /***********************************************************************/
172: /*
173: STCreate_User - This routine creates a user-defined
174: spectral transformation context.
176: Output Parameter:
177: . shell - user-defined spectral transformation context
178: */
179: PetscErrorCode STCreate_User(SampleShellST **shell)
180: {
181: SampleShellST *newctx;
185: PetscNew(&newctx);
186: KSPCreate(PETSC_COMM_WORLD,&newctx->ksp);
187: KSPAppendOptionsPrefix(newctx->ksp,"st_");
188: *shell = newctx;
189: return(0);
190: }
191: /* ------------------------------------------------------------------- */
194: /*
195: STSetUp_User - This routine sets up a user-defined
196: spectral transformation context.
198: Input Parameters:
199: . shell - user-defined spectral transformation context
200: . st - spectral transformation context containing the operator matrices
202: Output Parameter:
203: . shell - fully set up user-defined transformation context
205: Notes:
206: In this example, the user-defined transformation is simply OP=A^-1.
207: Therefore, the eigenpairs converge in reversed order. The KSP object
208: used for the solution of linear systems with A is handled via the
209: user-defined context SampleShellST.
210: */
211: PetscErrorCode STSetUp_User(SampleShellST *shell,ST st)
212: {
213: Mat A;
217: STGetOperators(st,0,&A);
218: KSPSetOperators(shell->ksp,A,A);
219: KSPSetFromOptions(shell->ksp);
220: return(0);
221: }
222: /* ------------------------------------------------------------------- */
225: /*
226: STApply_User - This routine demonstrates the use of a
227: user-provided spectral transformation.
229: Input Parameters:
230: . ctx - optional user-defined context, as set by STShellSetContext()
231: . x - input vector
233: Output Parameter:
234: . y - output vector
236: Notes:
237: The transformation implemented in this code is just OP=A^-1 and
238: therefore it is of little use, merely as an example of working with
239: a STSHELL.
240: */
241: PetscErrorCode STApply_User(ST st,Vec x,Vec y)
242: {
243: SampleShellST *shell;
247: STShellGetContext(st,(void**)&shell);
248: KSPSolve(shell->ksp,x,y);
249: return(0);
250: }
251: /* ------------------------------------------------------------------- */
254: /*
255: STBackTransform_User - This routine demonstrates the use of a
256: user-provided spectral transformation.
258: Input Parameters:
259: + ctx - optional user-defined context, as set by STShellSetContext()
260: . eigr - pointer to real part of eigenvalues
261: - eigi - pointer to imaginary part of eigenvalues
263: Output Parameters:
264: + eigr - modified real part of eigenvalues
265: - eigi - modified imaginary part of eigenvalues
267: Notes:
268: This code implements the back transformation of eigenvalues in
269: order to retrieve the eigenvalues of the original problem. In this
270: example, simply set k_i = 1/k_i.
271: */
272: PetscErrorCode STBackTransform_User(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
273: {
274: PetscInt j;
277: for (j=0;j<n;j++) {
278: eigr[j] = 1.0 / eigr[j];
279: }
280: return(0);
281: }
282: /* ------------------------------------------------------------------- */
285: /*
286: STDestroy_User - This routine destroys a user-defined
287: spectral transformation context.
289: Input Parameter:
290: . shell - user-defined spectral transformation context
291: */
292: PetscErrorCode STDestroy_User(SampleShellST *shell)
293: {
297: KSPDestroy(&shell->ksp);
298: PetscFree(shell);
299: return(0);
300: }