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[] = "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 | 22 | int main(int argc,char **argv) | |
37 | { | ||
38 | 22 | NEP nep; | |
39 | 22 | PEP pep; | |
40 | 22 | Mat Id,A,B; | |
41 | 22 | FN f1,f2,f3; | |
42 | 22 | RG rg; | |
43 | 22 | Mat mats[3]; | |
44 | 22 | FN funs[3]; | |
45 | 22 | PetscScalar coeffs[2],b; | |
46 | 22 | PetscInt n=128,nev,Istart,Iend,i,deg; | |
47 | 22 | PetscReal tau=0.001,h,a=20,xi,tol; | |
48 | 22 | PetscBool terse; | |
49 | |||
50 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
22 | PetscFunctionBeginUser; |
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.
|
22 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
52 |
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.
|
22 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
53 |
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.
|
22 | PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL)); |
54 |
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.
|
22 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau)); |
55 | 22 | h = PETSC_PI/(PetscReal)(n+1); | |
56 | |||
57 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
58 | Create a standalone PEP and RG objects with appropriate settings | ||
59 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
60 | |||
61 |
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.
|
22 | PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep)); |
62 |
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.
|
22 | PetscCall(PEPSetType(pep,PEPTOAR)); |
63 |
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.
|
22 | PetscCall(PEPSetFromOptions(pep)); |
64 | |||
65 |
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.
|
22 | PetscCall(RGCreate(PETSC_COMM_WORLD,&rg)); |
66 |
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.
|
22 | PetscCall(RGSetType(rg,RGINTERVAL)); |
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.
|
22 | PetscCall(RGIntervalSetEndpoints(rg,5,20,-0.1,0.1)); |
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.
|
22 | PetscCall(RGSetFromOptions(rg)); |
69 | |||
70 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
71 | Create nonlinear eigensolver context | ||
72 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
73 | |||
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.
|
22 | PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep)); |
75 | |||
76 | /* Identity matrix */ | ||
77 |
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.
|
22 | PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id)); |
78 |
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.
|
22 | PetscCall(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE)); |
79 | |||
80 | /* A = 1/h^2*tridiag(1,-2,1) + a*I */ | ||
81 |
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.
|
22 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
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.
|
22 | PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
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.
|
22 | PetscCall(MatSetFromOptions(A)); |
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.
|
22 | PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); |
85 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2838 | for (i=Istart;i<Iend;i++) { |
86 |
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.
|
2816 | if (i>0) PetscCall(MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES)); |
87 |
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.
|
2816 | if (i<n-1) PetscCall(MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES)); |
88 |
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.
|
2816 | PetscCall(MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES)); |
89 | } | ||
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.
|
22 | PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); |
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.
|
22 | PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); |
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.
|
22 | PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE)); |
93 | |||
94 | /* B = diag(b(xi)) */ | ||
95 |
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.
|
22 | PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); |
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.
|
22 | PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
97 |
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.
|
22 | PetscCall(MatSetFromOptions(B)); |
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.
|
22 | PetscCall(MatGetOwnershipRange(B,&Istart,&Iend)); |
99 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2838 | for (i=Istart;i<Iend;i++) { |
100 | 2816 | xi = (i+1)*h; | |
101 | 2816 | b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI)); | |
102 |
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.
|
2816 | PetscCall(MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES)); |
103 | } | ||
104 |
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.
|
22 | PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); |
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.
|
22 | PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); |
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.
|
22 | PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE)); |
107 | |||
108 | /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */ | ||
109 |
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.
|
22 | PetscCall(FNCreate(PETSC_COMM_WORLD,&f1)); |
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.
|
22 | PetscCall(FNSetType(f1,FNRATIONAL)); |
111 | 22 | coeffs[0] = -1.0; coeffs[1] = 0.0; | |
112 |
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.
|
22 | PetscCall(FNRationalSetNumerator(f1,2,coeffs)); |
113 | |||
114 |
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.
|
22 | PetscCall(FNCreate(PETSC_COMM_WORLD,&f2)); |
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.
|
22 | PetscCall(FNSetType(f2,FNRATIONAL)); |
116 | 22 | coeffs[0] = 1.0; | |
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.
|
22 | PetscCall(FNRationalSetNumerator(f2,1,coeffs)); |
118 | |||
119 |
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.
|
22 | PetscCall(FNCreate(PETSC_COMM_WORLD,&f3)); |
120 |
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.
|
22 | PetscCall(FNSetType(f3,FNEXP)); |
121 |
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.
|
22 | PetscCall(FNSetScale(f3,-tau,1.0)); |
122 | |||
123 | /* Set the split operator */ | ||
124 | 22 | mats[0] = A; funs[0] = f2; | |
125 | 22 | mats[1] = Id; funs[1] = f1; | |
126 | 22 | mats[2] = B; funs[2] = f3; | |
127 |
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.
|
22 | PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN)); |
128 | |||
129 | /* Customize nonlinear solver; set runtime options */ | ||
130 |
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.
|
22 | PetscCall(NEPSetType(nep,NEPINTERPOL)); |
131 |
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.
|
22 | PetscCall(NEPSetRG(nep,rg)); |
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.
|
22 | PetscCall(NEPInterpolSetPEP(nep,pep)); |
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.
|
22 | PetscCall(NEPInterpolGetInterpolation(nep,&tol,°)); |
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.
|
22 | PetscCall(NEPInterpolSetInterpolation(nep,tol,deg+2)); /* increase degree of interpolation */ |
135 |
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.
|
22 | PetscCall(NEPSetFromOptions(nep)); |
136 | |||
137 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
138 | Solve the eigensystem | ||
139 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
140 | |||
141 |
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.
|
22 | PetscCall(NEPSolve(nep)); |
142 |
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.
|
22 | PetscCall(NEPGetDimensions(nep,&nev,NULL,NULL)); |
143 |
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.
|
22 | 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 |
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.
|
22 | PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
151 |
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.
|
22 | if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL)); |
152 | else { | ||
153 | ✗ | PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); | |
154 | ✗ | PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD)); | |
155 | ✗ | PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD)); | |
156 | ✗ | PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); | |
157 | } | ||
158 |
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.
|
22 | PetscCall(NEPDestroy(&nep)); |
159 |
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.
|
22 | PetscCall(PEPDestroy(&pep)); |
160 |
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.
|
22 | PetscCall(RGDestroy(&rg)); |
161 |
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.
|
22 | PetscCall(MatDestroy(&Id)); |
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.
|
22 | PetscCall(MatDestroy(&A)); |
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.
|
22 | PetscCall(MatDestroy(&B)); |
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.
|
22 | PetscCall(FNDestroy(&f1)); |
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.
|
22 | PetscCall(FNDestroy(&f2)); |
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.
|
22 | PetscCall(FNDestroy(&f3)); |
167 |
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.
|
22 | 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*/ | ||
194 |