Actual source code: ex9.c
slepc-3.18.3 2023-03-24
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[] = "Solves a problem associated to the Brusselator wave model in chemical reactions, illustrating the use of shell matrices.\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
14: " -L <L>, where <L> = bifurcation parameter.\n"
15: " -alpha <alpha>, -beta <beta>, -delta1 <delta1>, -delta2 <delta2>,\n"
16: " where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";
18: #include <slepceps.h>
20: /*
21: This example computes the eigenvalues with largest real part of the
22: following matrix
24: A = [ tau1*T+(beta-1)*I alpha^2*I
25: -beta*I tau2*T-alpha^2*I ],
27: where
29: T = tridiag{1,-2,1}
30: h = 1/(n+1)
31: tau1 = delta1/(h*L)^2
32: tau2 = delta2/(h*L)^2
33: */
35: /*
36: Matrix operations
37: */
38: PetscErrorCode MatMult_Brussel(Mat,Vec,Vec);
39: PetscErrorCode MatMultTranspose_Brussel(Mat,Vec,Vec);
40: PetscErrorCode MatGetDiagonal_Brussel(Mat,Vec);
42: typedef struct {
43: Mat T;
44: Vec x1,x2,y1,y2;
45: PetscScalar alpha,beta,tau1,tau2,sigma;
46: } CTX_BRUSSEL;
48: int main(int argc,char **argv)
49: {
50: Mat A; /* eigenvalue problem matrix */
51: EPS eps; /* eigenproblem solver context */
52: EPSType type;
53: PetscScalar delta1,delta2,L,h;
54: PetscInt N=30,n,i,Istart,Iend,nev;
55: CTX_BRUSSEL *ctx;
56: PetscBool terse;
57: PetscViewer viewer;
60: SlepcInitialize(&argc,&argv,(char*)0,help);
62: PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
63: PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%" PetscInt_FMT "\n\n",N);
65: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: Generate the matrix
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: /*
70: Create shell matrix context and set default parameters
71: */
72: PetscNew(&ctx);
73: ctx->alpha = 2.0;
74: ctx->beta = 5.45;
75: delta1 = 0.008;
76: delta2 = 0.004;
77: L = 0.51302;
79: /*
80: Look the command line for user-provided parameters
81: */
82: PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL);
83: PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL);
84: PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL);
85: PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL);
86: PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL);
88: /*
89: Create matrix T
90: */
91: MatCreate(PETSC_COMM_WORLD,&ctx->T);
92: MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N);
93: MatSetFromOptions(ctx->T);
94: MatSetUp(ctx->T);
96: MatGetOwnershipRange(ctx->T,&Istart,&Iend);
97: for (i=Istart;i<Iend;i++) {
98: if (i>0) MatSetValue(ctx->T,i,i-1,1.0,INSERT_VALUES);
99: if (i<N-1) MatSetValue(ctx->T,i,i+1,1.0,INSERT_VALUES);
100: MatSetValue(ctx->T,i,i,-2.0,INSERT_VALUES);
101: }
102: MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY);
103: MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY);
104: MatGetLocalSize(ctx->T,&n,NULL);
106: /*
107: Fill the remaining information in the shell matrix context
108: and create auxiliary vectors
109: */
110: h = 1.0 / (PetscReal)(N+1);
111: ctx->tau1 = delta1 / ((h*L)*(h*L));
112: ctx->tau2 = delta2 / ((h*L)*(h*L));
113: ctx->sigma = 0.0;
114: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1);
115: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2);
116: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1);
117: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2);
119: /*
120: Create the shell matrix
121: */
122: MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);
123: MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_Brussel);
124: MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Brussel);
125: MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Brussel);
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Create the eigensolver and set various options
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131: /*
132: Create eigensolver context
133: */
134: EPSCreate(PETSC_COMM_WORLD,&eps);
136: /*
137: Set operators. In this case, it is a standard eigenvalue problem
138: */
139: EPSSetOperators(eps,A,NULL);
140: EPSSetProblemType(eps,EPS_NHEP);
142: /*
143: Ask for the rightmost eigenvalues
144: */
145: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
147: /*
148: Set other solver options at runtime
149: */
150: EPSSetFromOptions(eps);
152: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153: Solve the eigensystem
154: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156: EPSSolve(eps);
158: /*
159: Optional: Get some information from the solver and display it
160: */
161: EPSGetType(eps,&type);
162: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
163: EPSGetDimensions(eps,&nev,NULL,NULL);
164: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
166: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167: Display solution and clean up
168: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170: /* show detailed info unless -terse option is given by user */
171: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
172: if (terse) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
173: else {
174: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
175: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
176: EPSConvergedReasonView(eps,viewer);
177: EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer);
178: PetscViewerPopFormat(viewer);
179: }
180: EPSDestroy(&eps);
181: MatDestroy(&A);
182: MatDestroy(&ctx->T);
183: VecDestroy(&ctx->x1);
184: VecDestroy(&ctx->x2);
185: VecDestroy(&ctx->y1);
186: VecDestroy(&ctx->y2);
187: PetscFree(ctx);
188: SlepcFinalize();
189: return 0;
190: }
192: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
193: {
194: PetscInt n;
195: const PetscScalar *px;
196: PetscScalar *py;
197: CTX_BRUSSEL *ctx;
200: MatShellGetContext(A,&ctx);
201: MatGetLocalSize(ctx->T,&n,NULL);
202: VecGetArrayRead(x,&px);
203: VecGetArray(y,&py);
204: VecPlaceArray(ctx->x1,px);
205: VecPlaceArray(ctx->x2,px+n);
206: VecPlaceArray(ctx->y1,py);
207: VecPlaceArray(ctx->y2,py+n);
209: MatMult(ctx->T,ctx->x1,ctx->y1);
210: VecScale(ctx->y1,ctx->tau1);
211: VecAXPY(ctx->y1,ctx->beta-1.0+ctx->sigma,ctx->x1);
212: VecAXPY(ctx->y1,ctx->alpha*ctx->alpha,ctx->x2);
214: MatMult(ctx->T,ctx->x2,ctx->y2);
215: VecScale(ctx->y2,ctx->tau2);
216: VecAXPY(ctx->y2,-ctx->beta,ctx->x1);
217: VecAXPY(ctx->y2,-ctx->alpha*ctx->alpha+ctx->sigma,ctx->x2);
219: VecRestoreArrayRead(x,&px);
220: VecRestoreArray(y,&py);
221: VecResetArray(ctx->x1);
222: VecResetArray(ctx->x2);
223: VecResetArray(ctx->y1);
224: VecResetArray(ctx->y2);
225: return 0;
226: }
228: PetscErrorCode MatMultTranspose_Brussel(Mat A,Vec x,Vec y)
229: {
230: PetscInt n;
231: const PetscScalar *px;
232: PetscScalar *py;
233: CTX_BRUSSEL *ctx;
236: MatShellGetContext(A,&ctx);
237: MatGetLocalSize(ctx->T,&n,NULL);
238: VecGetArrayRead(x,&px);
239: VecGetArray(y,&py);
240: VecPlaceArray(ctx->x1,px);
241: VecPlaceArray(ctx->x2,px+n);
242: VecPlaceArray(ctx->y1,py);
243: VecPlaceArray(ctx->y2,py+n);
245: MatMultTranspose(ctx->T,ctx->x1,ctx->y1);
246: VecScale(ctx->y1,ctx->tau1);
247: VecAXPY(ctx->y1,ctx->beta-1.0+ctx->sigma,ctx->x1);
248: VecAXPY(ctx->y1,-ctx->beta,ctx->x2);
250: MatMultTranspose(ctx->T,ctx->x2,ctx->y2);
251: VecScale(ctx->y2,ctx->tau2);
252: VecAXPY(ctx->y2,ctx->alpha*ctx->alpha,ctx->x1);
253: VecAXPY(ctx->y2,-ctx->alpha*ctx->alpha+ctx->sigma,ctx->x2);
255: VecRestoreArrayRead(x,&px);
256: VecRestoreArray(y,&py);
257: VecResetArray(ctx->x1);
258: VecResetArray(ctx->x2);
259: VecResetArray(ctx->y1);
260: VecResetArray(ctx->y2);
261: return 0;
262: }
264: PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
265: {
266: Vec d1,d2;
267: PetscInt n;
268: PetscScalar *pd;
269: MPI_Comm comm;
270: CTX_BRUSSEL *ctx;
273: MatShellGetContext(A,&ctx);
274: PetscObjectGetComm((PetscObject)A,&comm);
275: MatGetLocalSize(ctx->T,&n,NULL);
276: VecGetArray(diag,&pd);
277: VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1);
278: VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2);
280: VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma);
281: VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma);
283: VecDestroy(&d1);
284: VecDestroy(&d2);
285: VecRestoreArray(diag,&pd);
286: return 0;
287: }
289: /*TEST
291: test:
292: suffix: 1
293: args: -n 50 -eps_nev 4 -eps_two_sided {{0 1}} -eps_type {{krylovschur lapack}} -terse
294: requires: !single
295: filter: grep -v method
297: test:
298: suffix: 2
299: args: -eps_nev 8 -eps_max_it 300 -eps_target -28 -rg_type interval -rg_interval_endpoints -40,-20,-.1,.1 -terse
300: requires: !single
302: test:
303: suffix: 3
304: args: -n 50 -eps_nev 4 -eps_balance twoside -terse
305: requires: double
306: filter: grep -v method
307: output_file: output/ex9_1.out
309: test:
310: suffix: 4
311: args: -eps_smallest_imaginary -eps_ncv 24 -terse
312: requires: !complex !single
314: test:
315: suffix: 4_complex
316: args: -eps_smallest_imaginary -eps_ncv 24 -terse
317: requires: complex !single
319: test:
320: suffix: 5
321: args: -eps_nev 4 -eps_target_real -eps_target -3 -terse
322: requires: !single
324: test:
325: suffix: 6
326: args: -eps_nev 2 -eps_target_imaginary -eps_target 3i -terse
327: requires: complex !single
329: test:
330: suffix: 7
331: args: -n 40 -eps_nev 1 -eps_type arnoldi -eps_smallest_real -eps_refined -eps_ncv 40 -eps_max_it 300 -terse
332: requires: double
334: test:
335: suffix: 8
336: args: -eps_nev 2 -eps_target -30 -eps_type jd -st_matmode shell -eps_jd_fix 0.0001 -eps_jd_const_correction_tol 0 -terse
337: requires: !single
338: filter: sed -e "s/[+-]0\.0*i//g"
339: timeoutfactor: 2
341: test:
342: suffix: 9
343: args: -eps_largest_imaginary -eps_ncv 24 -terse
344: requires: !single
346: TEST*/