Actual source code: ex9.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:       SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:       Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain

  6:       This file is part of SLEPc. See the README file for conditions of use
  7:       and additional information.
  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:  */


 36: /* 
 37:    Matrix operations
 38: */
 39: PetscErrorCode MatBrussel_Mult(Mat,Vec,Vec);
 40: PetscErrorCode MatBrussel_Shift(PetscScalar*,Mat);
 41: PetscErrorCode MatBrussel_GetDiagonal(Mat,Vec);

 43: typedef struct {
 44:   Mat         T;
 45:   Vec         x1, x2, y1, y2;
 46:   PetscScalar alpha, beta, tau1, tau2, sigma;
 47: } CTX_BRUSSEL;

 51: int main( int argc, char **argv )
 52: {
 53:   Mat                  A;                  /* eigenvalue problem matrix */
 54:   EPS                  eps;                  /* eigenproblem solver context */
 55:   EPSType              type;
 56:   PetscReal            error, tol, re, im;
 57:   PetscScalar          delta1, delta2, L, h, kr, ki, value[3];
 58:   PetscInt             N=30, n, i, col[3], Istart, Iend;
 59:   int                  nev, maxit, its, nconv;
 60:   PetscTruth     FirstBlock=PETSC_FALSE, LastBlock=PETSC_FALSE;
 62:   CTX_BRUSSEL    *ctx;

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

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

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 70:         Generate the matrix 
 71:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 73:   /* 
 74:      Create shell matrix context and set default parameters
 75:   */
 76:   PetscNew(CTX_BRUSSEL,&ctx);
 77:   ctx->alpha = 2.0;
 78:   ctx->beta  = 5.45;
 79:   delta1     = 0.008;
 80:   delta2     = 0.004;
 81:   L          = 0.51302;

 83:   /* 
 84:      Look the command line for user-provided parameters
 85:   */
 86:   PetscOptionsGetScalar(PETSC_NULL,"-L",&L,PETSC_NULL);
 87:   PetscOptionsGetScalar(PETSC_NULL,"-alpha",&ctx->alpha,PETSC_NULL);
 88:   PetscOptionsGetScalar(PETSC_NULL,"-beta",&ctx->beta,PETSC_NULL);
 89:   PetscOptionsGetScalar(PETSC_NULL,"-delta1",&delta1,PETSC_NULL);
 90:   PetscOptionsGetScalar(PETSC_NULL,"-delta2",&delta2,PETSC_NULL);

 92:   /* 
 93:      Create matrix T
 94:   */
 95:   MatCreate(PETSC_COMM_WORLD,&ctx->T);
 96:   MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N);
 97:   MatSetFromOptions(ctx->T);
 98: 
 99:   MatGetOwnershipRange(ctx->T,&Istart,&Iend);
100:   if (Istart==0) FirstBlock=PETSC_TRUE;
101:   if (Iend==N) LastBlock=PETSC_TRUE;
102:   value[0]=1.0; value[1]=-2.0; value[2]=1.0;
103:   for( i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++ ) {
104:     col[0]=i-1; col[1]=i; col[2]=i+1;
105:     MatSetValues(ctx->T,1,&i,3,col,value,INSERT_VALUES);
106:   }
107:   if (LastBlock) {
108:     i=N-1; col[0]=N-2; col[1]=N-1;
109:     MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);
110:   }
111:   if (FirstBlock) {
112:     i=0; col[0]=0; col[1]=1; value[0]=-2.0; value[1]=1.0;
113:     MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);
114:   }

116:   MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY);
117:   MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY);
118:   MatGetLocalSize(ctx->T,&n,PETSC_NULL);

120:   /* 
121:      Fill the remaining information in the shell matrix context
122:      and create auxiliary vectors
123:   */
124:   h = 1.0 / (double)(N+1);
125:   ctx->tau1 = delta1 / ((h*L)*(h*L));
126:   ctx->tau2 = delta2 / ((h*L)*(h*L));
127:   ctx->sigma = 0.0;
128:   VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->x1);
129:   VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->x2);
130:   VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->y1);
131:   VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->y2);

133:   /* 
134:      Create the shell matrix
135:   */
136:   MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);
137:   MatShellSetOperation(A,MATOP_MULT,(void(*)())MatBrussel_Mult);
138:   MatShellSetOperation(A,MATOP_SHIFT,(void(*)())MatBrussel_Shift);
139:   MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatBrussel_GetDiagonal);

141:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
142:                 Create the eigensolver and set various options
143:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

145:   /* 
146:      Create eigensolver context
147:   */
148:   EPSCreate(PETSC_COMM_WORLD,&eps);

150:   /* 
151:      Set operators. In this case, it is a standard eigenvalue problem
152:   */
153:   EPSSetOperators(eps,A,PETSC_NULL);
154:   EPSSetProblemType(eps,EPS_NHEP);

156:   /*
157:      Ask for the rightmost eigenvalues
158:   */
159:   EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);

161:   /*
162:      Set other solver options at runtime
163:   */
164:   EPSSetFromOptions(eps);

166:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
167:                       Solve the eigensystem
168:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

184:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
185:                     Display solution and clean up
186:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

188:   /* 
189:      Get number of converged eigenpairs
190:   */
191:   EPSGetConverged(eps,&nconv);
192:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);

194:   if (nconv>0) {
195:     /*
196:        Display eigenvalues and relative errors
197:     */
198:     PetscPrintf(PETSC_COMM_WORLD,
199:          "           k             ||Ax-kx||/||kx||\n"
200:          "  --------------------- ------------------\n" );
201:     for( i=0; i<nconv; i++ ) {
202:       /* 
203:          Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
204:          ki (imaginary part)
205:       */
206:       EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);

208:       /*
209:          Compute the relative error associated to each eigenpair
210:       */
211:       EPSComputeRelativeError(eps,i,&error);

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

247: PetscErrorCode MatBrussel_Mult(Mat A,Vec x,Vec y)
248: {
250:   PetscInt       n;
251:   PetscScalar    *px, *py;
252:   CTX_BRUSSEL    *ctx;

255:   MatShellGetContext(A,(void**)&ctx);
256:   MatGetLocalSize(ctx->T,&n,PETSC_NULL);
257:   VecGetArray(x,&px);
258:   VecGetArray(y,&py);
259:   VecPlaceArray(ctx->x1,px);
260:   VecPlaceArray(ctx->x2,px+n);
261:   VecPlaceArray(ctx->y1,py);
262:   VecPlaceArray(ctx->y2,py+n);

264:   MatMult(ctx->T,ctx->x1,ctx->y1);
265:   VecScale(ctx->y1,ctx->tau1);
266:   VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1);
267:   VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2);

269:   MatMult(ctx->T,ctx->x2,ctx->y2);
270:   VecScale(ctx->y2,ctx->tau2);
271:   VecAXPY(ctx->y2,-ctx->beta,ctx->x1);
272:   VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2);

274:   VecRestoreArray(x,&px);
275:   VecRestoreArray(y,&py);
276:   VecResetArray(ctx->x1);
277:   VecResetArray(ctx->x2);
278:   VecResetArray(ctx->y1);
279:   VecResetArray(ctx->y2);
280:   return(0);
281: }

285: PetscErrorCode MatBrussel_Shift( PetscScalar* a, Mat Y )
286: {
287:   CTX_BRUSSEL    *ctx;

291:   MatShellGetContext( Y, (void**)&ctx );
292:   ctx->sigma += *a;
293:   return(0);
294: }

298: int MatBrussel_GetDiagonal(Mat A,Vec diag)
299: {
300:   Vec            d1, d2;
302:   PetscInt       n;
303:   PetscScalar    *pd;
304:   MPI_Comm       comm;
305:   CTX_BRUSSEL    *ctx;

308:   MatShellGetContext(A,(void**)&ctx);
309:   PetscObjectGetComm((PetscObject)A,&comm);
310:   MatGetLocalSize(ctx->T,&n,PETSC_NULL);
311:   VecGetArray(diag,&pd);
312:   VecCreateMPIWithArray(comm,n,PETSC_DECIDE,pd,&d1);
313:   VecCreateMPIWithArray(comm,n,PETSC_DECIDE,pd+n,&d2);

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

318:   VecDestroy(d1);
319:   VecDestroy(d2);
320:   VecRestoreArray(diag,&pd);
321:   return(0);
322: }