Actual source code: ex9.c

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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);if (ierr) return ierr;

 62:   PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
 63:   PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%D\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: %D\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) {
173:     EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
174:   } else {
175:     PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
176:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
177:     EPSConvergedReasonView(eps,viewer);
178:     EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer);
179:     PetscViewerPopFormat(viewer);
180:   }
181:   EPSDestroy(&eps);
182:   MatDestroy(&A);
183:   MatDestroy(&ctx->T);
184:   VecDestroy(&ctx->x1);
185:   VecDestroy(&ctx->x2);
186:   VecDestroy(&ctx->y1);
187:   VecDestroy(&ctx->y2);
188:   PetscFree(ctx);
189:   SlepcFinalize();
190:   return ierr;
191: }

193: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
194: {
195:   PetscInt          n;
196:   const PetscScalar *px;
197:   PetscScalar       *py;
198:   CTX_BRUSSEL       *ctx;
199:   PetscErrorCode    ierr;

202:   MatShellGetContext(A,&ctx);
203:   MatGetLocalSize(ctx->T,&n,NULL);
204:   VecGetArrayRead(x,&px);
205:   VecGetArray(y,&py);
206:   VecPlaceArray(ctx->x1,px);
207:   VecPlaceArray(ctx->x2,px+n);
208:   VecPlaceArray(ctx->y1,py);
209:   VecPlaceArray(ctx->y2,py+n);

211:   MatMult(ctx->T,ctx->x1,ctx->y1);
212:   VecScale(ctx->y1,ctx->tau1);
213:   VecAXPY(ctx->y1,ctx->beta-1.0+ctx->sigma,ctx->x1);
214:   VecAXPY(ctx->y1,ctx->alpha*ctx->alpha,ctx->x2);

216:   MatMult(ctx->T,ctx->x2,ctx->y2);
217:   VecScale(ctx->y2,ctx->tau2);
218:   VecAXPY(ctx->y2,-ctx->beta,ctx->x1);
219:   VecAXPY(ctx->y2,-ctx->alpha*ctx->alpha+ctx->sigma,ctx->x2);

221:   VecRestoreArrayRead(x,&px);
222:   VecRestoreArray(y,&py);
223:   VecResetArray(ctx->x1);
224:   VecResetArray(ctx->x2);
225:   VecResetArray(ctx->y1);
226:   VecResetArray(ctx->y2);
227:   return(0);
228: }

230: PetscErrorCode MatMultTranspose_Brussel(Mat A,Vec x,Vec y)
231: {
232:   PetscInt          n;
233:   const PetscScalar *px;
234:   PetscScalar       *py;
235:   CTX_BRUSSEL       *ctx;
236:   PetscErrorCode    ierr;

239:   MatShellGetContext(A,&ctx);
240:   MatGetLocalSize(ctx->T,&n,NULL);
241:   VecGetArrayRead(x,&px);
242:   VecGetArray(y,&py);
243:   VecPlaceArray(ctx->x1,px);
244:   VecPlaceArray(ctx->x2,px+n);
245:   VecPlaceArray(ctx->y1,py);
246:   VecPlaceArray(ctx->y2,py+n);

248:   MatMultTranspose(ctx->T,ctx->x1,ctx->y1);
249:   VecScale(ctx->y1,ctx->tau1);
250:   VecAXPY(ctx->y1,ctx->beta-1.0+ctx->sigma,ctx->x1);
251:   VecAXPY(ctx->y1,-ctx->beta,ctx->x2);

253:   MatMultTranspose(ctx->T,ctx->x2,ctx->y2);
254:   VecScale(ctx->y2,ctx->tau2);
255:   VecAXPY(ctx->y2,ctx->alpha*ctx->alpha,ctx->x1);
256:   VecAXPY(ctx->y2,-ctx->alpha*ctx->alpha+ctx->sigma,ctx->x2);

258:   VecRestoreArrayRead(x,&px);
259:   VecRestoreArray(y,&py);
260:   VecResetArray(ctx->x1);
261:   VecResetArray(ctx->x2);
262:   VecResetArray(ctx->y1);
263:   VecResetArray(ctx->y2);
264:   return(0);
265: }

267: PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
268: {
269:   Vec            d1,d2;
270:   PetscInt       n;
271:   PetscScalar    *pd;
272:   MPI_Comm       comm;
273:   CTX_BRUSSEL    *ctx;

277:   MatShellGetContext(A,&ctx);
278:   PetscObjectGetComm((PetscObject)A,&comm);
279:   MatGetLocalSize(ctx->T,&n,NULL);
280:   VecGetArray(diag,&pd);
281:   VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1);
282:   VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2);

284:   VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma);
285:   VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma);

287:   VecDestroy(&d1);
288:   VecDestroy(&d2);
289:   VecRestoreArray(diag,&pd);
290:   return(0);
291: }

293: /*TEST

295:    test:
296:       suffix: 1
297:       args: -n 50 -eps_nev 4 -eps_two_sided {{0 1}} -eps_type {{krylovschur lapack}} -terse
298:       requires: !single
299:       filter: grep -v method

301:    test:
302:       suffix: 2
303:       args: -eps_nev 8 -eps_max_it 300 -eps_target -28 -rg_type interval -rg_interval_endpoints -40,-20,-.1,.1 -terse
304:       requires: !single

306:    test:
307:       suffix: 3
308:       args: -n 50 -eps_nev 4 -eps_balance twoside -terse
309:       requires: double
310:       filter: grep -v method
311:       output_file: output/ex9_1.out

313:    test:
314:       suffix: 4
315:       args: -eps_smallest_imaginary -eps_ncv 24 -terse
316:       requires: !complex !single

318:    test:
319:       suffix: 4_complex
320:       args: -eps_smallest_imaginary -eps_ncv 24 -terse
321:       requires: complex !single

323:    test:
324:       suffix: 5
325:       args: -eps_nev 4 -eps_target_real -eps_target -3 -terse
326:       requires: !single

328:    test:
329:       suffix: 6
330:       args: -eps_nev 2 -eps_target_imaginary -eps_target 3i -terse
331:       requires: complex !single

333:    test:
334:       suffix: 7
335:       args: -n 40 -eps_nev 1 -eps_type arnoldi -eps_smallest_real -eps_refined -eps_ncv 40 -eps_max_it 300 -terse
336:       requires: double

338:    test:
339:       suffix: 8
340:       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
341:       requires: !single
342:       filter: sed -e "s/[+-]0\.0*i//g"

344:    test:
345:       suffix: 9
346:       args: -eps_largest_imaginary -eps_ncv 24 -terse
347:       requires: !single

349: TEST*/