Actual source code: ex35.c
slepc-3.20.0 2023-09-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Shell spectral transformations with a non-injective mapping. "
12: "Implements spectrum folding for the 2-D Laplacian, as in ex24.c.\n\n"
13: "The command line options are:\n"
14: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
15: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n";
17: #include <slepceps.h>
19: /* Context for spectrum folding spectral transformation */
20: typedef struct {
21: Mat A;
22: Vec w;
23: PetscScalar target;
24: } FoldShellST;
26: /* Routines for shell spectral transformation */
27: PetscErrorCode STCreate_Fold(Mat,PetscScalar,FoldShellST**);
28: PetscErrorCode STApply_Fold(ST,Vec,Vec);
29: PetscErrorCode STDestroy_Fold(FoldShellST*);
31: int main (int argc,char **argv)
32: {
33: Mat A; /* operator matrix */
34: EPS eps; /* eigenproblem solver context */
35: ST st; /* spectral transformation context */
36: FoldShellST *fold; /* user-defined spectral transform context */
37: EPSType type;
38: PetscInt N,n=10,m,i,j,II,Istart,Iend,nev;
39: PetscBool isShell,terse,flag;
40: PetscScalar target=1.1;
42: PetscFunctionBeginUser;
43: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
45: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
46: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
47: if (!flag) m = n;
48: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-target",&target,NULL));
49: N = n*m;
50: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum Folding via shell ST, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid) target=%3.2f\n\n",N,n,m,(double)PetscRealPart(target)));
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Compute the 5-point stencil Laplacian
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
57: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
58: PetscCall(MatSetFromOptions(A));
59: PetscCall(MatSetUp(A));
61: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
62: for (II=Istart;II<Iend;II++) {
63: i = II/n; j = II-i*n;
64: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
65: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
66: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
67: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
68: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
69: }
70: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
71: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Create the eigensolver and set various options
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
78: PetscCall(EPSSetOperators(eps,A,NULL));
79: PetscCall(EPSSetProblemType(eps,EPS_HEP));
80: PetscCall(EPSSetTarget(eps,target));
81: PetscCall(EPSGetST(eps,&st));
82: PetscCall(STSetType(st,STSHELL));
83: PetscCall(EPSSetFromOptions(eps));
85: /*
86: Initialize shell spectral transformation
87: */
88: PetscCall(PetscObjectTypeCompare((PetscObject)st,STSHELL,&isShell));
89: if (isShell) {
90: /* Change sorting criterion since this shell ST computes eigenvalues
91: of the transformed operator closest to 0 */
92: PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
94: /* Create the context for the user-defined spectral transform */
95: PetscCall(STCreate_Fold(A,target,&fold));
96: PetscCall(STShellSetContext(st,fold));
98: /* Set callback function for applying the operator (in this case we do not
99: provide a back-transformation callback since the mapping is not one-to-one) */
100: PetscCall(STShellSetApply(st,STApply_Fold));
101: PetscCall(PetscObjectSetName((PetscObject)st,"STFOLD"));
102: }
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Solve the eigensystem
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: PetscCall(EPSSolve(eps));
109: PetscCall(EPSGetType(eps,&type));
110: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type));
111: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
112: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
114: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115: Display solution and clean up
116: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118: /* show detailed info unless -terse option is given by user */
119: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
120: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
121: else {
122: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
123: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
124: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
125: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
126: }
127: if (isShell) PetscCall(STDestroy_Fold(fold));
128: PetscCall(EPSDestroy(&eps));
129: PetscCall(MatDestroy(&A));
130: PetscCall(SlepcFinalize());
131: return 0;
132: }
134: /*
135: STCreate_Fold - Creates the spectrum folding ST context.
137: Input Parameter:
138: + A - problem matrix
139: - target - target value
141: Output Parameter:
142: . fold - user-defined spectral transformation context
143: */
144: PetscErrorCode STCreate_Fold(Mat A,PetscScalar target,FoldShellST **fold)
145: {
146: FoldShellST *newctx;
148: PetscFunctionBeginUser;
149: PetscCall(PetscNew(&newctx));
150: newctx->A = A;
151: PetscCall(PetscObjectReference((PetscObject)A));
152: newctx->target = target;
153: PetscCall(MatCreateVecs(A,&newctx->w,NULL));
154: *fold = newctx;
155: PetscFunctionReturn(PETSC_SUCCESS);
156: }
158: /*
159: STApply_Fold - Applies the operator (A-target*I)^2 to a given vector.
161: Input Parameters:
162: + st - spectral transformation context
163: - x - input vector
165: Output Parameter:
166: . y - output vector
167: */
168: PetscErrorCode STApply_Fold(ST st,Vec x,Vec y)
169: {
170: FoldShellST *fold;
171: PetscScalar sigma;
173: PetscFunctionBeginUser;
174: PetscCall(STShellGetContext(st,&fold));
175: sigma = -fold->target;
176: PetscCall(MatMult(fold->A,x,fold->w));
177: PetscCall(VecAXPY(fold->w,sigma,x));
178: PetscCall(MatMult(fold->A,fold->w,y));
179: PetscCall(VecAXPY(y,sigma,fold->w));
180: PetscFunctionReturn(PETSC_SUCCESS);
181: }
183: /*
184: STDestroy_Fold - This routine destroys the shell ST context.
186: Input Parameter:
187: . fold - user-defined spectral transformation context
188: */
189: PetscErrorCode STDestroy_Fold(FoldShellST *fold)
190: {
191: PetscFunctionBeginUser;
192: PetscCall(MatDestroy(&fold->A));
193: PetscCall(VecDestroy(&fold->w));
194: PetscCall(PetscFree(fold));
195: PetscFunctionReturn(PETSC_SUCCESS);
196: }
198: /*TEST
200: test:
201: args: -m 11 -eps_nev 4 -terse
202: suffix: 1
203: requires: !single
205: TEST*/