Actual source code: test13.c
slepc-3.21.1 2024-04-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 NEPProjectOperator() function.\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";
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: int main(int argc,char **argv)
37: {
38: NEP nep;
39: Mat Id,A,B,mats[3];
40: FN f1,f2,f3,funs[3];
41: BV V;
42: DS ds;
43: Vec v;
44: PetscScalar coeffs[2],b,*M;
45: PetscInt n=32,Istart,Iend,i,j,k,nc;
46: PetscReal tau=0.001,h,a=20,xi;
48: PetscFunctionBeginUser;
49: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
50: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
51: PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
52: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n",n,(double)tau));
53: h = PETSC_PI/(PetscReal)(n+1);
55: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56: Create nonlinear eigensolver context
57: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59: PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
61: /* Identity matrix */
62: PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id));
63: PetscCall(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE));
65: /* A = 1/h^2*tridiag(1,-2,1) + a*I */
66: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
67: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
68: PetscCall(MatSetFromOptions(A));
69: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
70: for (i=Istart;i<Iend;i++) {
71: if (i>0) PetscCall(MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES));
72: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES));
73: PetscCall(MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES));
74: }
75: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
76: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
77: PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
79: /* B = diag(b(xi)) */
80: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
81: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
82: PetscCall(MatSetFromOptions(B));
83: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
84: for (i=Istart;i<Iend;i++) {
85: xi = (i+1)*h;
86: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
87: PetscCall(MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES));
88: }
89: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
90: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
91: PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE));
93: /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
94: PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
95: PetscCall(FNSetType(f1,FNRATIONAL));
96: coeffs[0] = -1.0; coeffs[1] = 0.0;
97: PetscCall(FNRationalSetNumerator(f1,2,coeffs));
99: PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
100: PetscCall(FNSetType(f2,FNRATIONAL));
101: coeffs[0] = 1.0;
102: PetscCall(FNRationalSetNumerator(f2,1,coeffs));
104: PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
105: PetscCall(FNSetType(f3,FNEXP));
106: PetscCall(FNSetScale(f3,-tau,1.0));
108: /* Set the split operator */
109: mats[0] = A; funs[0] = f2;
110: mats[1] = Id; funs[1] = f1;
111: mats[2] = B; funs[2] = f3;
112: PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
113: PetscCall(NEPSetType(nep,NEPNARNOLDI));
114: PetscCall(NEPSetFromOptions(nep));
116: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: Project the NEP
118: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120: PetscCall(NEPSetUp(nep));
121: PetscCall(NEPGetBV(nep,&V));
122: PetscCall(BVGetSizes(V,NULL,NULL,&nc));
123: for (i=0;i<nc;i++) {
124: PetscCall(BVGetColumn(V,i,&v));
125: PetscCall(VecSetValue(v,i,1.0,INSERT_VALUES));
126: PetscCall(VecAssemblyBegin(v));
127: PetscCall(VecAssemblyEnd(v));
128: PetscCall(BVRestoreColumn(V,i,&v));
129: }
130: PetscCall(NEPGetDS(nep,&ds));
131: PetscCall(DSSetType(ds,DSNEP));
132: PetscCall(DSNEPSetFN(ds,3,funs));
133: PetscCall(DSAllocate(ds,nc));
134: PetscCall(DSSetDimensions(ds,nc,0,0));
135: PetscCall(NEPProjectOperator(nep,0,nc));
137: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138: Display projected matrices and clean up
139: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141: for (k=0;k<3;k++) {
142: PetscCall(DSGetArray(ds,DSMatExtra[k],&M));
143: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMatrix E%" PetscInt_FMT " = \n",k));
144: for (i=0;i<nc;i++) {
145: for (j=0;j<nc;j++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %.5g",(double)PetscRealPart(M[i+j*nc])));
146: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
147: }
148: PetscCall(DSRestoreArray(ds,DSMatExtra[k],&M));
149: }
151: PetscCall(NEPDestroy(&nep));
152: PetscCall(MatDestroy(&Id));
153: PetscCall(MatDestroy(&A));
154: PetscCall(MatDestroy(&B));
155: PetscCall(FNDestroy(&f1));
156: PetscCall(FNDestroy(&f2));
157: PetscCall(FNDestroy(&f3));
158: PetscCall(SlepcFinalize());
159: return 0;
160: }
162: /*TEST
164: test:
165: suffix: 1
166: args: -nep_ncv 5
168: TEST*/