Actual source code: ex30.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[] = "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 MatGetDiagonal_Brussel(Mat,Vec);
32: PetscErrorCode MyStoppingTest(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*);
34: typedef struct {
35: Mat T;
36: Vec x1,x2,y1,y2;
37: PetscScalar alpha,beta,tau1,tau2,sigma;
38: PetscInt lastnconv; /* last value of nconv; used in stopping test */
39: PetscInt nreps; /* number of repetitions of nconv; used in stopping test */
40: } CTX_BRUSSEL;
42: int main(int argc,char **argv)
43: {
44: Mat A; /* eigenvalue problem matrix */
45: EPS eps; /* eigenproblem solver context */
46: RG rg; /* region object */
47: PetscScalar delta1,delta2,L,h;
48: PetscInt N=30,n,i,Istart,Iend,mpd;
49: CTX_BRUSSEL *ctx;
50: PetscBool terse;
51: PetscViewer viewer;
53: PetscFunctionBeginUser;
54: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
56: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
57: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%" PetscInt_FMT "\n\n",N));
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Generate the matrix
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: /*
64: Create shell matrix context and set default parameters
65: */
66: PetscCall(PetscNew(&ctx));
67: ctx->alpha = 2.0;
68: ctx->beta = 5.45;
69: delta1 = 0.008;
70: delta2 = 0.004;
71: L = 0.51302;
73: /*
74: Look the command line for user-provided parameters
75: */
76: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL));
77: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL));
78: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL));
79: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL));
80: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL));
82: /*
83: Create matrix T
84: */
85: PetscCall(MatCreate(PETSC_COMM_WORLD,&ctx->T));
86: PetscCall(MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N));
87: PetscCall(MatSetFromOptions(ctx->T));
89: PetscCall(MatGetOwnershipRange(ctx->T,&Istart,&Iend));
90: for (i=Istart;i<Iend;i++) {
91: if (i>0) PetscCall(MatSetValue(ctx->T,i,i-1,1.0,INSERT_VALUES));
92: if (i<N-1) PetscCall(MatSetValue(ctx->T,i,i+1,1.0,INSERT_VALUES));
93: PetscCall(MatSetValue(ctx->T,i,i,-2.0,INSERT_VALUES));
94: }
95: PetscCall(MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY));
96: PetscCall(MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY));
97: PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
99: /*
100: Fill the remaining information in the shell matrix context
101: and create auxiliary vectors
102: */
103: h = 1.0 / (PetscReal)(N+1);
104: ctx->tau1 = delta1 / ((h*L)*(h*L));
105: ctx->tau2 = delta2 / ((h*L)*(h*L));
106: ctx->sigma = 0.0;
107: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1));
108: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2));
109: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1));
110: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2));
112: /*
113: Create the shell matrix
114: */
115: PetscCall(MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A));
116: PetscCall(MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_Brussel));
117: PetscCall(MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Brussel));
119: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: Create the eigensolver and configure the region
121: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
124: PetscCall(EPSSetOperators(eps,A,NULL));
125: PetscCall(EPSSetProblemType(eps,EPS_NHEP));
127: /*
128: Define the region containing the eigenvalues of interest
129: */
130: PetscCall(EPSGetRG(eps,&rg));
131: PetscCall(RGSetType(rg,RGINTERVAL));
132: PetscCall(RGIntervalSetEndpoints(rg,-PETSC_INFINITY,PETSC_INFINITY,-0.01,0.01));
133: PetscCall(RGSetComplement(rg,PETSC_TRUE));
134: /* sort eigenvalue approximations wrt a target, otherwise convergence will be erratic */
135: PetscCall(EPSSetTarget(eps,0.0));
136: PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
138: /*
139: Set solver options. In particular, we must allocate sufficient
140: storage for all eigenpairs that may converge (ncv). This is
141: application-dependent.
142: */
143: mpd = 40;
144: PetscCall(EPSSetDimensions(eps,2*mpd,3*mpd,mpd));
145: PetscCall(EPSSetTolerances(eps,1e-7,2000));
146: ctx->lastnconv = 0;
147: ctx->nreps = 0;
148: PetscCall(EPSSetStoppingTestFunction(eps,MyStoppingTest,(void*)ctx,NULL));
149: PetscCall(EPSSetFromOptions(eps));
151: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152: Solve the eigensystem and display solution
153: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155: PetscCall(EPSSolve(eps));
157: /* show detailed info unless -terse option is given by user */
158: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
159: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
160: PetscCall(EPSConvergedReasonView(eps,viewer));
161: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
162: if (!terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer));
163: PetscCall(PetscViewerPopFormat(viewer));
165: PetscCall(EPSDestroy(&eps));
166: PetscCall(MatDestroy(&A));
167: PetscCall(MatDestroy(&ctx->T));
168: PetscCall(VecDestroy(&ctx->x1));
169: PetscCall(VecDestroy(&ctx->x2));
170: PetscCall(VecDestroy(&ctx->y1));
171: PetscCall(VecDestroy(&ctx->y2));
172: PetscCall(PetscFree(ctx));
173: PetscCall(SlepcFinalize());
174: return 0;
175: }
177: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
178: {
179: PetscInt n;
180: const PetscScalar *px;
181: PetscScalar *py;
182: CTX_BRUSSEL *ctx;
184: PetscFunctionBeginUser;
185: PetscCall(MatShellGetContext(A,&ctx));
186: PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
187: PetscCall(VecGetArrayRead(x,&px));
188: PetscCall(VecGetArray(y,&py));
189: PetscCall(VecPlaceArray(ctx->x1,px));
190: PetscCall(VecPlaceArray(ctx->x2,px+n));
191: PetscCall(VecPlaceArray(ctx->y1,py));
192: PetscCall(VecPlaceArray(ctx->y2,py+n));
194: PetscCall(MatMult(ctx->T,ctx->x1,ctx->y1));
195: PetscCall(VecScale(ctx->y1,ctx->tau1));
196: PetscCall(VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1));
197: PetscCall(VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2));
199: PetscCall(MatMult(ctx->T,ctx->x2,ctx->y2));
200: PetscCall(VecScale(ctx->y2,ctx->tau2));
201: PetscCall(VecAXPY(ctx->y2,-ctx->beta,ctx->x1));
202: PetscCall(VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2));
204: PetscCall(VecRestoreArrayRead(x,&px));
205: PetscCall(VecRestoreArray(y,&py));
206: PetscCall(VecResetArray(ctx->x1));
207: PetscCall(VecResetArray(ctx->x2));
208: PetscCall(VecResetArray(ctx->y1));
209: PetscCall(VecResetArray(ctx->y2));
210: PetscFunctionReturn(PETSC_SUCCESS);
211: }
213: PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
214: {
215: Vec d1,d2;
216: PetscInt n;
217: PetscScalar *pd;
218: MPI_Comm comm;
219: CTX_BRUSSEL *ctx;
221: PetscFunctionBeginUser;
222: PetscCall(MatShellGetContext(A,&ctx));
223: PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
224: PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
225: PetscCall(VecGetArray(diag,&pd));
226: PetscCall(VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1));
227: PetscCall(VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2));
229: PetscCall(VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma));
230: PetscCall(VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma));
232: PetscCall(VecDestroy(&d1));
233: PetscCall(VecDestroy(&d2));
234: PetscCall(VecRestoreArray(diag,&pd));
235: PetscFunctionReturn(PETSC_SUCCESS);
236: }
238: /*
239: Function for user-defined stopping test.
241: Ignores the value of nev. It only takes into account the number of
242: eigenpairs that have converged in recent outer iterations (restarts);
243: if no new eigenvalues have converged in the last few restarts,
244: we stop the iteration, assuming that no more eigenvalues are present
245: inside the region.
246: */
247: PetscErrorCode MyStoppingTest(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ptr)
248: {
249: CTX_BRUSSEL *ctx = (CTX_BRUSSEL*)ptr;
251: PetscFunctionBeginUser;
252: /* check usual termination conditions, but ignoring the case nconv>=nev */
253: PetscCall(EPSStoppingBasic(eps,its,max_it,nconv,PETSC_INT_MAX,reason,NULL));
254: if (*reason==EPS_CONVERGED_ITERATING) {
255: /* check if nconv is the same as before */
256: if (nconv==ctx->lastnconv) ctx->nreps++;
257: else {
258: ctx->lastnconv = nconv;
259: ctx->nreps = 0;
260: }
261: /* check if no eigenvalues converged in last 10 restarts */
262: if (nconv && ctx->nreps>10) *reason = EPS_CONVERGED_USER;
263: }
264: PetscFunctionReturn(PETSC_SUCCESS);
265: }
267: /*TEST
269: test:
270: suffix: 1
271: args: -n 100 -terse
272: requires: !single
274: TEST*/