Actual source code: test10.c

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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() with different matrix size.\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"
 15:   "  -split <0/1>, to select the split form in the problem definition (enabled by default).\n";

 17: /* Based on ex22.c (delay) */

 19: #include <slepcnep.h>

 21: /*
 22:    User-defined application context
 23: */
 24: typedef struct {
 25:   PetscScalar tau;
 26:   PetscReal   a;
 27: } ApplicationCtx;

 29: /*
 30:    Create problem matrices in split form
 31: */
 32: PetscErrorCode BuildSplitMatrices(PetscInt n,PetscReal a,Mat *Id,Mat *A,Mat *B)
 33: {
 35:   PetscInt       i,Istart,Iend;
 36:   PetscReal      h,xi;
 37:   PetscScalar    b;

 40:   h = PETSC_PI/(PetscReal)(n+1);

 42:   /* Identity matrix */
 43:   MatCreate(PETSC_COMM_WORLD,Id);
 44:   MatSetSizes(*Id,PETSC_DECIDE,PETSC_DECIDE,n,n);
 45:   MatSetFromOptions(*Id);
 46:   MatSetUp(*Id);
 47:   MatGetOwnershipRange(*Id,&Istart,&Iend);
 48:   for (i=Istart;i<Iend;i++) {
 49:     MatSetValue(*Id,i,i,1.0,INSERT_VALUES);
 50:   }
 51:   MatAssemblyBegin(*Id,MAT_FINAL_ASSEMBLY);
 52:   MatAssemblyEnd(*Id,MAT_FINAL_ASSEMBLY);
 53:   MatSetOption(*Id,MAT_HERMITIAN,PETSC_TRUE);

 55:   /* A = 1/h^2*tridiag(1,-2,1) + a*I */
 56:   MatCreate(PETSC_COMM_WORLD,A);
 57:   MatSetSizes(*A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 58:   MatSetFromOptions(*A);
 59:   MatSetUp(*A);
 60:   MatGetOwnershipRange(*A,&Istart,&Iend);
 61:   for (i=Istart;i<Iend;i++) {
 62:     if (i>0) { MatSetValue(*A,i,i-1,1.0/(h*h),INSERT_VALUES); }
 63:     if (i<n-1) { MatSetValue(*A,i,i+1,1.0/(h*h),INSERT_VALUES); }
 64:     MatSetValue(*A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
 65:   }
 66:   MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
 67:   MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
 68:   MatSetOption(*A,MAT_HERMITIAN,PETSC_TRUE);

 70:   /* B = diag(b(xi)) */
 71:   MatCreate(PETSC_COMM_WORLD,B);
 72:   MatSetSizes(*B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 73:   MatSetFromOptions(*B);
 74:   MatSetUp(*B);
 75:   MatGetOwnershipRange(*B,&Istart,&Iend);
 76:   for (i=Istart;i<Iend;i++) {
 77:     xi = (i+1)*h;
 78:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
 79:     MatSetValue(*B,i,i,b,INSERT_VALUES);
 80:   }
 81:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
 82:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
 83:   MatSetOption(*B,MAT_HERMITIAN,PETSC_TRUE);
 84:   return(0);
 85: }

 87: /*
 88:    Compute Function matrix  T(lambda)
 89: */
 90: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
 91: {
 93:   ApplicationCtx *user = (ApplicationCtx*)ctx;
 94:   PetscInt       i,n,Istart,Iend;
 95:   PetscReal      h,xi;
 96:   PetscScalar    b;

 99:   MatGetSize(fun,&n,NULL);
100:   h = PETSC_PI/(PetscReal)(n+1);
101:   MatGetOwnershipRange(fun,&Istart,&Iend);
102:   for (i=Istart;i<Iend;i++) {
103:     if (i>0) { MatSetValue(fun,i,i-1,1.0/(h*h),INSERT_VALUES); }
104:     if (i<n-1) { MatSetValue(fun,i,i+1,1.0/(h*h),INSERT_VALUES); }
105:     xi = (i+1)*h;
106:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
107:     MatSetValue(fun,i,i,-lambda-2.0/(h*h)+user->a+PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES);
108:   }
109:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
110:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
111:   if (fun != B) {
112:     MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
113:     MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
114:   }
115:   return(0);
116: }

118: /*
119:    Compute Jacobian matrix  T'(lambda)
120: */
121: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
122: {
124:   ApplicationCtx *user = (ApplicationCtx*)ctx;
125:   PetscInt       i,n,Istart,Iend;
126:   PetscReal      h,xi;
127:   PetscScalar    b;

130:   MatGetSize(jac,&n,NULL);
131:   h = PETSC_PI/(PetscReal)(n+1);
132:   MatGetOwnershipRange(jac,&Istart,&Iend);
133:   for (i=Istart;i<Iend;i++) {
134:     xi = (i+1)*h;
135:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
136:     MatSetValue(jac,i,i,-1.0-user->tau*PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES);
137:   }
138:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
139:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
140:   return(0);
141: }

143: int main(int argc,char **argv)
144: {
145:   NEP            nep;             /* nonlinear eigensolver context */
146:   Mat            Id,A,B,J,F;      /* problem matrices */
147:   FN             f1,f2,f3;        /* functions to define the nonlinear operator */
148:   ApplicationCtx ctx;             /* user-defined context */
149:   Mat            mats[3];
150:   FN             funs[3];
151:   PetscScalar    coeffs[2];
152:   PetscInt       n=128;
153:   PetscReal      tau=0.001,a=20;
154:   PetscBool      split=PETSC_TRUE;

157:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
158:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
159:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
160:   PetscOptionsGetBool(NULL,NULL,"-split",&split,NULL);
161:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%D, tau=%g\n\n",n,(double)tau);

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

167:   NEPCreate(PETSC_COMM_WORLD,&nep);
168:   NEPSetTolerances(nep,1e-9,PETSC_DEFAULT);

170:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171:                       First solve
172:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

174:   if (split) {
175:     BuildSplitMatrices(n,a,&Id,&A,&B);
176:     /* f1=-lambda */
177:     FNCreate(PETSC_COMM_WORLD,&f1);
178:     FNSetType(f1,FNRATIONAL);
179:     coeffs[0] = -1.0; coeffs[1] = 0.0;
180:     FNRationalSetNumerator(f1,2,coeffs);
181:     /* f2=1.0 */
182:     FNCreate(PETSC_COMM_WORLD,&f2);
183:     FNSetType(f2,FNRATIONAL);
184:     coeffs[0] = 1.0;
185:     FNRationalSetNumerator(f2,1,coeffs);
186:     /* f3=exp(-tau*lambda) */
187:     FNCreate(PETSC_COMM_WORLD,&f3);
188:     FNSetType(f3,FNEXP);
189:     FNSetScale(f3,-tau,1.0);
190:     mats[0] = A;  funs[0] = f2;
191:     mats[1] = Id; funs[1] = f1;
192:     mats[2] = B;  funs[2] = f3;
193:     NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);
194:   } else {
195:     /* callback form  */
196:     ctx.tau = tau;
197:     ctx.a   = a;
198:     MatCreate(PETSC_COMM_WORLD,&F);
199:     MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
200:     MatSetFromOptions(F);
201:     MatSeqAIJSetPreallocation(F,3,NULL);
202:     MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
203:     MatSetUp(F);
204:     NEPSetFunction(nep,F,F,FormFunction,&ctx);
205:     MatCreate(PETSC_COMM_WORLD,&J);
206:     MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
207:     MatSetFromOptions(J);
208:     MatSeqAIJSetPreallocation(J,3,NULL);
209:     MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
210:     MatSetUp(J);
211:     NEPSetJacobian(nep,J,FormJacobian,&ctx);
212:   }

214:   /* Set solver parameters at runtime */
215:   NEPSetFromOptions(nep);

217:   /* Solve the eigensystem */
218:   NEPSolve(nep);
219:   NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);

221:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222:                Second solve, with problem matrices of size 2*n
223:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

225:   n *= 2;
226:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%D, tau=%g\n\n",n,(double)tau);
227:   if (split) {
228:     MatDestroy(&Id);
229:     MatDestroy(&A);
230:     MatDestroy(&B);
231:     BuildSplitMatrices(n,a,&Id,&A,&B);
232:     mats[0] = A;
233:     mats[1] = Id;
234:     mats[2] = B;
235:     NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);
236:   } else {
237:     /* callback form  */
238:     MatDestroy(&F);
239:     MatDestroy(&J);
240:     MatCreate(PETSC_COMM_WORLD,&F);
241:     MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
242:     MatSetFromOptions(F);
243:     MatSeqAIJSetPreallocation(F,3,NULL);
244:     MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
245:     MatSetUp(F);
246:     NEPSetFunction(nep,F,F,FormFunction,&ctx);
247:     MatCreate(PETSC_COMM_WORLD,&J);
248:     MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
249:     MatSetFromOptions(J);
250:     MatSeqAIJSetPreallocation(J,3,NULL);
251:     MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
252:     MatSetUp(J);
253:     NEPSetJacobian(nep,J,FormJacobian,&ctx);
254:   }

256:   /* Solve the eigensystem */
257:   NEPSolve(nep);
258:   NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);

260:   NEPDestroy(&nep);
261:   if (split) {
262:     MatDestroy(&Id);
263:     MatDestroy(&A);
264:     MatDestroy(&B);
265:     FNDestroy(&f1);
266:     FNDestroy(&f2);
267:     FNDestroy(&f3);
268:   } else {
269:     MatDestroy(&F);
270:     MatDestroy(&J);
271:   }
272:   SlepcFinalize();
273:   return ierr;
274: }

276: /*TEST

278:    testset:
279:       nsize: 2
280:       requires: !single
281:       output_file: output/test10_1.out
282:       test:
283:          suffix: 1
284:          args: -nep_type narnoldi -nep_target 0.55
285:       test:
286:          suffix: 1_rii
287:          args: -nep_type rii -nep_target 0.55 -nep_rii_hermitian -split {{0 1}}
288:       test:
289:          suffix: 1_narnoldi
290:          args: -nep_type narnoldi -nep_target 0.55 -nep_narnoldi_lag_preconditioner 2
291:       test:
292:          suffix: 1_slp
293:          args: -nep_type slp -nep_slp_st_pc_type redundant -split {{0 1}}
294:       test:
295:          suffix: 1_interpol
296:          args: -nep_type interpol -rg_type interval -rg_interval_endpoints .5,1,-.1,.1 -nep_target .7 -nep_interpol_st_pc_type redundant
297:       test:
298:          suffix: 1_narnoldi_sync
299:          args: -nep_type narnoldi -ds_parallel synchronized

301:    testset:
302:       args: -nep_nev 2 -rg_type interval -rg_interval_endpoints .5,15,-.1,.1 -nep_target .7
303:       requires: !single
304:       output_file: output/test10_2.out
305:       filter: sed -e "s/[+-]0\.0*i//g"
306:       test:
307:          suffix: 2_interpol
308:          args: -nep_type interpol -nep_interpol_pep_type jd -nep_interpol_st_pc_type sor
309:       test:
310:          suffix: 2_nleigs
311:          args: -nep_type nleigs -split {{0 1}}
312:          requires: complex
313:       test:
314:          suffix: 2_nleigs_real
315:          args: -nep_type nleigs -rg_interval_endpoints .5,15 -split {{0 1}}
316:          requires: !complex

318: TEST*/