Line | Branch | Exec | Source |
---|---|---|---|
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[] = "Delay differential equation.\n\n" | ||
12 | "The command line options are:\n" | ||
13 | " -n <n>, where <n> = number of grid subdivisions.\n" | ||
14 | " -tau <tau>, where <tau> is the delay parameter.\n\n"; | ||
15 | |||
16 | /* | ||
17 | Solve parabolic partial differential equation with time delay tau | ||
18 | |||
19 | u_t = u_xx + a*u(t) + b*u(t-tau) | ||
20 | u(0,t) = u(pi,t) = 0 | ||
21 | |||
22 | with a = 20 and b(x) = -4.1+x*(1-exp(x-pi)). | ||
23 | |||
24 | Discretization leads to a DDE of dimension n | ||
25 | |||
26 | -u' = A*u(t) + B*u(t-tau) | ||
27 | |||
28 | which results in the nonlinear eigenproblem | ||
29 | |||
30 | (-lambda*I + A + exp(-tau*lambda)*B)*u = 0 | ||
31 | */ | ||
32 | |||
33 | #include <slepcnep.h> | ||
34 | |||
35 | 163 | int main(int argc,char **argv) | |
36 | { | ||
37 | 163 | NEP nep; /* nonlinear eigensolver context */ | |
38 | 163 | Mat Id,A,B; /* problem matrices */ | |
39 | 163 | FN f1,f2,f3; /* functions to define the nonlinear operator */ | |
40 | 163 | Mat mats[3]; | |
41 | 163 | FN funs[3]; | |
42 | 163 | PetscScalar coeffs[2],b; | |
43 | 163 | PetscInt n=128,nev,Istart,Iend,i; | |
44 | 163 | PetscReal tau=0.001,h,a=20,xi; | |
45 | 163 | PetscBool terse; | |
46 | |||
47 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
163 | PetscFunctionBeginUser; |
48 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
50 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL)); |
51 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau)); |
52 | 163 | h = PETSC_PI/(PetscReal)(n+1); | |
53 | |||
54 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
55 | Create nonlinear eigensolver context | ||
56 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
57 | |||
58 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep)); |
59 | |||
60 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
61 | Create problem matrices and coefficient functions. Pass them to NEP | ||
62 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
63 | |||
64 | /* | ||
65 | Identity matrix | ||
66 | */ | ||
67 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id)); |
68 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE)); |
69 | |||
70 | /* | ||
71 | A = 1/h^2*tridiag(1,-2,1) + a*I | ||
72 | */ | ||
73 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
74 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
75 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatSetFromOptions(A)); |
76 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); |
77 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
51747 | for (i=Istart;i<Iend;i++) { |
78 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
51584 | if (i>0) PetscCall(MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES)); |
79 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
51584 | if (i<n-1) PetscCall(MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES)); |
80 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
51584 | PetscCall(MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES)); |
81 | } | ||
82 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); |
83 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); |
84 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE)); |
85 | |||
86 | /* | ||
87 | B = diag(b(xi)) | ||
88 | */ | ||
89 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); |
90 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
91 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatSetFromOptions(B)); |
92 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatGetOwnershipRange(B,&Istart,&Iend)); |
93 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
51747 | for (i=Istart;i<Iend;i++) { |
94 | 51584 | xi = (i+1)*h; | |
95 | 51584 | b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI)); | |
96 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
51584 | PetscCall(MatSetValue(B,i,i,b,INSERT_VALUES)); |
97 | } | ||
98 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); |
99 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); |
100 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE)); |
101 | |||
102 | /* | ||
103 | Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) | ||
104 | */ | ||
105 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNCreate(PETSC_COMM_WORLD,&f1)); |
106 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNSetType(f1,FNRATIONAL)); |
107 | 163 | coeffs[0] = -1.0; coeffs[1] = 0.0; | |
108 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNRationalSetNumerator(f1,2,coeffs)); |
109 | |||
110 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNCreate(PETSC_COMM_WORLD,&f2)); |
111 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNSetType(f2,FNRATIONAL)); |
112 | 163 | coeffs[0] = 1.0; | |
113 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNRationalSetNumerator(f2,1,coeffs)); |
114 | |||
115 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNCreate(PETSC_COMM_WORLD,&f3)); |
116 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNSetType(f3,FNEXP)); |
117 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNSetScale(f3,-tau,1.0)); |
118 | |||
119 | /* | ||
120 | Set the split operator. Note that A is passed first so that | ||
121 | SUBSET_NONZERO_PATTERN can be used | ||
122 | */ | ||
123 | 163 | mats[0] = A; funs[0] = f2; | |
124 | 163 | mats[1] = Id; funs[1] = f1; | |
125 | 163 | mats[2] = B; funs[2] = f3; | |
126 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN)); |
127 | |||
128 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
129 | Customize nonlinear solver; set runtime options | ||
130 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
131 | |||
132 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(NEPSetTolerances(nep,1e-9,PETSC_CURRENT)); |
133 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(NEPSetDimensions(nep,1,PETSC_DETERMINE,PETSC_DETERMINE)); |
134 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(NEPRIISetLagPreconditioner(nep,0)); |
135 | |||
136 | /* | ||
137 | Set solver parameters at runtime | ||
138 | */ | ||
139 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(NEPSetFromOptions(nep)); |
140 | |||
141 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
142 | Solve the eigensystem | ||
143 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
144 | |||
145 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(NEPSolve(nep)); |
146 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(NEPGetDimensions(nep,&nev,NULL,NULL)); |
147 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev)); |
148 | |||
149 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
150 | Display solution and clean up | ||
151 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
152 | |||
153 | /* show detailed info unless -terse option is given by user */ | ||
154 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
155 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
163 | if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL)); |
156 | else { | ||
157 | ✗ | PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); | |
158 | ✗ | PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD)); | |
159 | ✗ | PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD)); | |
160 | ✗ | PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); | |
161 | } | ||
162 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(NEPDestroy(&nep)); |
163 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatDestroy(&Id)); |
164 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatDestroy(&A)); |
165 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(MatDestroy(&B)); |
166 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNDestroy(&f1)); |
167 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNDestroy(&f2)); |
168 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
163 | PetscCall(FNDestroy(&f3)); |
169 |
3/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
163 | PetscCall(SlepcFinalize()); |
170 | return 0; | ||
171 | } | ||
172 | |||
173 | /*TEST | ||
174 | |||
175 | testset: | ||
176 | suffix: 1 | ||
177 | args: -nep_type {{rii slp narnoldi}} -terse | ||
178 | filter: sed -e "s/[+-]0\.0*i//g" | ||
179 | requires: !single | ||
180 | |||
181 | test: | ||
182 | suffix: 1_ciss | ||
183 | args: -nep_type ciss -nep_ciss_extraction {{ritz hankel caa}} -rg_type ellipse -rg_ellipse_center 10 -rg_ellipse_radius 9.5 -nep_ncv 24 -terse | ||
184 | requires: complex !single | ||
185 | |||
186 | test: | ||
187 | suffix: 2 | ||
188 | args: -nep_type interpol -nep_interpol_pep_extract {{none norm residual}} -rg_type interval -rg_interval_endpoints 5,20,-.1,.1 -nep_nev 3 -nep_target 5 -terse | ||
189 | filter: sed -e "s/[+-]0\.0*i//g" | ||
190 | requires: !single | ||
191 | |||
192 | testset: | ||
193 | args: -n 512 -nep_target 10 -nep_nev 3 -terse | ||
194 | filter: sed -e "s/[+-]0\.0*i//g" | ||
195 | requires: !single | ||
196 | output_file: output/ex22_3.out | ||
197 | test: | ||
198 | suffix: 3 | ||
199 | args: -nep_type {{rii slp narnoldi}} | ||
200 | test: | ||
201 | suffix: 3_simpleu | ||
202 | args: -nep_type {{rii slp narnoldi}} -nep_deflation_simpleu | ||
203 | test: | ||
204 | suffix: 3_slp_thres | ||
205 | args: -nep_type slp -nep_slp_deflation_threshold 1e-8 | ||
206 | test: | ||
207 | suffix: 3_rii_thres | ||
208 | args: -nep_type rii -nep_rii_deflation_threshold 1e-8 | ||
209 | |||
210 | test: | ||
211 | suffix: 4 | ||
212 | args: -nep_type interpol -rg_type interval -rg_interval_endpoints 5,20,-.1,.1 -nep_nev 3 -nep_target 5 -terse -nep_monitor draw::draw_lg | ||
213 | filter: sed -e "s/[+-]0\.0*i//g" | ||
214 | requires: x !single | ||
215 | output_file: output/ex22_2.out | ||
216 | |||
217 | TEST*/ | ||
218 |