Actual source code: ex9.c

  2: static char help[] = "Solves a problem associated to the Brusselator wave model in chemical reactions, illustrating the use of shell matrices.\n\n"
  3:   "The command line options are:\n\n"
  4:   "  -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
  5:   "  -L <L>, where <L> = bifurcation parameter.\n"
  6:   "  -alpha <alpha>, -beta <beta>, -delta1 <delta1>,  -delta2 <delta2>,\n"
  7:   "       where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";

 9:  #include slepceps.h

 11: /*
 12:    This example computes the eigenvalues with largest real part of the 
 13:    following matrix

 15:         A = [ tau1*T+(beta-1)*I     alpha^2*I
 16:                   -beta*I        tau2*T-alpha^2*I ],

 18:    where

 20:         T = tridiag{1,-2,1}
 21:         h = 1/(n+1)
 22:         tau1 = delta1/(h*L)^2
 23:         tau2 = delta2/(h*L)^2
 24:  */


 27: /* 
 28:    Matrix operations
 29: */
 30: extern int MatBrussel_Mult(Mat,Vec,Vec);
 31: extern int MatBrussel_Shift(PetscScalar*,Mat);
 32: extern int MatBrussel_GetDiagonal(Mat,Vec);

 34: typedef struct {
 35:   Mat         T;
 36:   Vec         x1, x2, y1, y2;
 37:   PetscScalar alpha, beta, tau1, tau2, sigma;
 38: } CTX_BRUSSEL;

 42: int main( int argc, char **argv )
 43: {
 44:   Mat         A;               /* eigenvalue problem matrix */
 45:   EPS         eps;             /* eigenproblem solver context */
 46:   EPSType     type;
 47:   PetscReal   error, tol, re, im;
 48:   PetscScalar delta1, delta2, L, h, kr, ki, value[3];
 49:   int         N=30, n, nev, ierr, maxit, i, its, nconv,
 50:               col[3], Istart, Iend, FirstBlock=0, LastBlock=0;
 51:   CTX_BRUSSEL *ctx;

 53:   SlepcInitialize(&argc,&argv,(char*)0,help);

 55:   PetscOptionsGetInt(PETSC_NULL,"-n",&N,PETSC_NULL);
 56:   PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%d\n\n",N);

 58:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 59:         Generate the matrix 
 60:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 62:   /* 
 63:      Create shell matrix context and set default parameters
 64:   */
 65:   PetscNew(CTX_BRUSSEL,&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:   PetscOptionsGetScalar(PETSC_NULL,"-L",&L,PETSC_NULL);
 76:   PetscOptionsGetScalar(PETSC_NULL,"-alpha",&ctx->alpha,PETSC_NULL);
 77:   PetscOptionsGetScalar(PETSC_NULL,"-beta",&ctx->beta,PETSC_NULL);
 78:   PetscOptionsGetScalar(PETSC_NULL,"-delta1",&delta1,PETSC_NULL);
 79:   PetscOptionsGetScalar(PETSC_NULL,"-delta2",&delta2,PETSC_NULL);

 81:   /* 
 82:      Create matrix T
 83:   */
 84:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&ctx->T);
 85:   MatSetFromOptions(ctx->T);
 86: 
 87:   MatGetOwnershipRange(ctx->T,&Istart,&Iend);
 88:   if (Istart==0) FirstBlock=PETSC_TRUE;
 89:   if (Iend==N) LastBlock=PETSC_TRUE;
 90:   value[0]=1.0; value[1]=-2.0; value[2]=1.0;
 91:   for( i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++ ) {
 92:     col[0]=i-1; col[1]=i; col[2]=i+1;
 93:     MatSetValues(ctx->T,1,&i,3,col,value,INSERT_VALUES);
 94:   }
 95:   if (LastBlock) {
 96:     i=N-1; col[0]=N-2; col[1]=N-1;
 97:     MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);
 98:   }
 99:   if (FirstBlock) {
100:     i=0; col[0]=0; col[1]=1; value[0]=-2.0; value[1]=1.0;
101:     MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);
102:   }

104:   MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY);
105:   MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY);
106:   MatGetLocalSize(ctx->T,&n,PETSC_NULL);

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

121:   /* 
122:      Create the shell matrix
123:   */
124:   MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);
125:   MatShellSetOperation(A,MATOP_MULT,(void(*)())MatBrussel_Mult);
126:   MatShellSetOperation(A,MATOP_SHIFT,(void(*)())MatBrussel_Shift);
127:   MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatBrussel_GetDiagonal);

129:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
130:                 Create the eigensolver and set various options
131:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

133:   /* 
134:      Create eigensolver context
135:   */
136:   EPSCreate(PETSC_COMM_WORLD,&eps);

138:   /* 
139:      Set operators. In this case, it is a standard eigenvalue problem
140:   */
141:   EPSSetOperators(eps,A,PETSC_NULL);

143:   /*
144:      Force to use ARPACK if it is installed and ask for the rightmost eigenvalues
145:   */
146: #if defined(SLEPC_HAVE_ARPACK)
147:   EPSSetType(eps,EPSARPACK);
148: #else
149:   EPSSetType(eps, EPSLAPACK);
150: #endif
151:   EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);

153:   /*
154:      Set other solver options at runtime
155:   */
156:   EPSSetFromOptions(eps);

158:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
159:                       Solve the eigensystem
160:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

162:   EPSSolve(eps);
163:   EPSGetIterationNumber(eps, &its);
164:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
165: 
166:   /*
167:      Optional: Get some information from the solver and display it
168:   */
169:   EPSGetType(eps,&type);
170:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
171:   EPSGetDimensions(eps,&nev,PETSC_NULL);
172:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
173:   EPSGetTolerances(eps,&tol,&maxit);
174:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);

176:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
177:                     Display solution and clean up
178:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

180:   /* 
181:      Get number of converged eigenpairs
182:   */
183:   EPSGetConverged(eps,&nconv);
184:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);

186:   if (nconv>0) {
187:     /*
188:        Display eigenvalues and relative errors
189:     */
190:     PetscPrintf(PETSC_COMM_WORLD,
191:          "           k             ||Ax-kx||/||kx||\n"
192:          "  --------------------- ------------------\n" );
193:     for( i=0; i<nconv; i++ ) {
194:       /* 
195:          Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
196:          ki (imaginary part)
197:       */
198:       EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);

200:       /*
201:          Compute the relative error associated to each eigenpair
202:       */
203:       EPSComputeRelativeError(eps,i,&error);

205: #if defined(PETSC_USE_COMPLEX)
206:       re = PetscRealPart(kr);
207:       im = PetscImaginaryPart(kr);
208: #else
209:       re = kr;
210:       im = ki;
211: #endif
212:       if( im != 0.0 ) {
213:         PetscPrintf(PETSC_COMM_WORLD," % 6f %+6f i",re,im);
214:       } else {
215:         PetscPrintf(PETSC_COMM_WORLD,"       % 6f      ",re);
216:       }
217:       PetscPrintf(PETSC_COMM_WORLD," % 12f\n",error);
218:     }
219:     PetscPrintf(PETSC_COMM_WORLD,"\n" );
220:   }
221: 
222:   /* 
223:      Free work space
224:   */
225:   EPSDestroy(eps);
226:   MatDestroy(A);
227:   MatDestroy(ctx->T);
228:   VecDestroy(ctx->x1);
229:   VecDestroy(ctx->x2);
230:   VecDestroy(ctx->y1);
231:   VecDestroy(ctx->y2);
232:   PetscFree(ctx);
233:   SlepcFinalize();
234:   return 0;
235: }

239: int MatBrussel_Mult(Mat A,Vec x,Vec y)
240: {
241:   int         n, ierr;
242:   PetscScalar alpha, *px, *py;
243:   MPI_Comm    comm;
244:   CTX_BRUSSEL *ctx;

246:   MatShellGetContext(A,(void**)&ctx);
247:   PetscObjectGetComm((PetscObject)A,&comm);
248:   MatGetLocalSize(ctx->T,&n,PETSC_NULL);
249:   VecGetArray(x,&px);
250:   VecGetArray(y,&py);
251:   VecPlaceArray(ctx->x1,px);
252:   VecPlaceArray(ctx->x2,px+n);
253:   VecPlaceArray(ctx->y1,py);
254:   VecPlaceArray(ctx->y2,py+n);

256:   MatMult(ctx->T,ctx->x1,ctx->y1);
257:   VecScale(&ctx->tau1,ctx->y1);
258:   alpha = ctx->beta - 1.0 - ctx->sigma;
259:   VecAXPY(&alpha,ctx->x1,ctx->y1);
260:   alpha = ctx->alpha * ctx->alpha;
261:   VecAXPY(&alpha,ctx->x2,ctx->y1);

263:   MatMult(ctx->T,ctx->x2,ctx->y2);
264:   VecScale(&ctx->tau2,ctx->y2);
265:   alpha = -ctx->beta;
266:   VecAXPY(&alpha,ctx->x1,ctx->y2);
267:   alpha = -ctx->alpha * ctx->alpha - ctx->sigma;
268:   VecAXPY(&alpha,ctx->x2,ctx->y2);

270:   VecRestoreArray(x,&px);
271:   VecRestoreArray(y,&py);

273:   return 0;
274: }

278: int MatBrussel_Shift( PetscScalar* a, Mat Y )
279: {
280:   CTX_BRUSSEL *ctx;
281:   int         ierr;

283:   MatShellGetContext( Y, (void**)&ctx );
284:   ctx->sigma += *a;
285:   return(0);
286: }

290: int MatBrussel_GetDiagonal(Mat A,Vec diag)
291: {
292:   Vec         d1, d2;
293:   int         n, ierr;
294:   PetscScalar alpha, *pd;
295:   MPI_Comm    comm;
296:   CTX_BRUSSEL *ctx;

298:   MatShellGetContext(A,(void**)&ctx);
299:   PetscObjectGetComm((PetscObject)A,&comm);
300:   MatGetLocalSize(ctx->T,&n,PETSC_NULL);
301:   VecGetArray(diag,&pd);
302:   VecCreateMPIWithArray(comm,n,PETSC_DECIDE,pd,&d1);
303:   VecCreateMPIWithArray(comm,n,PETSC_DECIDE,pd+n,&d2);

305:   alpha = -2.0*ctx->tau1 + ctx->beta - 1.0 - ctx->sigma;
306:   VecSet(&alpha,d1);
307:   alpha = -2.0*ctx->tau2 - ctx->alpha*ctx->alpha - ctx->sigma;
308:   VecSet(&alpha,d2);

310:   VecDestroy(d1);
311:   VecDestroy(d2);
312:   VecRestoreArray(diag,&pd);

314:   return 0;
315: }