Actual source code: ex35.c
slepc-main 2024-11-15
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,NULL,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));
60: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
61: for (II=Istart;II<Iend;II++) {
62: i = II/n; j = II-i*n;
63: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
64: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
65: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
66: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
67: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
68: }
69: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
70: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
72: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73: Create the eigensolver and set various options
74: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
77: PetscCall(EPSSetOperators(eps,A,NULL));
78: PetscCall(EPSSetProblemType(eps,EPS_HEP));
79: PetscCall(EPSSetTarget(eps,target));
80: PetscCall(EPSGetST(eps,&st));
81: PetscCall(STSetType(st,STSHELL));
82: PetscCall(EPSSetFromOptions(eps));
84: /*
85: Initialize shell spectral transformation
86: */
87: PetscCall(PetscObjectTypeCompare((PetscObject)st,STSHELL,&isShell));
88: if (isShell) {
89: /* Change sorting criterion since this shell ST computes eigenvalues
90: of the transformed operator closest to 0 */
91: PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
93: /* Create the context for the user-defined spectral transform */
94: PetscCall(STCreate_Fold(A,target,&fold));
95: PetscCall(STShellSetContext(st,fold));
97: /* Set callback function for applying the operator (in this case we do not
98: provide a back-transformation callback since the mapping is not one-to-one) */
99: PetscCall(STShellSetApply(st,STApply_Fold));
100: PetscCall(PetscObjectSetName((PetscObject)st,"STFOLD"));
101: }
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Solve the eigensystem
105: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: PetscCall(EPSSolve(eps));
108: PetscCall(EPSGetType(eps,&type));
109: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type));
110: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
111: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
113: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114: Display solution and clean up
115: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117: /* show detailed info unless -terse option is given by user */
118: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
119: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
120: else {
121: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
122: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
123: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
124: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
125: }
126: if (isShell) PetscCall(STDestroy_Fold(fold));
127: PetscCall(EPSDestroy(&eps));
128: PetscCall(MatDestroy(&A));
129: PetscCall(SlepcFinalize());
130: return 0;
131: }
133: /*
134: STCreate_Fold - Creates the spectrum folding ST context.
136: Input Parameter:
137: + A - problem matrix
138: - target - target value
140: Output Parameter:
141: . fold - user-defined spectral transformation context
142: */
143: PetscErrorCode STCreate_Fold(Mat A,PetscScalar target,FoldShellST **fold)
144: {
145: FoldShellST *newctx;
147: PetscFunctionBeginUser;
148: PetscCall(PetscNew(&newctx));
149: newctx->A = A;
150: PetscCall(PetscObjectReference((PetscObject)A));
151: newctx->target = target;
152: PetscCall(MatCreateVecs(A,&newctx->w,NULL));
153: *fold = newctx;
154: PetscFunctionReturn(PETSC_SUCCESS);
155: }
157: /*
158: STApply_Fold - Applies the operator (A-target*I)^2 to a given vector.
160: Input Parameters:
161: + st - spectral transformation context
162: - x - input vector
164: Output Parameter:
165: . y - output vector
166: */
167: PetscErrorCode STApply_Fold(ST st,Vec x,Vec y)
168: {
169: FoldShellST *fold;
170: PetscScalar sigma;
172: PetscFunctionBeginUser;
173: PetscCall(STShellGetContext(st,&fold));
174: sigma = -fold->target;
175: PetscCall(MatMult(fold->A,x,fold->w));
176: PetscCall(VecAXPY(fold->w,sigma,x));
177: PetscCall(MatMult(fold->A,fold->w,y));
178: PetscCall(VecAXPY(y,sigma,fold->w));
179: PetscFunctionReturn(PETSC_SUCCESS);
180: }
182: /*
183: STDestroy_Fold - This routine destroys the shell ST context.
185: Input Parameter:
186: . fold - user-defined spectral transformation context
187: */
188: PetscErrorCode STDestroy_Fold(FoldShellST *fold)
189: {
190: PetscFunctionBeginUser;
191: PetscCall(MatDestroy(&fold->A));
192: PetscCall(VecDestroy(&fold->w));
193: PetscCall(PetscFree(fold));
194: PetscFunctionReturn(PETSC_SUCCESS);
195: }
197: /*TEST
199: test:
200: args: -m 11 -eps_nev 4 -terse
201: suffix: 1
202: requires: !single
204: TEST*/