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