Actual source code: test15.c

slepc-3.21.1 2024-04-26
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 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;

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

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

 72:   PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));

 74:   /* Identity matrix */
 75:   PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id));
 76:   PetscCall(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE));

 78:   /* A = 1/h^2*tridiag(1,-2,1) + a*I */
 79:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 80:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
 81:   PetscCall(MatSetFromOptions(A));
 82:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 83:   for (i=Istart;i<Iend;i++) {
 84:     if (i>0) PetscCall(MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES));
 85:     if (i<n-1) PetscCall(MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES));
 86:     PetscCall(MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES));
 87:   }
 88:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 89:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
 90:   PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));

 92:   /* B = diag(b(xi)) */
 93:   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
 94:   PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
 95:   PetscCall(MatSetFromOptions(B));
 96:   PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
 97:   for (i=Istart;i<Iend;i++) {
 98:     xi = (i+1)*h;
 99:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
100:     PetscCall(MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES));
101:   }
102:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
103:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
104:   PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE));

106:   /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
107:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
108:   PetscCall(FNSetType(f1,FNRATIONAL));
109:   coeffs[0] = -1.0; coeffs[1] = 0.0;
110:   PetscCall(FNRationalSetNumerator(f1,2,coeffs));

112:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
113:   PetscCall(FNSetType(f2,FNRATIONAL));
114:   coeffs[0] = 1.0;
115:   PetscCall(FNRationalSetNumerator(f2,1,coeffs));

117:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
118:   PetscCall(FNSetType(f3,FNEXP));
119:   PetscCall(FNSetScale(f3,-tau,1.0));

121:   /* Set the split operator */
122:   mats[0] = A;  funs[0] = f2;
123:   mats[1] = Id; funs[1] = f1;
124:   mats[2] = B;  funs[2] = f3;
125:   PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:                 Customize nonlinear solver; set runtime options
129:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

131:   PetscCall(NEPSetType(nep,NEPNLEIGS));
132:   PetscCall(NEPGetRG(nep,&rg));
133:   PetscCall(RGSetType(rg,RGINTERVAL));
134: #if defined(PETSC_USE_COMPLEX)
135:   PetscCall(RGIntervalSetEndpoints(rg,5,20,-0.001,0.001));
136: #else
137:   PetscCall(RGIntervalSetEndpoints(rg,5,20,-0.0,0.0));
138: #endif
139:   PetscCall(NEPSetTarget(nep,15.0));
140:   PetscCall(NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE));

142:   /*
143:      Set solver options. In particular, we must allocate sufficient
144:      storage for all eigenpairs that may converge (ncv). This is
145:      application-dependent.
146:   */
147:   mpd = 40;
148:   PetscCall(NEPSetDimensions(nep,2*mpd,3*mpd,mpd));
149:   PetscCall(NEPSetTolerances(nep,PETSC_DEFAULT,2000));
150:   PetscCall(PetscNew(&ctx));
151:   ctx->lastnconv = 0;
152:   ctx->nreps     = 0;
153:   PetscCall(NEPSetStoppingTestFunction(nep,MyStoppingTest,(void*)ctx,NULL));

155:   PetscCall(NEPSetFromOptions(nep));

157:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158:                       Solve the eigensystem
159:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

161:   PetscCall(NEPSolve(nep));

163:   /* show detailed info unless -terse option is given by user */
164:   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
165:   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
166:   PetscCall(NEPConvergedReasonView(nep,viewer));
167:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
168:   if (!terse) PetscCall(NEPErrorView(nep,NEP_ERROR_BACKWARD,viewer));
169:   PetscCall(PetscViewerPopFormat(viewer));

171:   PetscCall(NEPDestroy(&nep));
172:   PetscCall(MatDestroy(&Id));
173:   PetscCall(MatDestroy(&A));
174:   PetscCall(MatDestroy(&B));
175:   PetscCall(FNDestroy(&f1));
176:   PetscCall(FNDestroy(&f2));
177:   PetscCall(FNDestroy(&f3));
178:   PetscCall(PetscFree(ctx));
179:   PetscCall(SlepcFinalize());
180:   return 0;
181: }

183: /*
184:     Function for user-defined stopping test.

186:     Ignores the value of nev. It only takes into account the number of
187:     eigenpairs that have converged in recent outer iterations (restarts);
188:     if no new eigenvalues have converged in the last few restarts,
189:     we stop the iteration, assuming that no more eigenvalues are present
190:     inside the region.
191: */
192: PetscErrorCode MyStoppingTest(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ptr)
193: {
194:   CTX_DELAY      *ctx = (CTX_DELAY*)ptr;

196:   PetscFunctionBeginUser;
197:   /* check usual termination conditions, but ignoring the case nconv>=nev */
198:   PetscCall(NEPStoppingBasic(nep,its,max_it,nconv,PETSC_MAX_INT,reason,NULL));
199:   if (*reason==NEP_CONVERGED_ITERATING) {
200:     /* check if nconv is the same as before */
201:     if (nconv==ctx->lastnconv) ctx->nreps++;
202:     else {
203:       ctx->lastnconv = nconv;
204:       ctx->nreps     = 0;
205:     }
206:     /* check if no eigenvalues converged in last 10 restarts */
207:     if (nconv && ctx->nreps>10) *reason = NEP_CONVERGED_USER;
208:   }
209:   PetscFunctionReturn(PETSC_SUCCESS);
210: }

212: /*TEST

214:    test:
215:       suffix: 1
216:       args: -terse

218: TEST*/