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);
142:   EPSSetProblemType(eps,EPS_NHEP);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

274:   return 0;
275: }

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

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

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

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

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

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

315:   return 0;
316: }