Actual source code: ex30.c
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[] = "Illustrates the use of a region for filtering; the number of wanted eigenvalues is not known a priori.\n\n"
12: "The problem is the Brusselator wave model as in ex9.c.\n"
13: "The command line options are:\n"
14: " -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
15: " -L <L>, where <L> = bifurcation parameter.\n"
16: " -alpha <alpha>, -beta <beta>, -delta1 <delta1>, -delta2 <delta2>,\n"
17: " where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";
19: #include <slepceps.h>
21: /*
22: This example tries to compute all eigenvalues lying outside the real axis.
23: This could be achieved by computing LARGEST_IMAGINARY eigenvalues, but
24: here we take a different route: define a region of the complex plane where
25: eigenvalues must be emphasized (eigenvalues outside the region are filtered
26: out). In this case, we select the region as the complement of a thin stripe
27: around the real axis.
28: */
30: PetscErrorCode MatMult_Brussel(Mat,Vec,Vec);
31: PetscErrorCode MyStoppingTest(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*);
33: typedef struct {
34: Mat T;
35: Vec x1,x2,y1,y2;
36: PetscScalar alpha,beta,tau1,tau2,sigma;
37: PetscInt lastnconv; /* last value of nconv; used in stopping test */
38: PetscInt nreps; /* number of repetitions of nconv; used in stopping test */
39: } CTX_BRUSSEL;
41: int main(int argc,char **argv)
42: {
43: Mat A; /* eigenvalue problem matrix */
44: EPS eps; /* eigenproblem solver context */
45: RG rg; /* region object */
46: PetscScalar delta1,delta2,L,h;
47: PetscInt N=30,n,i,Istart,Iend,mpd;
48: CTX_BRUSSEL *ctx;
49: PetscBool terse;
50: PetscViewer viewer;
52: PetscFunctionBeginUser;
53: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
55: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
56: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%" PetscInt_FMT "\n\n",N));
58: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59: Generate the matrix
60: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
62: /*
63: Create shell matrix context and set default parameters
64: */
65: PetscCall(PetscNew(&ctx));
66: ctx->alpha = 2.0;
67: ctx->beta = 5.45;
68: delta1 = 0.008;
69: delta2 = 0.004;
70: L = 0.51302;
72: /*
73: Look the command line for user-provided parameters
74: */
75: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL));
76: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL));
77: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL));
78: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL));
79: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL));
81: /*
82: Create matrix T
83: */
84: PetscCall(MatCreate(PETSC_COMM_WORLD,&ctx->T));
85: PetscCall(MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N));
86: PetscCall(MatSetFromOptions(ctx->T));
88: PetscCall(MatGetOwnershipRange(ctx->T,&Istart,&Iend));
89: for (i=Istart;i<Iend;i++) {
90: if (i>0) PetscCall(MatSetValue(ctx->T,i,i-1,1.0,INSERT_VALUES));
91: if (i<N-1) PetscCall(MatSetValue(ctx->T,i,i+1,1.0,INSERT_VALUES));
92: PetscCall(MatSetValue(ctx->T,i,i,-2.0,INSERT_VALUES));
93: }
94: PetscCall(MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY));
95: PetscCall(MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY));
96: PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
98: /*
99: Fill the remaining information in the shell matrix context
100: and create auxiliary vectors
101: */
102: h = 1.0 / (PetscReal)(N+1);
103: ctx->tau1 = delta1 / ((h*L)*(h*L));
104: ctx->tau2 = delta2 / ((h*L)*(h*L));
105: ctx->sigma = 0.0;
106: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1));
107: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2));
108: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1));
109: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2));
111: /*
112: Create the shell matrix
113: */
114: PetscCall(MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A));
115: PetscCall(MatShellSetOperation(A,MATOP_MULT,(PetscErrorCodeFn*)MatMult_Brussel));
117: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118: Create the eigensolver and configure the region
119: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
121: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
122: PetscCall(EPSSetOperators(eps,A,NULL));
123: PetscCall(EPSSetProblemType(eps,EPS_NHEP));
125: /*
126: Define the region containing the eigenvalues of interest
127: */
128: PetscCall(EPSGetRG(eps,&rg));
129: PetscCall(RGSetType(rg,RGINTERVAL));
130: PetscCall(RGIntervalSetEndpoints(rg,-PETSC_INFINITY,PETSC_INFINITY,-0.01,0.01));
131: PetscCall(RGSetComplement(rg,PETSC_TRUE));
132: /* sort eigenvalue approximations wrt a target, otherwise convergence will be erratic */
133: PetscCall(EPSSetTarget(eps,0.0));
134: PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
136: /*
137: Set solver options. In particular, we must allocate sufficient
138: storage for all eigenpairs that may converge (ncv). This is
139: application-dependent.
140: */
141: mpd = 40;
142: PetscCall(EPSSetDimensions(eps,2*mpd,3*mpd,mpd));
143: PetscCall(EPSSetTolerances(eps,1e-7,2000));
144: ctx->lastnconv = 0;
145: ctx->nreps = 0;
146: PetscCall(EPSSetStoppingTestFunction(eps,MyStoppingTest,(void*)ctx,NULL));
147: PetscCall(EPSSetFromOptions(eps));
149: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Solve the eigensystem and display solution
151: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153: PetscCall(EPSSolve(eps));
155: /* show detailed info unless -terse option is given by user */
156: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
157: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
158: PetscCall(EPSConvergedReasonView(eps,viewer));
159: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
160: if (!terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer));
161: PetscCall(PetscViewerPopFormat(viewer));
163: PetscCall(EPSDestroy(&eps));
164: PetscCall(MatDestroy(&A));
165: PetscCall(MatDestroy(&ctx->T));
166: PetscCall(VecDestroy(&ctx->x1));
167: PetscCall(VecDestroy(&ctx->x2));
168: PetscCall(VecDestroy(&ctx->y1));
169: PetscCall(VecDestroy(&ctx->y2));
170: PetscCall(PetscFree(ctx));
171: PetscCall(SlepcFinalize());
172: return 0;
173: }
175: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
176: {
177: PetscInt n;
178: const PetscScalar *px;
179: PetscScalar *py;
180: CTX_BRUSSEL *ctx;
182: PetscFunctionBeginUser;
183: PetscCall(MatShellGetContext(A,&ctx));
184: PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
185: PetscCall(VecGetArrayRead(x,&px));
186: PetscCall(VecGetArray(y,&py));
187: PetscCall(VecPlaceArray(ctx->x1,px));
188: PetscCall(VecPlaceArray(ctx->x2,px+n));
189: PetscCall(VecPlaceArray(ctx->y1,py));
190: PetscCall(VecPlaceArray(ctx->y2,py+n));
192: PetscCall(MatMult(ctx->T,ctx->x1,ctx->y1));
193: PetscCall(VecScale(ctx->y1,ctx->tau1));
194: PetscCall(VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1));
195: PetscCall(VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2));
197: PetscCall(MatMult(ctx->T,ctx->x2,ctx->y2));
198: PetscCall(VecScale(ctx->y2,ctx->tau2));
199: PetscCall(VecAXPY(ctx->y2,-ctx->beta,ctx->x1));
200: PetscCall(VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2));
202: PetscCall(VecRestoreArrayRead(x,&px));
203: PetscCall(VecRestoreArray(y,&py));
204: PetscCall(VecResetArray(ctx->x1));
205: PetscCall(VecResetArray(ctx->x2));
206: PetscCall(VecResetArray(ctx->y1));
207: PetscCall(VecResetArray(ctx->y2));
208: PetscFunctionReturn(PETSC_SUCCESS);
209: }
211: /*
212: Function for user-defined stopping test.
214: Ignores the value of nev. It only takes into account the number of
215: eigenpairs that have converged in recent outer iterations (restarts);
216: if no new eigenvalues have converged in the last few restarts,
217: we stop the iteration, assuming that no more eigenvalues are present
218: inside the region.
219: */
220: PetscErrorCode MyStoppingTest(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ptr)
221: {
222: CTX_BRUSSEL *ctx = (CTX_BRUSSEL*)ptr;
224: PetscFunctionBeginUser;
225: /* check usual termination conditions, but ignoring the case nconv>=nev */
226: PetscCall(EPSStoppingBasic(eps,its,max_it,nconv,PETSC_INT_MAX,reason,NULL));
227: if (*reason==EPS_CONVERGED_ITERATING) {
228: /* check if nconv is the same as before */
229: if (nconv==ctx->lastnconv) ctx->nreps++;
230: else {
231: ctx->lastnconv = nconv;
232: ctx->nreps = 0;
233: }
234: /* check if no eigenvalues converged in last 10 restarts */
235: if (nconv && ctx->nreps>10) *reason = EPS_CONVERGED_USER;
236: }
237: PetscFunctionReturn(PETSC_SUCCESS);
238: }
240: /*TEST
242: test:
243: suffix: 1
244: args: -n 100 -terse
245: requires: !single
247: TEST*/