Actual source code: test15.c
slepc-3.17.0 2022-03-31
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: SlepcInitialize(&argc,&argv,(char*)0,help);
62: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
63: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
64: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau);
65: h = PETSC_PI/(PetscReal)(n+1);
67: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: Create nonlinear eigensolver context
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: NEPCreate(PETSC_COMM_WORLD,&nep);
73: /* Identity matrix */
74: MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id);
75: MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);
77: /* A = 1/h^2*tridiag(1,-2,1) + a*I */
78: MatCreate(PETSC_COMM_WORLD,&A);
79: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
80: MatSetFromOptions(A);
81: MatSetUp(A);
82: MatGetOwnershipRange(A,&Istart,&Iend);
83: for (i=Istart;i<Iend;i++) {
84: if (i>0) MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES);
85: if (i<n-1) MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES);
86: MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
87: }
88: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
89: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
90: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
92: /* B = diag(b(xi)) */
93: MatCreate(PETSC_COMM_WORLD,&B);
94: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
95: MatSetFromOptions(B);
96: MatSetUp(B);
97: MatGetOwnershipRange(B,&Istart,&Iend);
98: for (i=Istart;i<Iend;i++) {
99: xi = (i+1)*h;
100: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
101: MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES);
102: }
103: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
104: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
105: MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);
107: /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
108: FNCreate(PETSC_COMM_WORLD,&f1);
109: FNSetType(f1,FNRATIONAL);
110: coeffs[0] = -1.0; coeffs[1] = 0.0;
111: FNRationalSetNumerator(f1,2,coeffs);
113: FNCreate(PETSC_COMM_WORLD,&f2);
114: FNSetType(f2,FNRATIONAL);
115: coeffs[0] = 1.0;
116: FNRationalSetNumerator(f2,1,coeffs);
118: FNCreate(PETSC_COMM_WORLD,&f3);
119: FNSetType(f3,FNEXP);
120: FNSetScale(f3,-tau,1.0);
122: /* Set the split operator */
123: mats[0] = A; funs[0] = f2;
124: mats[1] = Id; funs[1] = f1;
125: mats[2] = B; funs[2] = f3;
126: NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);
128: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129: Customize nonlinear solver; set runtime options
130: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132: NEPSetType(nep,NEPNLEIGS);
133: NEPGetRG(nep,&rg);
134: RGSetType(rg,RGINTERVAL);
135: #if defined(PETSC_USE_COMPLEX)
136: RGIntervalSetEndpoints(rg,5,20,-0.001,0.001);
137: #else
138: RGIntervalSetEndpoints(rg,5,20,-0.0,0.0);
139: #endif
140: NEPSetTarget(nep,15.0);
141: NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE);
143: /*
144: Set solver options. In particular, we must allocate sufficient
145: storage for all eigenpairs that may converge (ncv). This is
146: application-dependent.
147: */
148: mpd = 40;
149: NEPSetDimensions(nep,2*mpd,3*mpd,mpd);
150: NEPSetTolerances(nep,PETSC_DEFAULT,2000);
151: PetscNew(&ctx);
152: ctx->lastnconv = 0;
153: ctx->nreps = 0;
154: NEPSetStoppingTestFunction(nep,MyStoppingTest,(void*)ctx,NULL);
156: NEPSetFromOptions(nep);
158: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159: Solve the eigensystem
160: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
162: NEPSolve(nep);
164: /* show detailed info unless -terse option is given by user */
165: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
166: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
167: NEPConvergedReasonView(nep,viewer);
168: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
169: if (!terse) NEPErrorView(nep,NEP_ERROR_BACKWARD,viewer);
170: PetscViewerPopFormat(viewer);
172: NEPDestroy(&nep);
173: MatDestroy(&Id);
174: MatDestroy(&A);
175: MatDestroy(&B);
176: FNDestroy(&f1);
177: FNDestroy(&f2);
178: FNDestroy(&f3);
179: PetscFree(ctx);
180: SlepcFinalize();
181: return 0;
182: }
184: /*
185: Function for user-defined stopping test.
187: Ignores the value of nev. It only takes into account the number of
188: eigenpairs that have converged in recent outer iterations (restarts);
189: if no new eigenvalues have converged in the last few restarts,
190: we stop the iteration, assuming that no more eigenvalues are present
191: inside the region.
192: */
193: PetscErrorCode MyStoppingTest(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ptr)
194: {
195: CTX_DELAY *ctx = (CTX_DELAY*)ptr;
198: /* check usual termination conditions, but ignoring the case nconv>=nev */
199: NEPStoppingBasic(nep,its,max_it,nconv,PETSC_MAX_INT,reason,NULL);
200: if (*reason==NEP_CONVERGED_ITERATING) {
201: /* check if nconv is the same as before */
202: if (nconv==ctx->lastnconv) ctx->nreps++;
203: else {
204: ctx->lastnconv = nconv;
205: ctx->nreps = 0;
206: }
207: /* check if no eigenvalues converged in last 10 restarts */
208: if (nconv && ctx->nreps>10) *reason = NEP_CONVERGED_USER;
209: }
210: PetscFunctionReturn(0);
211: }
213: /*TEST
215: test:
216: suffix: 1
217: args: -terse
219: TEST*/