Actual source code: test15.c

slepc-3.14.1 2020-12-08
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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 user-defined stopping test.\n\n"
 12:   "This is based on ex22.\n"
 13:   "The command line options are:\n"
 14:   "  -n <n>, where <n> = number of grid subdivisions.\n"
 15:   "  -tau <tau>, where <tau> is the delay parameter.\n\n";

 17: /*
 18:    Solve parabolic partial differential equation with time delay tau

 20:             u_t = u_xx + a*u(t) + b*u(t-tau)
 21:             u(0,t) = u(pi,t) = 0

 23:    with a = 20 and b(x) = -4.1+x*(1-exp(x-pi)).

 25:    Discretization leads to a DDE of dimension n

 27:             -u' = A*u(t) + B*u(t-tau)

 29:    which results in the nonlinear eigenproblem

 31:             (-lambda*I + A + exp(-tau*lambda)*B)*u = 0
 32: */

 34: #include <slepcnep.h>

 36: /*
 37:    User-defined routines
 38: */
 39: PetscErrorCode MyStoppingTest(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*);

 41: typedef struct {
 42:   PetscInt    lastnconv;      /* last value of nconv; used in stopping test */
 43:   PetscInt    nreps;          /* number of repetitions of nconv; used in stopping test */
 44: } CTX_DELAY;

 46: int main(int argc,char **argv)
 47: {
 48:   NEP            nep;
 49:   Mat            Id,A,B;
 50:   FN             f1,f2,f3;
 51:   RG             rg;
 52:   CTX_DELAY      *ctx;
 53:   Mat            mats[3];
 54:   FN             funs[3];
 55:   PetscScalar    coeffs[2],b;
 56:   PetscInt       n=128,Istart,Iend,i,mpd;
 57:   PetscReal      tau=0.001,h,a=20,xi;
 58:   PetscBool      terse;
 59:   PetscViewer    viewer;

 62:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 63:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 64:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 65:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%D, tau=%g\n\n",n,(double)tau);
 66:   h = PETSC_PI/(PetscReal)(n+1);

 68:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 69:      Create nonlinear eigensolver context
 70:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 72:   NEPCreate(PETSC_COMM_WORLD,&nep);

 74:   /* Identity matrix */
 75:   MatCreate(PETSC_COMM_WORLD,&Id);
 76:   MatSetSizes(Id,PETSC_DECIDE,PETSC_DECIDE,n,n);
 77:   MatSetFromOptions(Id);
 78:   MatSetUp(Id);
 79:   MatGetOwnershipRange(Id,&Istart,&Iend);
 80:   for (i=Istart;i<Iend;i++) {
 81:     MatSetValue(Id,i,i,1.0,INSERT_VALUES);
 82:   }
 83:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
 84:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
 85:   MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);

 87:   /* A = 1/h^2*tridiag(1,-2,1) + a*I */
 88:   MatCreate(PETSC_COMM_WORLD,&A);
 89:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 90:   MatSetFromOptions(A);
 91:   MatSetUp(A);
 92:   MatGetOwnershipRange(A,&Istart,&Iend);
 93:   for (i=Istart;i<Iend;i++) {
 94:     if (i>0) { MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES); }
 95:     if (i<n-1) { MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES); }
 96:     MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
 97:   }
 98:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 99:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
100:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);

102:   /* B = diag(b(xi)) */
103:   MatCreate(PETSC_COMM_WORLD,&B);
104:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
105:   MatSetFromOptions(B);
106:   MatSetUp(B);
107:   MatGetOwnershipRange(B,&Istart,&Iend);
108:   for (i=Istart;i<Iend;i++) {
109:     xi = (i+1)*h;
110:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
111:     MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES);
112:   }
113:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
114:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
115:   MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);

117:   /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
118:   FNCreate(PETSC_COMM_WORLD,&f1);
119:   FNSetType(f1,FNRATIONAL);
120:   coeffs[0] = -1.0; coeffs[1] = 0.0;
121:   FNRationalSetNumerator(f1,2,coeffs);

123:   FNCreate(PETSC_COMM_WORLD,&f2);
124:   FNSetType(f2,FNRATIONAL);
125:   coeffs[0] = 1.0;
126:   FNRationalSetNumerator(f2,1,coeffs);

128:   FNCreate(PETSC_COMM_WORLD,&f3);
129:   FNSetType(f3,FNEXP);
130:   FNSetScale(f3,-tau,1.0);

132:   /* Set the split operator */
133:   mats[0] = A;  funs[0] = f2;
134:   mats[1] = Id; funs[1] = f1;
135:   mats[2] = B;  funs[2] = f3;
136:   NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);

138:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139:                 Customize nonlinear solver; set runtime options
140:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

142:   NEPSetType(nep,NEPNLEIGS);
143:   NEPGetRG(nep,&rg);
144:   RGSetType(rg,RGINTERVAL);
145: #if defined(PETSC_USE_COMPLEX)
146:   RGIntervalSetEndpoints(rg,5,20,-0.001,0.001);
147: #else
148:   RGIntervalSetEndpoints(rg,5,20,-0.0,0.0);
149: #endif
150:   NEPSetTarget(nep,15.0);
151:   NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE);

153:   /*
154:      Set solver options. In particular, we must allocate sufficient
155:      storage for all eigenpairs that may converge (ncv). This is
156:      application-dependent.
157:   */
158:   mpd = 40;
159:   NEPSetDimensions(nep,2*mpd,3*mpd,mpd);
160:   NEPSetTolerances(nep,PETSC_DEFAULT,2000);
161:   PetscNew(&ctx);
162:   ctx->lastnconv = 0;
163:   ctx->nreps     = 0;
164:   NEPSetStoppingTestFunction(nep,MyStoppingTest,(void*)ctx,NULL);

166:   NEPSetFromOptions(nep);

168:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169:                       Solve the eigensystem
170:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

172:   NEPSolve(nep);

174:   /* show detailed info unless -terse option is given by user */
175:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
176:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
177:   NEPConvergedReasonView(nep,viewer);
178:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
179:   if (!terse) {
180:     NEPErrorView(nep,NEP_ERROR_BACKWARD,viewer);
181:   }
182:   PetscViewerPopFormat(viewer);

184:   NEPDestroy(&nep);
185:   MatDestroy(&Id);
186:   MatDestroy(&A);
187:   MatDestroy(&B);
188:   FNDestroy(&f1);
189:   FNDestroy(&f2);
190:   FNDestroy(&f3);
191:   PetscFree(ctx);
192:   SlepcFinalize();
193:   return ierr;
194: }

196: /*
197:     Function for user-defined stopping test.

199:     Ignores the value of nev. It only takes into account the number of
200:     eigenpairs that have converged in recent outer iterations (restarts);
201:     if no new eigenvalues have converged in the last few restarts,
202:     we stop the iteration, assuming that no more eigenvalues are present
203:     inside the region.
204: */
205: PetscErrorCode MyStoppingTest(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ptr)
206: {
208:   CTX_DELAY      *ctx = (CTX_DELAY*)ptr;

211:   /* check usual termination conditions, but ignoring the case nconv>=nev */
212:   NEPStoppingBasic(nep,its,max_it,nconv,PETSC_MAX_INT,reason,NULL);
213:   if (*reason==NEP_CONVERGED_ITERATING) {
214:     /* check if nconv is the same as before */
215:     if (nconv==ctx->lastnconv) ctx->nreps++;
216:     else {
217:       ctx->lastnconv = nconv;
218:       ctx->nreps     = 0;
219:     }
220:     /* check if no eigenvalues converged in last 10 restarts */
221:     if (nconv && ctx->nreps>10) *reason = NEP_CONVERGED_USER;
222:   }
223:   return(0);
224: }

226: /*TEST

228:    test:
229:       suffix: 1
230:       args: -terse

232: TEST*/