Actual source code: ex30.c

slepc-3.20.2 2024-03-15
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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 is 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 MatGetDiagonal_Brussel(Mat,Vec);
 32: PetscErrorCode MyStoppingTest(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*);

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

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

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

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

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

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

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

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

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

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

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

120:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121:                 Create the eigensolver and configure the region
122:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

124:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
125:   PetscCall(EPSSetOperators(eps,A,NULL));
126:   PetscCall(EPSSetProblemType(eps,EPS_NHEP));

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

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

152:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153:                 Solve the eigensystem and display solution
154:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

156:   PetscCall(EPSSolve(eps));

158:   /* show detailed info unless -terse option is given by user */
159:   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
160:   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
161:   PetscCall(EPSConvergedReasonView(eps,viewer));
162:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
163:   if (!terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer));
164:   PetscCall(PetscViewerPopFormat(viewer));

166:   PetscCall(EPSDestroy(&eps));
167:   PetscCall(MatDestroy(&A));
168:   PetscCall(MatDestroy(&ctx->T));
169:   PetscCall(VecDestroy(&ctx->x1));
170:   PetscCall(VecDestroy(&ctx->x2));
171:   PetscCall(VecDestroy(&ctx->y1));
172:   PetscCall(VecDestroy(&ctx->y2));
173:   PetscCall(PetscFree(ctx));
174:   PetscCall(SlepcFinalize());
175:   return 0;
176: }

178: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
179: {
180:   PetscInt          n;
181:   const PetscScalar *px;
182:   PetscScalar       *py;
183:   CTX_BRUSSEL       *ctx;

185:   PetscFunctionBeginUser;
186:   PetscCall(MatShellGetContext(A,&ctx));
187:   PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
188:   PetscCall(VecGetArrayRead(x,&px));
189:   PetscCall(VecGetArray(y,&py));
190:   PetscCall(VecPlaceArray(ctx->x1,px));
191:   PetscCall(VecPlaceArray(ctx->x2,px+n));
192:   PetscCall(VecPlaceArray(ctx->y1,py));
193:   PetscCall(VecPlaceArray(ctx->y2,py+n));

195:   PetscCall(MatMult(ctx->T,ctx->x1,ctx->y1));
196:   PetscCall(VecScale(ctx->y1,ctx->tau1));
197:   PetscCall(VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1));
198:   PetscCall(VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2));

200:   PetscCall(MatMult(ctx->T,ctx->x2,ctx->y2));
201:   PetscCall(VecScale(ctx->y2,ctx->tau2));
202:   PetscCall(VecAXPY(ctx->y2,-ctx->beta,ctx->x1));
203:   PetscCall(VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2));

205:   PetscCall(VecRestoreArrayRead(x,&px));
206:   PetscCall(VecRestoreArray(y,&py));
207:   PetscCall(VecResetArray(ctx->x1));
208:   PetscCall(VecResetArray(ctx->x2));
209:   PetscCall(VecResetArray(ctx->y1));
210:   PetscCall(VecResetArray(ctx->y2));
211:   PetscFunctionReturn(PETSC_SUCCESS);
212: }

214: PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
215: {
216:   Vec            d1,d2;
217:   PetscInt       n;
218:   PetscScalar    *pd;
219:   MPI_Comm       comm;
220:   CTX_BRUSSEL    *ctx;

222:   PetscFunctionBeginUser;
223:   PetscCall(MatShellGetContext(A,&ctx));
224:   PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
225:   PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
226:   PetscCall(VecGetArray(diag,&pd));
227:   PetscCall(VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1));
228:   PetscCall(VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2));

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

233:   PetscCall(VecDestroy(&d1));
234:   PetscCall(VecDestroy(&d2));
235:   PetscCall(VecRestoreArray(diag,&pd));
236:   PetscFunctionReturn(PETSC_SUCCESS);
237: }

239: /*
240:     Function for user-defined stopping test.

242:     Ignores the value of nev. It only takes into account the number of
243:     eigenpairs that have converged in recent outer iterations (restarts);
244:     if no new eigenvalues have converged in the last few restarts,
245:     we stop the iteration, assuming that no more eigenvalues are present
246:     inside the region.
247: */
248: PetscErrorCode MyStoppingTest(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ptr)
249: {
250:   CTX_BRUSSEL    *ctx = (CTX_BRUSSEL*)ptr;

252:   PetscFunctionBeginUser;
253:   /* check usual termination conditions, but ignoring the case nconv>=nev */
254:   PetscCall(EPSStoppingBasic(eps,its,max_it,nconv,PETSC_MAX_INT,reason,NULL));
255:   if (*reason==EPS_CONVERGED_ITERATING) {
256:     /* check if nconv is the same as before */
257:     if (nconv==ctx->lastnconv) ctx->nreps++;
258:     else {
259:       ctx->lastnconv = nconv;
260:       ctx->nreps     = 0;
261:     }
262:     /* check if no eigenvalues converged in last 10 restarts */
263:     if (nconv && ctx->nreps>10) *reason = EPS_CONVERGED_USER;
264:   }
265:   PetscFunctionReturn(PETSC_SUCCESS);
266: }

268: /*TEST

270:    test:
271:       suffix: 1
272:       args: -n 100 -terse
273:       requires: !single

275: TEST*/