Actual source code: test10.c

slepc-3.14.2 2021-02-01
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[] = "Tests multiple calls to NEPSolve(). Based on ex22.c.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions.\n"
 14:   "  -tau <tau>, where <tau> is the delay parameter.\n\n";

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

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

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

 24:    Discretization leads to a DDE of dimension n

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

 28:    which results in the nonlinear eigenproblem

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

 33: #include <slepcnep.h>

 35: /*
 36:    Check if computed eigenvectors have unit norm
 37: */
 38: PetscErrorCode CheckNormalizedVectors(NEP nep)
 39: {
 41:   PetscInt       i,nconv;
 42:   Mat            A;
 43:   Vec            xr,xi;
 44:   PetscReal      error=0.0,normr;
 45: #if !defined(PETSC_USE_COMPLEX)
 46:   PetscReal      normi;
 47: #endif

 50:   NEPGetConverged(nep,&nconv);
 51:   if (nconv>0) {
 52:     NEPGetSplitOperatorTerm(nep,0,&A,NULL);
 53:     MatCreateVecs(A,&xr,&xi);
 54:     for (i=0;i<nconv;i++) {
 55:       NEPGetEigenpair(nep,i,NULL,NULL,xr,xi);
 56: #if defined(PETSC_USE_COMPLEX)
 57:       VecNorm(xr,NORM_2,&normr);
 58:       error = PetscMax(error,PetscAbsReal(normr-PetscRealConstant(1.0)));
 59: #else
 60:       VecNormBegin(xr,NORM_2,&normr);
 61:       VecNormBegin(xi,NORM_2,&normi);
 62:       VecNormEnd(xr,NORM_2,&normr);
 63:       VecNormEnd(xi,NORM_2,&normi);
 64:       error = PetscMax(error,PetscAbsReal(SlepcAbsEigenvalue(normr,normi)-PetscRealConstant(1.0)));
 65: #endif
 66:     }
 67:     VecDestroy(&xr);
 68:     VecDestroy(&xi);
 69:     if (error>100*PETSC_MACHINE_EPSILON) {
 70:       PetscPrintf(PETSC_COMM_WORLD,"Vectors are not normalized. Error=%g\n",(double)error);
 71:     }
 72:   }
 73:   return(0);
 74: }

 76: int main(int argc,char **argv)
 77: {
 78:   NEP            nep;             /* nonlinear eigensolver context */
 79:   Mat            Id,A,B;          /* problem matrices */
 80:   FN             f1,f2,f3;        /* functions to define the nonlinear operator */
 81:   Mat            mats[3];
 82:   FN             funs[3];
 83:   PetscScalar    coeffs[2],b;
 84:   PetscInt       n=128,Istart,Iend,i;
 85:   PetscReal      tau=0.001,h,a=20,xi;
 86:   PetscBool      skipnorm=PETSC_FALSE;

 89:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 90:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 91:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 92:   PetscOptionsGetBool(NULL,NULL,"-skipnorm",&skipnorm,NULL);
 93:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%D, tau=%g\n\n",n,(double)tau);
 94:   h = PETSC_PI/(PetscReal)(n+1);

 96:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97:              Create functions that define the split operator
 98:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

100:   /* f1=-lambda */
101:   FNCreate(PETSC_COMM_WORLD,&f1);
102:   FNSetType(f1,FNRATIONAL);
103:   coeffs[0] = -1.0; coeffs[1] = 0.0;
104:   FNRationalSetNumerator(f1,2,coeffs);

106:   /* f2=1.0 */
107:   FNCreate(PETSC_COMM_WORLD,&f2);
108:   FNSetType(f2,FNRATIONAL);
109:   coeffs[0] = 1.0;
110:   FNRationalSetNumerator(f2,1,coeffs);

112:   /* f3=exp(-tau*lambda) */
113:   FNCreate(PETSC_COMM_WORLD,&f3);
114:   FNSetType(f3,FNEXP);
115:   FNSetScale(f3,-tau,1.0);

117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118:                          Create problem matrices
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

121:   /* Identity matrix */
122:   MatCreate(PETSC_COMM_WORLD,&Id);
123:   MatSetSizes(Id,PETSC_DECIDE,PETSC_DECIDE,n,n);
124:   MatSetFromOptions(Id);
125:   MatSetUp(Id);
126:   MatGetOwnershipRange(Id,&Istart,&Iend);
127:   for (i=Istart;i<Iend;i++) {
128:     MatSetValue(Id,i,i,1.0,INSERT_VALUES);
129:   }
130:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
131:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
132:   MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);

134:   /* A = 1/h^2*tridiag(1,-2,1) + a*I */
135:   MatCreate(PETSC_COMM_WORLD,&A);
136:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
137:   MatSetFromOptions(A);
138:   MatSetUp(A);
139:   MatGetOwnershipRange(A,&Istart,&Iend);
140:   for (i=Istart;i<Iend;i++) {
141:     if (i>0) { MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES); }
142:     if (i<n-1) { MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES); }
143:     MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
144:   }
145:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
146:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
147:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);

149:   /* B = diag(b(xi)) */
150:   MatCreate(PETSC_COMM_WORLD,&B);
151:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
152:   MatSetFromOptions(B);
153:   MatSetUp(B);
154:   MatGetOwnershipRange(B,&Istart,&Iend);
155:   for (i=Istart;i<Iend;i++) {
156:     xi = (i+1)*h;
157:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
158:     MatSetValue(B,i,i,b,INSERT_VALUES);
159:   }
160:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
161:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
162:   MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);

164:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165:               Create nonlinear eigensolver and set options
166:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

168:   NEPCreate(PETSC_COMM_WORLD,&nep);
169:   mats[0] = A;  funs[0] = f2;
170:   mats[1] = Id; funs[1] = f1;
171:   mats[2] = B;  funs[2] = f3;
172:   NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);
173:   NEPSetTolerances(nep,1e-9,PETSC_DEFAULT);
174:   NEPSetFromOptions(nep);

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

180:   NEPSolve(nep);
181:   NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
182:   if (!skipnorm) { CheckNormalizedVectors(nep); }

184:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185:                    Create problem matrices of size 2*n
186:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

188:   MatDestroy(&Id);
189:   MatDestroy(&A);
190:   MatDestroy(&B);
191:   n *= 2;
192:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%D, tau=%g\n\n",n,(double)tau);
193:   h = PETSC_PI/(PetscReal)(n+1);

195:   /* Identity matrix */
196:   MatCreate(PETSC_COMM_WORLD,&Id);
197:   MatSetSizes(Id,PETSC_DECIDE,PETSC_DECIDE,n,n);
198:   MatSetFromOptions(Id);
199:   MatSetUp(Id);
200:   MatGetOwnershipRange(Id,&Istart,&Iend);
201:   for (i=Istart;i<Iend;i++) {
202:     MatSetValue(Id,i,i,1.0,INSERT_VALUES);
203:   }
204:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
205:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
206:   MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);

208:   /* A = 1/h^2*tridiag(1,-2,1) + a*I */
209:   MatCreate(PETSC_COMM_WORLD,&A);
210:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
211:   MatSetFromOptions(A);
212:   MatSetUp(A);
213:   MatGetOwnershipRange(A,&Istart,&Iend);
214:   for (i=Istart;i<Iend;i++) {
215:     if (i>0) { MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES); }
216:     if (i<n-1) { MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES); }
217:     MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
218:   }
219:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
220:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
221:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);

223:   /* B = diag(b(xi)) */
224:   MatCreate(PETSC_COMM_WORLD,&B);
225:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
226:   MatSetFromOptions(B);
227:   MatSetUp(B);
228:   MatGetOwnershipRange(B,&Istart,&Iend);
229:   for (i=Istart;i<Iend;i++) {
230:     xi = (i+1)*h;
231:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
232:     MatSetValue(B,i,i,b,INSERT_VALUES);
233:   }
234:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
235:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
236:   MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);

238:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
239:        Solve again, calling NEPReset() since matrix size has changed
240:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

242:   NEPReset(nep);  /* if this is omitted, it will be called in NEPSetSplitOperators() */
243:   mats[0] = A;  funs[0] = f2;
244:   mats[1] = Id; funs[1] = f1;
245:   mats[2] = B;  funs[2] = f3;
246:   NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);
247:   NEPSolve(nep);
248:   NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
249:   if (!skipnorm) { CheckNormalizedVectors(nep); }

251:   NEPDestroy(&nep);
252:   MatDestroy(&Id);
253:   MatDestroy(&A);
254:   MatDestroy(&B);
255:   FNDestroy(&f1);
256:   FNDestroy(&f2);
257:   FNDestroy(&f3);
258:   SlepcFinalize();
259:   return ierr;
260: }

262: /*TEST

264:    testset:
265:       nsize: 2
266:       requires: !single
267:       output_file: output/test10_1.out
268:       test:
269:          suffix: 1
270:          args: -nep_type narnoldi -nep_target 0.55
271:       test:
272:          suffix: 1_rii
273:          args: -nep_type rii -nep_target 0.55 -nep_rii_hermitian
274:       test:
275:          suffix: 1_narnoldi
276:          args: -nep_type narnoldi -nep_target 0.55 -nep_narnoldi_lag_preconditioner 2
277:       test:
278:          suffix: 1_slp
279:          args: -nep_type slp -nep_slp_st_pc_type redundant
280:       test:
281:          suffix: 1_interpol
282:          args: -nep_type interpol -rg_type interval -rg_interval_endpoints .5,1,-.1,.1 -nep_target .7 -nep_interpol_st_pc_type redundant
283:       test:
284:          suffix: 1_narnoldi_sync
285:          args: -nep_type narnoldi -ds_parallel synchronized

287:    testset:
288:       args: -nep_nev 2 -rg_type interval -rg_interval_endpoints .5,15,-.1,.1 -nep_target .7
289:       requires: !single
290:       output_file: output/test10_2.out
291:       test:
292:          suffix: 2_interpol
293:          args: -nep_type interpol -nep_interpol_pep_type jd -nep_interpol_st_pc_type sor
294:       test:
295:          suffix: 2_nleigs
296:          args: -nep_type nleigs
297:          requires: complex
298:       test:
299:          suffix: 2_nleigs_real
300:          args: -nep_type nleigs -rg_interval_endpoints .5,15
301:          requires: !complex

303: TEST*/