Actual source code: test6.c
slepc-3.21.2 2024-09-25
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;
51: PetscFunctionBeginUser;
52: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
53: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
54: PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
55: PetscCall(PetscOptionsGetBool(NULL,NULL,"-initv",&initv,NULL));
56: PetscCall(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: PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp));
64: PetscCall(KSPSetType(ksp,KSPBCGS));
65: PetscCall(KSPGetPC(ksp,&pc));
66: PetscCall(PCSetType(pc,PCBJACOBI));
67: PetscCall(KSPSetFromOptions(ksp));
69: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
70: Create nonlinear eigensolver context
71: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73: PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
75: /* Identity matrix */
76: PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id));
77: PetscCall(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE));
79: /* A = 1/h^2*tridiag(1,-2,1) + a*I */
80: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
81: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
82: PetscCall(MatSetFromOptions(A));
83: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
84: for (i=Istart;i<Iend;i++) {
85: if (i>0) PetscCall(MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES));
86: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES));
87: PetscCall(MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES));
88: }
89: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
90: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
91: PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
93: /* B = diag(b(xi)) */
94: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
95: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
96: PetscCall(MatSetFromOptions(B));
97: PetscCall(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: PetscCall(MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES));
102: }
103: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
104: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
105: PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE));
107: /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
108: PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
109: PetscCall(FNSetType(f1,FNRATIONAL));
110: coeffs[0] = -1.0; coeffs[1] = 0.0;
111: PetscCall(FNRationalSetNumerator(f1,2,coeffs));
113: PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
114: PetscCall(FNSetType(f2,FNRATIONAL));
115: coeffs[0] = 1.0;
116: PetscCall(FNRationalSetNumerator(f2,1,coeffs));
118: PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
119: PetscCall(FNSetType(f3,FNEXP));
120: PetscCall(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: PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
128: /* Customize nonlinear solver; set runtime options */
129: PetscCall(NEPSetOptionsPrefix(nep,"check_"));
130: PetscCall(NEPAppendOptionsPrefix(nep,"myprefix_"));
131: PetscCall(NEPGetOptionsPrefix(nep,&prefix));
132: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"NEP prefix is currently: %s\n\n",prefix));
133: PetscCall(NEPSetType(nep,NEPNARNOLDI));
134: PetscCall(NEPNArnoldiSetKSP(nep,ksp));
135: if (initv) { /* initial vector */
136: PetscCall(MatCreateVecs(A,&v0,NULL));
137: PetscCall(VecGetArray(v0,&pv));
138: for (i=Istart;i<Iend;i++) pv[i-Istart] = PetscSinReal((4.0*PETSC_PI*i)/n);
139: PetscCall(VecRestoreArray(v0,&pv));
140: PetscCall(NEPSetInitialSpace(nep,1,&v0));
141: PetscCall(VecDestroy(&v0));
142: }
143: PetscCall(NEPSetFromOptions(nep));
145: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: Solve the eigensystem
147: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149: PetscCall(NEPSolve(nep));
150: PetscCall(NEPGetDimensions(nep,&nev,NULL,NULL));
151: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
152: PetscCall(NEPNArnoldiGetLagPreconditioner(nep,&lag));
153: PetscCall(PetscPrintf(PETSC_COMM_WORLD," N-Arnoldi lag parameter: %" PetscInt_FMT "\n",lag));
155: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156: Display solution and clean up
157: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159: /* show detailed info unless -terse option is given by user */
160: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
161: if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
162: else {
163: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
164: PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
165: PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
166: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
167: }
168: PetscCall(NEPDestroy(&nep));
169: PetscCall(KSPDestroy(&ksp));
170: PetscCall(MatDestroy(&Id));
171: PetscCall(MatDestroy(&A));
172: PetscCall(MatDestroy(&B));
173: PetscCall(FNDestroy(&f1));
174: PetscCall(FNDestroy(&f2));
175: PetscCall(FNDestroy(&f3));
176: PetscCall(SlepcFinalize());
177: return 0;
178: }
180: /*TEST
182: test:
183: suffix: 1
184: args: -check_myprefix_nep_view -check_myprefix_nep_monitor_conv -initv -terse
185: filter: grep -v "tolerance" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g" -e "s/+0i//g"
186: requires: double
188: TEST*/