Actual source code: ex9.c

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

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

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

 29: #include <slepceps.h>

 31: /*
 32:    This example computes the eigenvalues with largest real part of the
 33:    following matrix

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

 38:    where

 40:         T = tridiag{1,-2,1}
 41:         h = 1/(n+1)
 42:         tau1 = delta1/(h*L)^2
 43:         tau2 = delta2/(h*L)^2
 44:  */


 47: /*
 48:    Matrix operations
 49: */
 50: PetscErrorCode MatMult_Brussel(Mat,Vec,Vec);
 51: PetscErrorCode MatShift_Brussel(PetscScalar*,Mat);
 52: PetscErrorCode MatGetDiagonal_Brussel(Mat,Vec);

 54: typedef struct {
 55:   Mat         T;
 56:   Vec         x1,x2,y1,y2;
 57:   PetscScalar alpha,beta,tau1,tau2,sigma;
 58: } CTX_BRUSSEL;

 62: int main(int argc,char **argv)
 63: {
 64:   Mat            A;               /* eigenvalue problem matrix */
 65:   EPS            eps;             /* eigenproblem solver context */
 66:   EPSType        type;
 67:   PetscScalar    delta1,delta2,L,h,value[3];
 68:   PetscInt       N=30,n,i,col[3],Istart,Iend,nev;
 69:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
 70:   CTX_BRUSSEL    *ctx;

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

 75:   PetscOptionsGetInt(NULL,"-n",&N,NULL);
 76:   PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%D\n\n",N);

 78:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 79:         Generate the matrix
 80:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 82:   /*
 83:      Create shell matrix context and set default parameters
 84:   */
 85:   PetscNew(CTX_BRUSSEL,&ctx);
 86:   ctx->alpha = 2.0;
 87:   ctx->beta  = 5.45;
 88:   delta1     = 0.008;
 89:   delta2     = 0.004;
 90:   L          = 0.51302;

 92:   /*
 93:      Look the command line for user-provided parameters
 94:   */
 95:   PetscOptionsGetScalar(NULL,"-L",&L,NULL);
 96:   PetscOptionsGetScalar(NULL,"-alpha",&ctx->alpha,NULL);
 97:   PetscOptionsGetScalar(NULL,"-beta",&ctx->beta,NULL);
 98:   PetscOptionsGetScalar(NULL,"-delta1",&delta1,NULL);
 99:   PetscOptionsGetScalar(NULL,"-delta2",&delta2,NULL);

101:   /*
102:      Create matrix T
103:   */
104:   MatCreate(PETSC_COMM_WORLD,&ctx->T);
105:   MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N);
106:   MatSetFromOptions(ctx->T);
107:   MatSetUp(ctx->T);

109:   MatGetOwnershipRange(ctx->T,&Istart,&Iend);
110:   if (Istart==0) FirstBlock=PETSC_TRUE;
111:   if (Iend==N) LastBlock=PETSC_TRUE;
112:   value[0]=1.0; value[1]=-2.0; value[2]=1.0;
113:   for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
114:     col[0]=i-1; col[1]=i; col[2]=i+1;
115:     MatSetValues(ctx->T,1,&i,3,col,value,INSERT_VALUES);
116:   }
117:   if (LastBlock) {
118:     i=N-1; col[0]=N-2; col[1]=N-1;
119:     MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);
120:   }
121:   if (FirstBlock) {
122:     i=0; col[0]=0; col[1]=1; value[0]=-2.0; value[1]=1.0;
123:     MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);
124:   }

126:   MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY);
127:   MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY);
128:   MatGetLocalSize(ctx->T,&n,NULL);

130:   /*
131:      Fill the remaining information in the shell matrix context
132:      and create auxiliary vectors
133:   */
134:   h = 1.0 / (PetscReal)(N+1);
135:   ctx->tau1 = delta1 / ((h*L)*(h*L));
136:   ctx->tau2 = delta2 / ((h*L)*(h*L));
137:   ctx->sigma = 0.0;
138:   VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1);
139:   VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2);
140:   VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1);
141:   VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2);

143:   /*
144:      Create the shell matrix
145:   */
146:   MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);
147:   MatShellSetOperation(A,MATOP_MULT,(void(*)())MatMult_Brussel);
148:   MatShellSetOperation(A,MATOP_SHIFT,(void(*)())MatShift_Brussel);
149:   MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Brussel);

151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152:                 Create the eigensolver and set various options
153:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

155:   /*
156:      Create eigensolver context
157:   */
158:   EPSCreate(PETSC_COMM_WORLD,&eps);

160:   /*
161:      Set operators. In this case, it is a standard eigenvalue problem
162:   */
163:   EPSSetOperators(eps,A,NULL);
164:   EPSSetProblemType(eps,EPS_NHEP);

166:   /*
167:      Ask for the rightmost eigenvalues
168:   */
169:   EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);

171:   /*
172:      Set other solver options at runtime
173:   */
174:   EPSSetFromOptions(eps);

176:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177:                       Solve the eigensystem
178:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

180:   EPSSolve(eps);

182:   /*
183:      Optional: Get some information from the solver and display it
184:   */
185:   EPSGetType(eps,&type);
186:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
187:   EPSGetDimensions(eps,&nev,NULL,NULL);
188:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

190:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191:                     Display solution and clean up
192:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

194:   EPSPrintSolution(eps,NULL);
195:   EPSDestroy(&eps);
196:   MatDestroy(&A);
197:   MatDestroy(&ctx->T);
198:   VecDestroy(&ctx->x1);
199:   VecDestroy(&ctx->x2);
200:   VecDestroy(&ctx->y1);
201:   VecDestroy(&ctx->y2);
202:   PetscFree(ctx);
203:   SlepcFinalize();
204:   return 0;
205: }

209: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
210: {
211:   PetscInt          n;
212:   const PetscScalar *px;
213:   PetscScalar       *py;
214:   CTX_BRUSSEL       *ctx;
215:   PetscErrorCode    ierr;

218:   MatShellGetContext(A,(void**)&ctx);
219:   MatGetLocalSize(ctx->T,&n,NULL);
220:   VecGetArrayRead(x,&px);
221:   VecGetArray(y,&py);
222:   VecPlaceArray(ctx->x1,px);
223:   VecPlaceArray(ctx->x2,px+n);
224:   VecPlaceArray(ctx->y1,py);
225:   VecPlaceArray(ctx->y2,py+n);

227:   MatMult(ctx->T,ctx->x1,ctx->y1);
228:   VecScale(ctx->y1,ctx->tau1);
229:   VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1);
230:   VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2);

232:   MatMult(ctx->T,ctx->x2,ctx->y2);
233:   VecScale(ctx->y2,ctx->tau2);
234:   VecAXPY(ctx->y2,-ctx->beta,ctx->x1);
235:   VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2);

237:   VecRestoreArrayRead(x,&px);
238:   VecRestoreArray(y,&py);
239:   VecResetArray(ctx->x1);
240:   VecResetArray(ctx->x2);
241:   VecResetArray(ctx->y1);
242:   VecResetArray(ctx->y2);
243:   return(0);
244: }

248: PetscErrorCode MatShift_Brussel(PetscScalar* a,Mat Y)
249: {
250:   CTX_BRUSSEL    *ctx;

254:   MatShellGetContext(Y,(void**)&ctx);
255:   ctx->sigma += *a;
256:   return(0);
257: }

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

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

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

281:   VecDestroy(&d1);
282:   VecDestroy(&d2);
283:   VecRestoreArray(diag,&pd);
284:   return(0);
285: }