Actual source code: ex9.c

slepc-main 2024-11-15
Report Typos and Errors
  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;

 59:   PetscFunctionBeginUser;
 60:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

 62:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
 63:   PetscCall(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:   PetscCall(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:   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL));
 83:   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL));
 84:   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL));
 85:   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL));
 86:   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL));

 88:   /*
 89:      Create matrix T
 90:   */
 91:   PetscCall(MatCreate(PETSC_COMM_WORLD,&ctx->T));
 92:   PetscCall(MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N));
 93:   PetscCall(MatSetFromOptions(ctx->T));

 95:   PetscCall(MatGetOwnershipRange(ctx->T,&Istart,&Iend));
 96:   for (i=Istart;i<Iend;i++) {
 97:     if (i>0) PetscCall(MatSetValue(ctx->T,i,i-1,1.0,INSERT_VALUES));
 98:     if (i<N-1) PetscCall(MatSetValue(ctx->T,i,i+1,1.0,INSERT_VALUES));
 99:     PetscCall(MatSetValue(ctx->T,i,i,-2.0,INSERT_VALUES));
100:   }
101:   PetscCall(MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY));
102:   PetscCall(MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY));
103:   PetscCall(MatGetLocalSize(ctx->T,&n,NULL));

105:   /*
106:      Fill the remaining information in the shell matrix context
107:      and create auxiliary vectors
108:   */
109:   h = 1.0 / (PetscReal)(N+1);
110:   ctx->tau1 = delta1 / ((h*L)*(h*L));
111:   ctx->tau2 = delta2 / ((h*L)*(h*L));
112:   ctx->sigma = 0.0;
113:   PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1));
114:   PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2));
115:   PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1));
116:   PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2));

118:   /*
119:      Create the shell matrix
120:   */
121:   PetscCall(MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A));
122:   PetscCall(MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_Brussel));
123:   PetscCall(MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Brussel));
124:   PetscCall(MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Brussel));

126:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127:                 Create the eigensolver and set various options
128:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

130:   /*
131:      Create eigensolver context
132:   */
133:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));

135:   /*
136:      Set operators. In this case, it is a standard eigenvalue problem
137:   */
138:   PetscCall(EPSSetOperators(eps,A,NULL));
139:   PetscCall(EPSSetProblemType(eps,EPS_NHEP));

141:   /*
142:      Ask for the rightmost eigenvalues
143:   */
144:   PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));

146:   /*
147:      Set other solver options at runtime
148:   */
149:   PetscCall(EPSSetFromOptions(eps));

151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152:                       Solve the eigensystem
153:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

155:   PetscCall(EPSSolve(eps));

157:   /*
158:      Optional: Get some information from the solver and display it
159:   */
160:   PetscCall(EPSGetType(eps,&type));
161:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type));
162:   PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
163:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));

165:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166:                     Display solution and clean up
167:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

169:   /* show detailed info unless -terse option is given by user */
170:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
171:   if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
172:   else {
173:     PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
174:     PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
175:     PetscCall(EPSConvergedReasonView(eps,viewer));
176:     PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer));
177:     PetscCall(PetscViewerPopFormat(viewer));
178:   }
179:   PetscCall(EPSDestroy(&eps));
180:   PetscCall(MatDestroy(&A));
181:   PetscCall(MatDestroy(&ctx->T));
182:   PetscCall(VecDestroy(&ctx->x1));
183:   PetscCall(VecDestroy(&ctx->x2));
184:   PetscCall(VecDestroy(&ctx->y1));
185:   PetscCall(VecDestroy(&ctx->y2));
186:   PetscCall(PetscFree(ctx));
187:   PetscCall(SlepcFinalize());
188:   return 0;
189: }

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

198:   PetscFunctionBeginUser;
199:   PetscCall(MatShellGetContext(A,&ctx));
200:   PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
201:   PetscCall(VecGetArrayRead(x,&px));
202:   PetscCall(VecGetArray(y,&py));
203:   PetscCall(VecPlaceArray(ctx->x1,px));
204:   PetscCall(VecPlaceArray(ctx->x2,px+n));
205:   PetscCall(VecPlaceArray(ctx->y1,py));
206:   PetscCall(VecPlaceArray(ctx->y2,py+n));

208:   PetscCall(MatMult(ctx->T,ctx->x1,ctx->y1));
209:   PetscCall(VecScale(ctx->y1,ctx->tau1));
210:   PetscCall(VecAXPY(ctx->y1,ctx->beta-1.0+ctx->sigma,ctx->x1));
211:   PetscCall(VecAXPY(ctx->y1,ctx->alpha*ctx->alpha,ctx->x2));

213:   PetscCall(MatMult(ctx->T,ctx->x2,ctx->y2));
214:   PetscCall(VecScale(ctx->y2,ctx->tau2));
215:   PetscCall(VecAXPY(ctx->y2,-ctx->beta,ctx->x1));
216:   PetscCall(VecAXPY(ctx->y2,-ctx->alpha*ctx->alpha+ctx->sigma,ctx->x2));

218:   PetscCall(VecRestoreArrayRead(x,&px));
219:   PetscCall(VecRestoreArray(y,&py));
220:   PetscCall(VecResetArray(ctx->x1));
221:   PetscCall(VecResetArray(ctx->x2));
222:   PetscCall(VecResetArray(ctx->y1));
223:   PetscCall(VecResetArray(ctx->y2));
224:   PetscFunctionReturn(PETSC_SUCCESS);
225: }

227: PetscErrorCode MatMultTranspose_Brussel(Mat A,Vec x,Vec y)
228: {
229:   PetscInt          n;
230:   const PetscScalar *px;
231:   PetscScalar       *py;
232:   CTX_BRUSSEL       *ctx;

234:   PetscFunctionBeginUser;
235:   PetscCall(MatShellGetContext(A,&ctx));
236:   PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
237:   PetscCall(VecGetArrayRead(x,&px));
238:   PetscCall(VecGetArray(y,&py));
239:   PetscCall(VecPlaceArray(ctx->x1,px));
240:   PetscCall(VecPlaceArray(ctx->x2,px+n));
241:   PetscCall(VecPlaceArray(ctx->y1,py));
242:   PetscCall(VecPlaceArray(ctx->y2,py+n));

244:   PetscCall(MatMultTranspose(ctx->T,ctx->x1,ctx->y1));
245:   PetscCall(VecScale(ctx->y1,ctx->tau1));
246:   PetscCall(VecAXPY(ctx->y1,ctx->beta-1.0+ctx->sigma,ctx->x1));
247:   PetscCall(VecAXPY(ctx->y1,-ctx->beta,ctx->x2));

249:   PetscCall(MatMultTranspose(ctx->T,ctx->x2,ctx->y2));
250:   PetscCall(VecScale(ctx->y2,ctx->tau2));
251:   PetscCall(VecAXPY(ctx->y2,ctx->alpha*ctx->alpha,ctx->x1));
252:   PetscCall(VecAXPY(ctx->y2,-ctx->alpha*ctx->alpha+ctx->sigma,ctx->x2));

254:   PetscCall(VecRestoreArrayRead(x,&px));
255:   PetscCall(VecRestoreArray(y,&py));
256:   PetscCall(VecResetArray(ctx->x1));
257:   PetscCall(VecResetArray(ctx->x2));
258:   PetscCall(VecResetArray(ctx->y1));
259:   PetscCall(VecResetArray(ctx->y2));
260:   PetscFunctionReturn(PETSC_SUCCESS);
261: }

263: PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
264: {
265:   Vec            d1,d2;
266:   PetscInt       n;
267:   PetscScalar    *pd;
268:   MPI_Comm       comm;
269:   CTX_BRUSSEL    *ctx;

271:   PetscFunctionBeginUser;
272:   PetscCall(MatShellGetContext(A,&ctx));
273:   PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
274:   PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
275:   PetscCall(VecGetArray(diag,&pd));
276:   PetscCall(VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1));
277:   PetscCall(VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2));

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

282:   PetscCall(VecDestroy(&d1));
283:   PetscCall(VecDestroy(&d2));
284:   PetscCall(VecRestoreArray(diag,&pd));
285:   PetscFunctionReturn(PETSC_SUCCESS);
286: }

288: /*TEST

290:    test:
291:       suffix: 1
292:       args: -n 50 -eps_nev 4 -eps_two_sided {{0 1}} -eps_type {{krylovschur lapack}} -terse
293:       requires: !single
294:       filter: grep -v method

296:    test:
297:       suffix: 2
298:       args: -eps_nev 8 -eps_max_it 300 -eps_target -28 -rg_type interval -rg_interval_endpoints -40,-20,-.1,.1 -terse
299:       requires: !single

301:    test:
302:       suffix: 3
303:       args: -n 50 -eps_nev 4 -eps_balance twoside -terse
304:       requires: double
305:       filter: grep -v method
306:       output_file: output/ex9_1.out

308:    test:
309:       suffix: 4
310:       args: -eps_smallest_imaginary -eps_ncv 24 -terse
311:       requires: !complex !single

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

318:    test:
319:       suffix: 5
320:       args: -eps_nev 4 -eps_target_real -eps_target -3 -terse
321:       requires: !single

323:    test:
324:       suffix: 6
325:       args: -eps_nev 2 -eps_target_imaginary -eps_target 3i -terse
326:       requires: complex !single

328:    test:
329:       suffix: 7
330:       args: -n 40 -eps_nev 1 -eps_type arnoldi -eps_smallest_real -eps_refined -eps_ncv 42 -eps_max_it 300 -terse
331:       requires: double

333:    test:
334:       suffix: 8
335:       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
336:       requires: !single
337:       filter: sed -e "s/[+-]0\.0*i//g"
338:       timeoutfactor: 2

340:    test:
341:       suffix: 9
342:       args: -eps_largest_imaginary -eps_ncv 24 -terse
343:       requires: !single

345: TEST*/