Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Test the INTERPOL solver with a user-provided PEP.\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";
16 :
17 : /*
18 : Solve parabolic partial differential equation with time delay tau
19 :
20 : u_t = u_xx + a*u(t) + b*u(t-tau)
21 : u(0,t) = u(pi,t) = 0
22 :
23 : with a = 20 and b(x) = -4.1+x*(1-exp(x-pi)).
24 :
25 : Discretization leads to a DDE of dimension n
26 :
27 : -u' = A*u(t) + B*u(t-tau)
28 :
29 : which results in the nonlinear eigenproblem
30 :
31 : (-lambda*I + A + exp(-tau*lambda)*B)*u = 0
32 : */
33 :
34 : #include <slepcnep.h>
35 :
36 1 : int main(int argc,char **argv)
37 : {
38 1 : NEP nep;
39 1 : PEP pep;
40 1 : Mat Id,A,B;
41 1 : FN f1,f2,f3;
42 1 : RG rg;
43 1 : Mat mats[3];
44 1 : FN funs[3];
45 1 : PetscScalar coeffs[2],b;
46 1 : PetscInt n=128,nev,Istart,Iend,i,deg;
47 1 : PetscReal tau=0.001,h,a=20,xi,tol;
48 1 : PetscBool terse;
49 :
50 1 : PetscFunctionBeginUser;
51 1 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
52 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
53 1 : PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
54 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau));
55 1 : h = PETSC_PI/(PetscReal)(n+1);
56 :
57 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58 : Create a standalone PEP and RG objects with appropriate settings
59 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60 :
61 1 : PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
62 1 : PetscCall(PEPSetType(pep,PEPTOAR));
63 1 : PetscCall(PEPSetFromOptions(pep));
64 :
65 1 : PetscCall(RGCreate(PETSC_COMM_WORLD,&rg));
66 1 : PetscCall(RGSetType(rg,RGINTERVAL));
67 1 : PetscCall(RGIntervalSetEndpoints(rg,5,20,-0.1,0.1));
68 1 : PetscCall(RGSetFromOptions(rg));
69 :
70 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71 : Create nonlinear eigensolver context
72 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73 :
74 1 : PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
75 :
76 : /* Identity matrix */
77 1 : PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id));
78 1 : PetscCall(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE));
79 :
80 : /* A = 1/h^2*tridiag(1,-2,1) + a*I */
81 1 : PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
82 1 : PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
83 1 : PetscCall(MatSetFromOptions(A));
84 1 : PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
85 129 : for (i=Istart;i<Iend;i++) {
86 128 : if (i>0) PetscCall(MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES));
87 128 : if (i<n-1) PetscCall(MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES));
88 128 : PetscCall(MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES));
89 : }
90 1 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
91 1 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
92 1 : PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
93 :
94 : /* B = diag(b(xi)) */
95 1 : PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
96 1 : PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
97 1 : PetscCall(MatSetFromOptions(B));
98 1 : PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
99 129 : for (i=Istart;i<Iend;i++) {
100 128 : xi = (i+1)*h;
101 128 : b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
102 128 : PetscCall(MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES));
103 : }
104 1 : PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
105 1 : PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
106 1 : PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE));
107 :
108 : /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
109 1 : PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
110 1 : PetscCall(FNSetType(f1,FNRATIONAL));
111 1 : coeffs[0] = -1.0; coeffs[1] = 0.0;
112 1 : PetscCall(FNRationalSetNumerator(f1,2,coeffs));
113 :
114 1 : PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
115 1 : PetscCall(FNSetType(f2,FNRATIONAL));
116 1 : coeffs[0] = 1.0;
117 1 : PetscCall(FNRationalSetNumerator(f2,1,coeffs));
118 :
119 1 : PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
120 1 : PetscCall(FNSetType(f3,FNEXP));
121 1 : PetscCall(FNSetScale(f3,-tau,1.0));
122 :
123 : /* Set the split operator */
124 1 : mats[0] = A; funs[0] = f2;
125 1 : mats[1] = Id; funs[1] = f1;
126 1 : mats[2] = B; funs[2] = f3;
127 1 : PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
128 :
129 : /* Customize nonlinear solver; set runtime options */
130 1 : PetscCall(NEPSetType(nep,NEPINTERPOL));
131 1 : PetscCall(NEPSetRG(nep,rg));
132 1 : PetscCall(NEPInterpolSetPEP(nep,pep));
133 1 : PetscCall(NEPInterpolGetInterpolation(nep,&tol,°));
134 1 : PetscCall(NEPInterpolSetInterpolation(nep,tol,deg+2)); /* increase degree of interpolation */
135 1 : PetscCall(NEPSetFromOptions(nep));
136 :
137 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138 : Solve the eigensystem
139 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140 :
141 1 : PetscCall(NEPSolve(nep));
142 1 : PetscCall(NEPGetDimensions(nep,&nev,NULL,NULL));
143 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
144 :
145 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146 : Display solution and clean up
147 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148 :
149 : /* show detailed info unless -terse option is given by user */
150 1 : PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
151 1 : if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
152 : else {
153 0 : PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
154 0 : PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
155 0 : PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
156 0 : PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
157 : }
158 1 : PetscCall(NEPDestroy(&nep));
159 1 : PetscCall(PEPDestroy(&pep));
160 1 : PetscCall(RGDestroy(&rg));
161 1 : PetscCall(MatDestroy(&Id));
162 1 : PetscCall(MatDestroy(&A));
163 1 : PetscCall(MatDestroy(&B));
164 1 : PetscCall(FNDestroy(&f1));
165 1 : PetscCall(FNDestroy(&f2));
166 1 : PetscCall(FNDestroy(&f3));
167 1 : PetscCall(SlepcFinalize());
168 : return 0;
169 : }
170 :
171 : /*TEST
172 :
173 : testset:
174 : args: -nep_nev 3 -nep_target 5 -terse
175 : output_file: output/test5_1.out
176 : filter: sed -e "s/[+-]0\.0*i//g"
177 : requires: !single
178 : test:
179 : suffix: 1
180 : test:
181 : suffix: 2_cuda
182 : args: -mat_type aijcusparse
183 : requires: cuda
184 : test:
185 : suffix: 2_hip
186 : args: -mat_type aijhipsparse
187 : requires: hip
188 : test:
189 : suffix: 3
190 : args: -nep_view_values draw
191 : requires: x
192 :
193 : TEST*/
|