Actual source code: ex30.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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[] = "Illustrates the use of a region for filtering; the number of wanted eigenvalues in not known a priori.\n\n"
 12:   "The problem is the Brusselator wave model as in ex9.c.\n"
 13:   "The command line options are:\n"
 14:   "  -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
 15:   "  -L <L>, where <L> = bifurcation parameter.\n"
 16:   "  -alpha <alpha>, -beta <beta>, -delta1 <delta1>,  -delta2 <delta2>,\n"
 17:   "       where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";

 19: #include <slepceps.h>

 21: /*
 22:    This example tries to compute all eigenvalues lying outside the real axis.
 23:    This could be achieved by computing LARGEST_IMAGINARY eigenvalues, but
 24:    here we take a different route: define a region of the complex plane where
 25:    eigenvalues must be emphasized (eigenvalues outside the region are filtered
 26:    out). In this case, we select the region as the complement of a thin stripe
 27:    around the real axis.
 28:  */

 30: PetscErrorCode MatMult_Brussel(Mat,Vec,Vec);
 31: PetscErrorCode MatShift_Brussel(PetscScalar*,Mat);
 32: PetscErrorCode MatGetDiagonal_Brussel(Mat,Vec);
 33: PetscErrorCode MyStoppingTest(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*);

 35: typedef struct {
 36:   Mat         T;
 37:   Vec         x1,x2,y1,y2;
 38:   PetscScalar alpha,beta,tau1,tau2,sigma;
 39:   PetscInt    lastnconv;      /* last value of nconv; used in stopping test */
 40:   PetscInt    nreps;          /* number of repetitions of nconv; used in stopping test */
 41: } CTX_BRUSSEL;

 43: int main(int argc,char **argv)
 44: {
 45:   Mat            A;               /* eigenvalue problem matrix */
 46:   EPS            eps;             /* eigenproblem solver context */
 47:   RG             rg;              /* region object */
 48:   PetscScalar    delta1,delta2,L,h;
 49:   PetscInt       N=30,n,i,Istart,Iend,mpd;
 50:   CTX_BRUSSEL    *ctx;
 51:   PetscBool      terse;
 52:   PetscViewer    viewer;

 55:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

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

 60:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 61:         Generate the matrix
 62:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 64:   /*
 65:      Create shell matrix context and set default parameters
 66:   */
 67:   PetscNew(&ctx);
 68:   ctx->alpha = 2.0;
 69:   ctx->beta  = 5.45;
 70:   delta1     = 0.008;
 71:   delta2     = 0.004;
 72:   L          = 0.51302;

 74:   /*
 75:      Look the command line for user-provided parameters
 76:   */
 77:   PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL);
 78:   PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL);
 79:   PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL);
 80:   PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL);
 81:   PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL);

 83:   /*
 84:      Create matrix T
 85:   */
 86:   MatCreate(PETSC_COMM_WORLD,&ctx->T);
 87:   MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N);
 88:   MatSetFromOptions(ctx->T);
 89:   MatSetUp(ctx->T);

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

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

114:   /*
115:      Create the shell matrix
116:   */
117:   MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);
118:   MatShellSetOperation(A,MATOP_MULT,(void(*)())MatMult_Brussel);
119:   MatShellSetOperation(A,MATOP_SHIFT,(void(*)())MatShift_Brussel);
120:   MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Brussel);

122:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123:                 Create the eigensolver and configure the region
124:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

126:   EPSCreate(PETSC_COMM_WORLD,&eps);
127:   EPSSetOperators(eps,A,NULL);
128:   EPSSetProblemType(eps,EPS_NHEP);

130:   /*
131:      Define the region containing the eigenvalues of interest
132:   */
133:   EPSGetRG(eps,&rg);
134:   RGSetType(rg,RGINTERVAL);
135:   RGIntervalSetEndpoints(rg,-PETSC_INFINITY,PETSC_INFINITY,-0.01,0.01);
136:   RGSetComplement(rg,PETSC_TRUE);
137:   /* sort eigenvalue approximations wrt a target, otherwise convergence will be erratic */
138:   EPSSetTarget(eps,0.0);
139:   EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);

141:   /*
142:      Set solver options. In particular, we must allocate sufficient
143:      storage for all eigenpairs that may converge (ncv). This is
144:      application-dependent.
145:   */
146:   mpd = 40;
147:   EPSSetDimensions(eps,2*mpd,3*mpd,mpd);
148:   EPSSetTolerances(eps,1e-7,2000);
149:   ctx->lastnconv = 0;
150:   ctx->nreps     = 0;
151:   EPSSetStoppingTestFunction(eps,MyStoppingTest,(void*)ctx,NULL);
152:   EPSSetFromOptions(eps);

154:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155:                 Solve the eigensystem and display solution
156:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

158:   EPSSolve(eps);

160:   /* show detailed info unless -terse option is given by user */
161:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
162:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
163:   EPSReasonView(eps,viewer);
164:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
165:   if (!terse) {
166:     EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer);
167:   }
168:   PetscViewerPopFormat(viewer);

170:   EPSDestroy(&eps);
171:   MatDestroy(&A);
172:   MatDestroy(&ctx->T);
173:   VecDestroy(&ctx->x1);
174:   VecDestroy(&ctx->x2);
175:   VecDestroy(&ctx->y1);
176:   VecDestroy(&ctx->y2);
177:   PetscFree(ctx);
178:   SlepcFinalize();
179:   return ierr;
180: }

182: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
183: {
184:   PetscInt          n;
185:   const PetscScalar *px;
186:   PetscScalar       *py;
187:   CTX_BRUSSEL       *ctx;
188:   PetscErrorCode    ierr;

191:   MatShellGetContext(A,(void**)&ctx);
192:   MatGetLocalSize(ctx->T,&n,NULL);
193:   VecGetArrayRead(x,&px);
194:   VecGetArray(y,&py);
195:   VecPlaceArray(ctx->x1,px);
196:   VecPlaceArray(ctx->x2,px+n);
197:   VecPlaceArray(ctx->y1,py);
198:   VecPlaceArray(ctx->y2,py+n);

200:   MatMult(ctx->T,ctx->x1,ctx->y1);
201:   VecScale(ctx->y1,ctx->tau1);
202:   VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1);
203:   VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2);

205:   MatMult(ctx->T,ctx->x2,ctx->y2);
206:   VecScale(ctx->y2,ctx->tau2);
207:   VecAXPY(ctx->y2,-ctx->beta,ctx->x1);
208:   VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2);

210:   VecRestoreArrayRead(x,&px);
211:   VecRestoreArray(y,&py);
212:   VecResetArray(ctx->x1);
213:   VecResetArray(ctx->x2);
214:   VecResetArray(ctx->y1);
215:   VecResetArray(ctx->y2);
216:   return(0);
217: }

219: PetscErrorCode MatShift_Brussel(PetscScalar* a,Mat Y)
220: {
221:   CTX_BRUSSEL    *ctx;

225:   MatShellGetContext(Y,(void**)&ctx);
226:   ctx->sigma += *a;
227:   return(0);
228: }

230: PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
231: {
232:   Vec            d1,d2;
233:   PetscInt       n;
234:   PetscScalar    *pd;
235:   MPI_Comm       comm;
236:   CTX_BRUSSEL    *ctx;

240:   MatShellGetContext(A,(void**)&ctx);
241:   PetscObjectGetComm((PetscObject)A,&comm);
242:   MatGetLocalSize(ctx->T,&n,NULL);
243:   VecGetArray(diag,&pd);
244:   VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1);
245:   VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2);

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

250:   VecDestroy(&d1);
251:   VecDestroy(&d2);
252:   VecRestoreArray(diag,&pd);
253:   return(0);
254: }

256: /*
257:     Function for user-defined stopping test.

259:     Ignores the value of nev. It only takes into account the number of
260:     eigenpairs that have converged in recent outer iterations (restarts);
261:     if no new eigenvalus have converged in the last few restarts,
262:     we stop the iteration, assuming that no more eigenvalues are present
263:     inside the region.
264: */
265: PetscErrorCode MyStoppingTest(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ptr)
266: {
268:   CTX_BRUSSEL    *ctx = (CTX_BRUSSEL*)ptr;

271:   /* check usual termination conditions, but ignoring the case nconv>=nev */
272:   EPSStoppingBasic(eps,its,max_it,nconv,PETSC_MAX_INT,reason,NULL);
273:   if (*reason==EPS_CONVERGED_ITERATING) {
274:     /* check if nconv is the same as before */
275:     if (nconv==ctx->lastnconv) ctx->nreps++;
276:     else {
277:       ctx->lastnconv = nconv;
278:       ctx->nreps     = 0;
279:     }
280:     /* check if no eigenvalues converged in last 10 restarts */
281:     if (nconv && ctx->nreps>10) *reason = EPS_CONVERGED_USER;
282:   }
283:   return(0);
284: }