Actual source code: test17.c

slepc-3.17.0 2022-03-31
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[] = "Tests a user-provided preconditioner.\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:   "  -a <a>, where <a> is the coefficient that multiplies u in the equation.\n"
 16:   "  -split <0/1>, to select the split form in the problem definition (enabled by default).\n";

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

 20: #include <slepcnep.h>

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

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

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

 42:   /* Identity matrix */
 43:   MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,Id);
 44:   MatSetOption(*Id,MAT_HERMITIAN,PETSC_TRUE);

 46:   /* A = 1/h^2*tridiag(1,-2,1) + a*I */
 47:   MatCreate(PETSC_COMM_WORLD,A);
 48:   MatSetSizes(*A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 49:   MatSetFromOptions(*A);
 50:   MatSetUp(*A);
 51:   MatGetOwnershipRange(*A,&Istart,&Iend);
 52:   for (i=Istart;i<Iend;i++) {
 53:     if (i>0) MatSetValue(*A,i,i-1,1.0/(h*h),INSERT_VALUES);
 54:     if (i<n-1) MatSetValue(*A,i,i+1,1.0/(h*h),INSERT_VALUES);
 55:     MatSetValue(*A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
 56:   }
 57:   MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
 58:   MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
 59:   MatSetOption(*A,MAT_HERMITIAN,PETSC_TRUE);

 61:   /* B = diag(b(xi)) */
 62:   MatCreate(PETSC_COMM_WORLD,B);
 63:   MatSetSizes(*B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 64:   MatSetFromOptions(*B);
 65:   MatSetUp(*B);
 66:   MatGetOwnershipRange(*B,&Istart,&Iend);
 67:   for (i=Istart;i<Iend;i++) {
 68:     xi = (i+1)*h;
 69:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
 70:     MatSetValue(*B,i,i,b,INSERT_VALUES);
 71:   }
 72:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
 73:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
 74:   MatSetOption(*B,MAT_HERMITIAN,PETSC_TRUE);
 75:   PetscFunctionReturn(0);
 76: }

 78: /*
 79:    Create preconditioner matrices (only Ap=diag(A))
 80: */
 81: PetscErrorCode BuildSplitPreconditioner(PetscInt n,PetscReal a,Mat *Ap)
 82: {
 83:   PetscInt       i,Istart,Iend;
 84:   PetscReal      h;

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

 89:   /* Ap = diag(A) */
 90:   MatCreate(PETSC_COMM_WORLD,Ap);
 91:   MatSetSizes(*Ap,PETSC_DECIDE,PETSC_DECIDE,n,n);
 92:   MatSetFromOptions(*Ap);
 93:   MatSetUp(*Ap);
 94:   MatGetOwnershipRange(*Ap,&Istart,&Iend);
 95:   for (i=Istart;i<Iend;i++) MatSetValue(*Ap,i,i,-2.0/(h*h)+a,INSERT_VALUES);
 96:   MatAssemblyBegin(*Ap,MAT_FINAL_ASSEMBLY);
 97:   MatAssemblyEnd(*Ap,MAT_FINAL_ASSEMBLY);
 98:   MatSetOption(*Ap,MAT_HERMITIAN,PETSC_TRUE);
 99:   PetscFunctionReturn(0);
100: }

102: /*
103:    Compute Function matrix  T(lambda)
104: */
105: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
106: {
107:   ApplicationCtx *user = (ApplicationCtx*)ctx;
108:   PetscInt       i,n,Istart,Iend;
109:   PetscReal      h,xi;
110:   PetscScalar    b;

113:   MatGetSize(fun,&n,NULL);
114:   h = PETSC_PI/(PetscReal)(n+1);
115:   MatGetOwnershipRange(fun,&Istart,&Iend);
116:   for (i=Istart;i<Iend;i++) {
117:     if (i>0) MatSetValue(fun,i,i-1,1.0/(h*h),INSERT_VALUES);
118:     if (i<n-1) MatSetValue(fun,i,i+1,1.0/(h*h),INSERT_VALUES);
119:     xi = (i+1)*h;
120:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
121:     MatSetValue(fun,i,i,-lambda-2.0/(h*h)+user->a+PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES);
122:     if (B!=fun) MatSetValue(B,i,i,-lambda-2.0/(h*h)+user->a+PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES);
123:   }
124:   MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
125:   MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
126:   if (fun != B) {
127:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
128:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
129:   }
130:   PetscFunctionReturn(0);
131: }

133: /*
134:    Compute Jacobian matrix  T'(lambda)
135: */
136: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
137: {
138:   ApplicationCtx *user = (ApplicationCtx*)ctx;
139:   PetscInt       i,n,Istart,Iend;
140:   PetscReal      h,xi;
141:   PetscScalar    b;

144:   MatGetSize(jac,&n,NULL);
145:   h = PETSC_PI/(PetscReal)(n+1);
146:   MatGetOwnershipRange(jac,&Istart,&Iend);
147:   for (i=Istart;i<Iend;i++) {
148:     xi = (i+1)*h;
149:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
150:     MatSetValue(jac,i,i,-1.0-user->tau*PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES);
151:   }
152:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
153:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
154:   PetscFunctionReturn(0);
155: }

157: int main(int argc,char **argv)
158: {
159:   NEP            nep;             /* nonlinear eigensolver context */
160:   Mat            Id,A,B,Ap,J,F,P; /* problem matrices */
161:   FN             f1,f2,f3;        /* functions to define the nonlinear operator */
162:   ApplicationCtx ctx;             /* user-defined context */
163:   Mat            mats[3];
164:   FN             funs[3];
165:   PetscScalar    coeffs[2];
166:   PetscInt       n=128;
167:   PetscReal      tau=0.001,a=20;
168:   PetscBool      split=PETSC_TRUE;

170:   SlepcInitialize(&argc,&argv,(char*)0,help);
171:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
172:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
173:   PetscOptionsGetReal(NULL,NULL,"-a",&a,NULL);
174:   PetscOptionsGetBool(NULL,NULL,"-split",&split,NULL);
175:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g, a=%g\n\n",n,(double)tau,(double)a);

177:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178:               Create nonlinear eigensolver and solve the problem
179:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

181:   NEPCreate(PETSC_COMM_WORLD,&nep);
182:   if (split) {
183:     BuildSplitMatrices(n,a,&Id,&A,&B);
184:     /* f1=-lambda */
185:     FNCreate(PETSC_COMM_WORLD,&f1);
186:     FNSetType(f1,FNRATIONAL);
187:     coeffs[0] = -1.0; coeffs[1] = 0.0;
188:     FNRationalSetNumerator(f1,2,coeffs);
189:     /* f2=1.0 */
190:     FNCreate(PETSC_COMM_WORLD,&f2);
191:     FNSetType(f2,FNRATIONAL);
192:     coeffs[0] = 1.0;
193:     FNRationalSetNumerator(f2,1,coeffs);
194:     /* f3=exp(-tau*lambda) */
195:     FNCreate(PETSC_COMM_WORLD,&f3);
196:     FNSetType(f3,FNEXP);
197:     FNSetScale(f3,-tau,1.0);
198:     mats[0] = A;  funs[0] = f2;
199:     mats[1] = Id; funs[1] = f1;
200:     mats[2] = B;  funs[2] = f3;
201:     NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);
202:     BuildSplitPreconditioner(n,a,&Ap);
203:     mats[0] = Ap;
204:     mats[1] = Id;
205:     mats[2] = B;
206:     NEPSetSplitPreconditioner(nep,3,mats,SAME_NONZERO_PATTERN);
207:   } else {
208:     /* callback form  */
209:     ctx.tau = tau;
210:     ctx.a   = a;
211:     MatCreate(PETSC_COMM_WORLD,&F);
212:     MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
213:     MatSetFromOptions(F);
214:     MatSeqAIJSetPreallocation(F,3,NULL);
215:     MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
216:     MatSetUp(F);
217:     MatDuplicate(F,MAT_DO_NOT_COPY_VALUES,&P);
218:     NEPSetFunction(nep,F,P,FormFunction,&ctx);
219:     MatCreate(PETSC_COMM_WORLD,&J);
220:     MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
221:     MatSetFromOptions(J);
222:     MatSeqAIJSetPreallocation(J,3,NULL);
223:     MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
224:     MatSetUp(J);
225:     NEPSetJacobian(nep,J,FormJacobian,&ctx);
226:   }

228:   /* Set solver parameters at runtime */
229:   NEPSetFromOptions(nep);

231:   /* Solve the eigensystem */
232:   NEPSolve(nep);
233:   NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);

235:   NEPDestroy(&nep);
236:   if (split) {
237:     MatDestroy(&Id);
238:     MatDestroy(&A);
239:     MatDestroy(&B);
240:     MatDestroy(&Ap);
241:     FNDestroy(&f1);
242:     FNDestroy(&f2);
243:     FNDestroy(&f3);
244:   } else {
245:     MatDestroy(&F);
246:     MatDestroy(&P);
247:     MatDestroy(&J);
248:   }
249:   SlepcFinalize();
250:   return 0;
251: }

253: /*TEST

255:    testset:
256:       args: -a 90000 -nep_nev 2
257:       requires: double
258:       output_file: output/test17_1.out
259:       test:
260:          suffix: 1
261:          args: -nep_type slp -nep_two_sided {{0 1}} -split {{0 1}}

263:    testset:
264:       args: -nep_nev 2 -rg_type interval -rg_interval_endpoints .5,15,-.1,.1 -nep_target .7
265:       requires: !single
266:       output_file: output/test17_2.out
267:       filter: sed -e "s/[+-]0\.0*i//g"
268:       test:
269:          suffix: 2_interpol
270:          args: -nep_type interpol -nep_interpol_st_ksp_type bcgs -nep_interpol_st_pc_type sor -nep_tol 1e-6 -nep_interpol_st_ksp_rtol 1e-7
271:       test:
272:          suffix: 2_nleigs
273:          args: -nep_type nleigs -split {{0 1}}
274:          requires: complex
275:       test:
276:          suffix: 2_nleigs_real
277:          args: -nep_type nleigs -rg_interval_endpoints .5,15 -split {{0 1}} -nep_nleigs_ksp_type tfqmr
278:          requires: !complex

280:    testset:
281:       args: -nep_type ciss -rg_type ellipse -rg_ellipse_center 10 -rg_ellipse_radius 9.5 -rg_ellipse_vscale 0.1 -nep_ciss_ksp_type bcgs -nep_ciss_pc_type sor
282:       output_file: output/test17_3.out
283:       requires: complex !single
284:       test:
285:          suffix: 3
286:          args: -split {{0 1}}
287:       test:
288:          suffix: 3_par
289:          nsize: 2
290:          args: -nep_ciss_partitions 2

292: TEST*/