Actual source code: test6.c
slepc-3.18.2 2023-01-26
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[] = "Test the NArnoldi solver with a user-provided KSP.\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"
16: " -initv ... set an initial vector.\n\n";
18: /*
19: Solve parabolic partial differential equation with time delay tau
21: u_t = u_xx + a*u(t) + b*u(t-tau)
22: u(0,t) = u(pi,t) = 0
24: with a = 20 and b(x) = -4.1+x*(1-exp(x-pi)).
26: Discretization leads to a DDE of dimension n
28: -u' = A*u(t) + B*u(t-tau)
30: which results in the nonlinear eigenproblem
32: (-lambda*I + A + exp(-tau*lambda)*B)*u = 0
33: */
35: #include <slepcnep.h>
37: int main(int argc,char **argv)
38: {
39: NEP nep;
40: KSP ksp;
41: PC pc;
42: Mat Id,A,B,mats[3];
43: FN f1,f2,f3,funs[3];
44: Vec v0;
45: PetscScalar coeffs[2],b,*pv;
46: PetscInt n=128,nev,Istart,Iend,i,lag;
47: PetscReal tau=0.001,h,a=20,xi;
48: PetscBool terse,initv=PETSC_FALSE;
49: const char *prefix;
52: SlepcInitialize(&argc,&argv,(char*)0,help);
53: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
54: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
55: PetscOptionsGetBool(NULL,NULL,"-initv",&initv,NULL);
56: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau);
57: h = PETSC_PI/(PetscReal)(n+1);
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Create a standalone KSP with appropriate settings
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: KSPCreate(PETSC_COMM_WORLD,&ksp);
64: KSPSetType(ksp,KSPBCGS);
65: KSPGetPC(ksp,&pc);
66: PCSetType(pc,PCBJACOBI);
67: KSPSetFromOptions(ksp);
69: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
70: Create nonlinear eigensolver context
71: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73: NEPCreate(PETSC_COMM_WORLD,&nep);
75: /* Identity matrix */
76: MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id);
77: MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);
79: /* A = 1/h^2*tridiag(1,-2,1) + a*I */
80: MatCreate(PETSC_COMM_WORLD,&A);
81: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
82: MatSetFromOptions(A);
83: MatSetUp(A);
84: MatGetOwnershipRange(A,&Istart,&Iend);
85: for (i=Istart;i<Iend;i++) {
86: if (i>0) MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES);
87: if (i<n-1) MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES);
88: MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
89: }
90: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
91: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
92: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
94: /* B = diag(b(xi)) */
95: MatCreate(PETSC_COMM_WORLD,&B);
96: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
97: MatSetFromOptions(B);
98: MatSetUp(B);
99: MatGetOwnershipRange(B,&Istart,&Iend);
100: for (i=Istart;i<Iend;i++) {
101: xi = (i+1)*h;
102: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
103: MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES);
104: }
105: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
106: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
107: MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);
109: /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
110: FNCreate(PETSC_COMM_WORLD,&f1);
111: FNSetType(f1,FNRATIONAL);
112: coeffs[0] = -1.0; coeffs[1] = 0.0;
113: FNRationalSetNumerator(f1,2,coeffs);
115: FNCreate(PETSC_COMM_WORLD,&f2);
116: FNSetType(f2,FNRATIONAL);
117: coeffs[0] = 1.0;
118: FNRationalSetNumerator(f2,1,coeffs);
120: FNCreate(PETSC_COMM_WORLD,&f3);
121: FNSetType(f3,FNEXP);
122: FNSetScale(f3,-tau,1.0);
124: /* Set the split operator */
125: mats[0] = A; funs[0] = f2;
126: mats[1] = Id; funs[1] = f1;
127: mats[2] = B; funs[2] = f3;
128: NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);
130: /* Customize nonlinear solver; set runtime options */
131: NEPSetOptionsPrefix(nep,"check_");
132: NEPAppendOptionsPrefix(nep,"myprefix_");
133: NEPGetOptionsPrefix(nep,&prefix);
134: PetscPrintf(PETSC_COMM_WORLD,"NEP prefix is currently: %s\n\n",prefix);
135: NEPSetType(nep,NEPNARNOLDI);
136: NEPNArnoldiSetKSP(nep,ksp);
137: if (initv) { /* initial vector */
138: MatCreateVecs(A,&v0,NULL);
139: VecGetArray(v0,&pv);
140: for (i=Istart;i<Iend;i++) pv[i-Istart] = PetscSinReal((4.0*PETSC_PI*i)/n);
141: VecRestoreArray(v0,&pv);
142: NEPSetInitialSpace(nep,1,&v0);
143: VecDestroy(&v0);
144: }
145: NEPSetFromOptions(nep);
147: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: Solve the eigensystem
149: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151: NEPSolve(nep);
152: NEPGetDimensions(nep,&nev,NULL,NULL);
153: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
154: NEPNArnoldiGetLagPreconditioner(nep,&lag);
155: PetscPrintf(PETSC_COMM_WORLD," N-Arnoldi lag parameter: %" PetscInt_FMT "\n",lag);
157: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158: Display solution and clean up
159: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161: /* show detailed info unless -terse option is given by user */
162: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
163: if (terse) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
164: else {
165: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
166: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
167: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
168: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
169: }
170: NEPDestroy(&nep);
171: KSPDestroy(&ksp);
172: MatDestroy(&Id);
173: MatDestroy(&A);
174: MatDestroy(&B);
175: FNDestroy(&f1);
176: FNDestroy(&f2);
177: FNDestroy(&f3);
178: SlepcFinalize();
179: return 0;
180: }
182: /*TEST
184: test:
185: suffix: 1
186: args: -check_myprefix_nep_view -check_myprefix_nep_monitor_conv -initv -terse
187: filter: grep -v "tolerance" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g" -e "s/+0i//g"
188: requires: double
190: TEST*/