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[] = "Illustrates the use of a user-defined stopping test.\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 | /* | ||
37 | User-defined routines | ||
38 | */ | ||
39 | PetscErrorCode MyStoppingTest(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*); | ||
40 | |||
41 | typedef struct { | ||
42 | PetscInt lastnconv; /* last value of nconv; used in stopping test */ | ||
43 | PetscInt nreps; /* number of repetitions of nconv; used in stopping test */ | ||
44 | } CTX_DELAY; | ||
45 | |||
46 | 10 | int main(int argc,char **argv) | |
47 | { | ||
48 | 10 | NEP nep; | |
49 | 10 | Mat Id,A,B; | |
50 | 10 | FN f1,f2,f3; | |
51 | 10 | RG rg; | |
52 | 10 | CTX_DELAY *ctx; | |
53 | 10 | Mat mats[3]; | |
54 | 10 | FN funs[3]; | |
55 | 10 | PetscScalar coeffs[2],b; | |
56 | 10 | PetscInt n=128,Istart,Iend,i,mpd; | |
57 | 10 | PetscReal tau=0.001,h,a=20,xi; | |
58 | 10 | PetscBool terse; | |
59 | 10 | PetscViewer viewer; | |
60 | |||
61 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBeginUser; |
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.
|
10 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
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.
|
10 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
64 |
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.
|
10 | PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL)); |
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.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau)); |
66 | 10 | h = PETSC_PI/(PetscReal)(n+1); | |
67 | |||
68 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
69 | Create nonlinear eigensolver context | ||
70 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
71 | |||
72 |
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.
|
10 | PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep)); |
73 | |||
74 | /* Identity matrix */ | ||
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.
|
10 | PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id)); |
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.
|
10 | PetscCall(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE)); |
77 | |||
78 | /* A = 1/h^2*tridiag(1,-2,1) + a*I */ | ||
79 |
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.
|
10 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
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.
|
10 | PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
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.
|
10 | PetscCall(MatSetFromOptions(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.
|
10 | PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); |
83 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1290 | for (i=Istart;i<Iend;i++) { |
84 |
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.
|
1280 | if (i>0) PetscCall(MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES)); |
85 |
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.
|
1280 | if (i<n-1) PetscCall(MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES)); |
86 |
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.
|
1280 | PetscCall(MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES)); |
87 | } | ||
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.
|
10 | PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); |
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.
|
10 | PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); |
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.
|
10 | PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE)); |
91 | |||
92 | /* B = diag(b(xi)) */ | ||
93 |
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.
|
10 | PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); |
94 |
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.
|
10 | PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
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.
|
10 | PetscCall(MatSetFromOptions(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.
|
10 | PetscCall(MatGetOwnershipRange(B,&Istart,&Iend)); |
97 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1290 | for (i=Istart;i<Iend;i++) { |
98 | 1280 | xi = (i+1)*h; | |
99 | 1280 | b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI)); | |
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.
|
1280 | PetscCall(MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES)); |
101 | } | ||
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.
|
10 | PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); |
103 |
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.
|
10 | PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); |
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.
|
10 | PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE)); |
105 | |||
106 | /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */ | ||
107 |
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.
|
10 | PetscCall(FNCreate(PETSC_COMM_WORLD,&f1)); |
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.
|
10 | PetscCall(FNSetType(f1,FNRATIONAL)); |
109 | 10 | coeffs[0] = -1.0; coeffs[1] = 0.0; | |
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.
|
10 | PetscCall(FNRationalSetNumerator(f1,2,coeffs)); |
111 | |||
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.
|
10 | PetscCall(FNCreate(PETSC_COMM_WORLD,&f2)); |
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.
|
10 | PetscCall(FNSetType(f2,FNRATIONAL)); |
114 | 10 | coeffs[0] = 1.0; | |
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.
|
10 | PetscCall(FNRationalSetNumerator(f2,1,coeffs)); |
116 | |||
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.
|
10 | PetscCall(FNCreate(PETSC_COMM_WORLD,&f3)); |
118 |
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.
|
10 | PetscCall(FNSetType(f3,FNEXP)); |
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.
|
10 | PetscCall(FNSetScale(f3,-tau,1.0)); |
120 | |||
121 | /* Set the split operator */ | ||
122 | 10 | mats[0] = A; funs[0] = f2; | |
123 | 10 | mats[1] = Id; funs[1] = f1; | |
124 | 10 | mats[2] = B; funs[2] = f3; | |
125 |
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.
|
10 | PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN)); |
126 | |||
127 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
128 | Customize nonlinear solver; set runtime options | ||
129 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
130 | |||
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.
|
10 | PetscCall(NEPSetType(nep,NEPNLEIGS)); |
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.
|
10 | PetscCall(NEPGetRG(nep,&rg)); |
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.
|
10 | PetscCall(RGSetType(rg,RGINTERVAL)); |
134 | #if defined(PETSC_USE_COMPLEX) | ||
135 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(RGIntervalSetEndpoints(rg,5,20,-0.001,0.001)); |
136 | #else | ||
137 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(RGIntervalSetEndpoints(rg,5,20,-0.0,0.0)); |
138 | #endif | ||
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.
|
10 | PetscCall(NEPSetTarget(nep,15.0)); |
140 |
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.
|
10 | PetscCall(NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE)); |
141 | |||
142 | /* | ||
143 | Set solver options. In particular, we must allocate sufficient | ||
144 | storage for all eigenpairs that may converge (ncv). This is | ||
145 | application-dependent. | ||
146 | */ | ||
147 | 10 | mpd = 40; | |
148 |
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.
|
10 | PetscCall(NEPSetDimensions(nep,2*mpd,3*mpd,mpd)); |
149 |
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.
|
10 | PetscCall(NEPSetTolerances(nep,PETSC_CURRENT,2000)); |
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.
|
10 | PetscCall(PetscNew(&ctx)); |
151 | 10 | ctx->lastnconv = 0; | |
152 | 10 | ctx->nreps = 0; | |
153 |
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.
|
10 | PetscCall(NEPSetStoppingTestFunction(nep,MyStoppingTest,(void*)ctx,NULL)); |
154 | |||
155 |
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.
|
10 | PetscCall(NEPSetFromOptions(nep)); |
156 | |||
157 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
158 | Solve the eigensystem | ||
159 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
160 | |||
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.
|
10 | PetscCall(NEPSolve(nep)); |
162 | |||
163 | /* show detailed info unless -terse option is given by user */ | ||
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.
|
10 | PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer)); |
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.
|
10 | PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL)); |
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.
|
10 | PetscCall(NEPConvergedReasonView(nep,viewer)); |
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.
|
10 | PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
168 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | if (!terse) PetscCall(NEPErrorView(nep,NEP_ERROR_BACKWARD,viewer)); |
169 |
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.
|
10 | PetscCall(PetscViewerPopFormat(viewer)); |
170 | |||
171 |
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.
|
10 | PetscCall(NEPDestroy(&nep)); |
172 |
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.
|
10 | PetscCall(MatDestroy(&Id)); |
173 |
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.
|
10 | PetscCall(MatDestroy(&A)); |
174 |
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.
|
10 | PetscCall(MatDestroy(&B)); |
175 |
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.
|
10 | PetscCall(FNDestroy(&f1)); |
176 |
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.
|
10 | PetscCall(FNDestroy(&f2)); |
177 |
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.
|
10 | PetscCall(FNDestroy(&f3)); |
178 |
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.
|
10 | PetscCall(PetscFree(ctx)); |
179 |
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.
|
10 | PetscCall(SlepcFinalize()); |
180 | return 0; | ||
181 | } | ||
182 | |||
183 | /* | ||
184 | Function for user-defined stopping test. | ||
185 | |||
186 | Ignores the value of nev. It only takes into account the number of | ||
187 | eigenpairs that have converged in recent outer iterations (restarts); | ||
188 | if no new eigenvalues have converged in the last few restarts, | ||
189 | we stop the iteration, assuming that no more eigenvalues are present | ||
190 | inside the region. | ||
191 | */ | ||
192 | 120 | PetscErrorCode MyStoppingTest(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ptr) | |
193 | { | ||
194 | 120 | CTX_DELAY *ctx = (CTX_DELAY*)ptr; | |
195 | |||
196 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
120 | PetscFunctionBeginUser; |
197 | /* check usual termination conditions, but ignoring the case nconv>=nev */ | ||
198 |
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.
|
120 | PetscCall(NEPStoppingBasic(nep,its,max_it,nconv,PETSC_INT_MAX,reason,NULL)); |
199 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
120 | if (*reason==NEP_CONVERGED_ITERATING) { |
200 | /* check if nconv is the same as before */ | ||
201 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
120 | if (nconv==ctx->lastnconv) ctx->nreps++; |
202 | else { | ||
203 | 10 | ctx->lastnconv = nconv; | |
204 | 10 | ctx->nreps = 0; | |
205 | } | ||
206 | /* check if no eigenvalues converged in last 10 restarts */ | ||
207 |
3/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
120 | if (nconv && ctx->nreps>10) *reason = NEP_CONVERGED_USER; |
208 | } | ||
209 |
5/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
24 | PetscFunctionReturn(PETSC_SUCCESS); |
210 | } | ||
211 | |||
212 | /*TEST | ||
213 | |||
214 | test: | ||
215 | suffix: 1 | ||
216 | args: -terse | ||
217 | |||
218 | TEST*/ | ||
219 |