GCC Code Coverage Report


Directory: ./
File: src/pep/interface/pepsolve.c
Date: 2025-11-19 04:19:03
Exec Total Coverage
Lines: 188 194 96.9%
Functions: 10 10 100.0%
Branches: 480 1024 46.9%

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 PEP routines related to the solution process
12
13 References:
14
15 [1] C. Campos and J.E. Roman, "Parallel iterative refinement in
16 polynomial eigenvalue problems", Numer. Linear Algebra Appl.
17 23(4):730-745, 2016.
18 */
19
20 #include <slepc/private/pepimpl.h> /*I "slepcpep.h" I*/
21 #include <slepc/private/bvimpl.h>
22 #include <petscdraw.h>
23
24 static PetscBool cited = PETSC_FALSE;
25 static const char citation[] =
26 "@Article{slepc-pep-refine,\n"
27 " author = \"C. Campos and J. E. Roman\",\n"
28 " title = \"Parallel iterative refinement in polynomial eigenvalue problems\",\n"
29 " journal = \"Numer. Linear Algebra Appl.\",\n"
30 " volume = \"23\",\n"
31 " number = \"4\",\n"
32 " pages = \"730--745\",\n"
33 " year = \"2016,\"\n"
34 " doi = \"https://doi.org/10.1002/nla.2052\"\n"
35 "}\n";
36
37 16727 PetscErrorCode PEPComputeVectors(PEP pep)
38 {
39
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
16727 PetscFunctionBegin;
40
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
16727 PEPCheckSolved(pep,1);
41
8/10
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
16727 if (pep->state==PEP_STATE_SOLVED) PetscTryTypeMethod(pep,computevectors);
42 16727 pep->state = PEP_STATE_EIGENVECTORS;
43
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
16727 PetscFunctionReturn(PETSC_SUCCESS);
44 }
45
46 1756 PetscErrorCode PEPExtractVectors(PEP pep)
47 {
48
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1756 PetscFunctionBegin;
49
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
1756 PEPCheckSolved(pep,1);
50
7/10
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
1756 if (pep->state==PEP_STATE_SOLVED) PetscTryTypeMethod(pep,extractvectors);
51
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
342 PetscFunctionReturn(PETSC_SUCCESS);
52 }
53
54 /*@
55 PEPSolve - Solves the polynomial eigenproblem.
56
57 Collective
58
59 Input Parameter:
60 . pep - the polynomial eigensolver context
61
62 Options Database Keys:
63 + -pep_view - print information about the solver once the solve is complete
64 . -pep_view_pre - print information about the solver before the solve starts
65 . -pep_view_matk - view the coefficient matrix $A_k$ (replace `k` by an integer from 0 to `nmat`-1)
66 . -pep_view_vectors - view the computed eigenvectors
67 . -pep_view_values - view the computed eigenvalues
68 . -pep_converged_reason - print reason for convergence/divergence, and number of iterations
69 . -pep_error_absolute - print absolute errors of each eigenpair
70 . -pep_error_relative - print relative errors of each eigenpair
71 - -pep_error_backward - print backward errors of each eigenpair
72
73 Notes:
74 The problem matrices are specified with `PEPSetOperators()`.
75
76 `PEPSolve()` will return without generating an error regardless of whether
77 all requested solutions were computed or not. Call `PEPGetConverged()` to get the
78 actual number of computed solutions, and `PEPGetConvergedReason()` to determine if
79 the solver converged or failed and why.
80
81 All the command-line options listed above admit an optional argument specifying
82 the viewer type and options. For instance, use `-pep_view_mat0 binary:matrix0.bin`
83 to save the $A_0$ matrix to a binary file, `-pep_view_values draw` to draw the computed
84 eigenvalues graphically, or `-pep_error_relative :myerr.m:ascii_matlab` to save
85 the errors in a file that can be executed in Matlab.
86
87 Level: beginner
88
89 .seealso: [](ch:pep), `PEPCreate()`, `PEPSetUp()`, `PEPDestroy()`, `PEPSetTolerances()`, `PEPGetConverged()`, `PEPGetConvergedReason()`
90 @*/
91 1956 PetscErrorCode PEPSolve(PEP pep)
92 {
93 1956 PetscInt i,k;
94 1956 PetscBool flg,islinear;
95 1956 char str[16];
96
97
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1956 PetscFunctionBegin;
98
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
1956 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
99
2/14
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 2 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.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
1956 if (pep->state>=PEP_STATE_SOLVED) PetscFunctionReturn(PETSC_SUCCESS);
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.
1956 PetscCall(PetscLogEventBegin(PEP_Solve,pep,0,0,0));
101
102 /* call setup */
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.
1956 PetscCall(PEPSetUp(pep));
104 1956 pep->nconv = 0;
105 1956 pep->its = 0;
106
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1956 k = pep->lineariz? pep->ncv: pep->ncv*(pep->nmat-1);
107
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
133657 for (i=0;i<k;i++) {
108 131701 pep->eigr[i] = 0.0;
109 131701 pep->eigi[i] = 0.0;
110 131701 pep->errest[i] = 0.0;
111 131701 pep->perm[i] = i;
112 }
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.
1956 PetscCall(PEPViewFromOptions(pep,NULL,"-pep_view_pre"));
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.
1956 PetscCall(RGViewFromOptions(pep->rg,NULL,"-rg_view"));
115
116 /* Call solver */
117
5/10
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
1956 PetscUseTypeMethod(pep,solve);
118
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
1956 PetscCheck(pep->reason,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
119 1956 pep->state = PEP_STATE_SOLVED;
120
121 /* Only the first nconv columns contain useful information */
122
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.
1956 PetscCall(BVSetActiveColumns(pep->V,0,pep->nconv));
123
124
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.
1956 PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&islinear));
125
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1956 if (!islinear) {
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.
1543 PetscCall(STPostSolve(pep->st));
127 /* Map eigenvalues back to the original problem */
128
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.
1543 PetscCall(STGetTransform(pep->st,&flg));
129
9/10
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 6 times.
✓ Branch 5 taken 4 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✓ Branch 9 taken 1 times.
1543 if (flg) PetscTryTypeMethod(pep,backtransform);
130 }
131
132 #if !defined(PETSC_USE_COMPLEX)
133 /* reorder conjugate eigenvalues (positive imaginary first) */
134
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
5606 for (i=0;i<pep->nconv-1;i++) {
135
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
4664 if (pep->eigi[i] != 0) {
136
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
704 if (pep->eigi[i] < 0) {
137 175 pep->eigi[i] = -pep->eigi[i];
138 175 pep->eigi[i+1] = -pep->eigi[i+1];
139 /* the next correction only works with eigenvectors */
140
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.
175 PetscCall(PEPComputeVectors(pep));
141
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.
175 PetscCall(BVScaleColumn(pep->V,i+1,-1.0));
142 }
143 704 i++;
144 }
145 }
146 #endif
147
148
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.
1956 if (pep->refine!=PEP_REFINE_NONE) PetscCall(PetscCitationsRegister(citation,&cited));
149
150
4/6
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 10 times.
✗ Branch 5 not taken.
1956 if (pep->refine==PEP_REFINE_SIMPLE && pep->rits>0 && pep->nconv>0) {
151
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.
156 PetscCall(PEPComputeVectors(pep));
152
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.
156 PetscCall(PEPNewtonRefinementSimple(pep,&pep->rits,pep->rtol,pep->nconv));
153 }
154
155 /* sort eigenvalues according to pep->which parameter */
156
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.
1956 PetscCall(SlepcSortEigenvalues(pep->sc,pep->nconv,pep->eigr,pep->eigi,pep->perm));
157
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.
1956 PetscCall(PetscLogEventEnd(PEP_Solve,pep,0,0,0));
158
159 /* various viewers */
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.
1956 PetscCall(PEPViewFromOptions(pep,NULL,"-pep_view"));
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.
1956 PetscCall(PEPConvergedReasonViewFromOptions(pep));
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.
1956 PetscCall(PEPErrorViewFromOptions(pep));
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.
1956 PetscCall(PEPValuesViewFromOptions(pep));
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.
1956 PetscCall(PEPVectorsViewFromOptions(pep));
165
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
9302 for (i=0;i<pep->nmat;i++) {
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.
7346 PetscCall(PetscSNPrintf(str,sizeof(str),"-pep_view_mat%" PetscInt_FMT,i));
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.
7346 PetscCall(MatViewFromOptions(pep->A[i],(PetscObject)pep,str));
168 }
169
170 /* Remove the initial subspace */
171 1956 pep->nini = 0;
172
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
1956 PetscFunctionReturn(PETSC_SUCCESS);
173 }
174
175 /*@
176 PEPGetIterationNumber - Gets the current iteration number. If the
177 call to `PEPSolve()` is complete, then it returns the number of iterations
178 carried out by the solution method.
179
180 Not Collective
181
182 Input Parameter:
183 . pep - the polynomial eigensolver context
184
185 Output Parameter:
186 . its - number of iterations
187
188 Note:
189 During the $i$-th iteration this call returns $i-1$. If `PEPSolve()` is
190 complete, then parameter `its` contains either the iteration number at
191 which convergence was successfully reached, or failure was detected.
192 Call `PEPGetConvergedReason()` to determine if the solver converged or
193 failed and why.
194
195 Level: intermediate
196
197 .seealso: [](ch:pep), `PEPGetConvergedReason()`, `PEPSetTolerances()`
198 @*/
199 310 PetscErrorCode PEPGetIterationNumber(PEP pep,PetscInt *its)
200 {
201
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
310 PetscFunctionBegin;
202
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
310 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
203
2/8
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
310 PetscAssertPointer(its,2);
204 310 *its = pep->its;
205
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
310 PetscFunctionReturn(PETSC_SUCCESS);
206 }
207
208 /*@
209 PEPGetConverged - Gets the number of converged eigenpairs.
210
211 Not Collective
212
213 Input Parameter:
214 . pep - the polynomial eigensolver context
215
216 Output Parameter:
217 . nconv - number of converged eigenpairs
218
219 Notes:
220 This function should be called after `PEPSolve()` has finished.
221
222 The value `nconv` may be different from the number of requested solutions
223 `nev`, but not larger than `ncv`, see `PEPSetDimensions()`.
224
225 Level: beginner
226
227 .seealso: [](ch:pep), `PEPSetDimensions()`, `PEPSolve()`, `PEPGetEigenpair()`
228 @*/
229 924 PetscErrorCode PEPGetConverged(PEP pep,PetscInt *nconv)
230 {
231
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
924 PetscFunctionBegin;
232
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
924 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
233
2/8
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
924 PetscAssertPointer(nconv,2);
234
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
924 PEPCheckSolved(pep,1);
235 924 *nconv = pep->nconv;
236
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
924 PetscFunctionReturn(PETSC_SUCCESS);
237 }
238
239 /*@
240 PEPGetConvergedReason - Gets the reason why the `PEPSolve()` iteration was
241 stopped.
242
243 Not Collective
244
245 Input Parameter:
246 . pep - the polynomial eigensolver context
247
248 Output Parameter:
249 . reason - negative value indicates diverged, positive value converged, see
250 `PEPConvergedReason` for the possible values
251
252 Options Database Key:
253 . -pep_converged_reason - print reason for convergence/divergence, and number of iterations
254
255 Note:
256 If this routine is called before or doing the `PEPSolve()` the value of
257 `PEP_CONVERGED_ITERATING` is returned.
258
259 Level: intermediate
260
261 .seealso: [](ch:pep), `PEPSetTolerances()`, `PEPSolve()`, `PEPConvergedReason`
262 @*/
263 266 PetscErrorCode PEPGetConvergedReason(PEP pep,PEPConvergedReason *reason)
264 {
265
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
266 PetscFunctionBegin;
266
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
266 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
267
2/8
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
266 PetscAssertPointer(reason,2);
268
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
266 PEPCheckSolved(pep,1);
269 266 *reason = pep->reason;
270
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
266 PetscFunctionReturn(PETSC_SUCCESS);
271 }
272
273 /*@
274 PEPGetEigenpair - Gets the `i`-th solution of the eigenproblem as computed by
275 `PEPSolve()`. The solution consists in both the eigenvalue and the eigenvector.
276
277 Collective
278
279 Input Parameters:
280 + pep - the polynomial eigensolver context
281 - i - index of the solution
282
283 Output Parameters:
284 + eigr - real part of eigenvalue
285 . eigi - imaginary part of eigenvalue
286 . Vr - real part of eigenvector
287 - Vi - imaginary part of eigenvector
288
289 Notes:
290 It is allowed to pass `NULL` for `Vr` and `Vi`, if the eigenvector is not
291 required. Otherwise, the caller must provide valid `Vec` objects, i.e.,
292 they must be created by the calling program with e.g. `MatCreateVecs()`.
293
294 If the eigenvalue is real, then `eigi` and `Vi` are set to zero. If PETSc is
295 configured with complex scalars the eigenvalue is stored
296 directly in `eigr` (`eigi` is set to zero) and the eigenvector in `Vr` (`Vi` is
297 set to zero). In any case, the user can pass `NULL` in `Vr` or `Vi` if one of
298 them is not required.
299
300 The index `i` should be a value between 0 and `nconv`-1 (see `PEPGetConverged()`).
301 Eigenpairs are indexed according to the ordering criterion established
302 with `PEPSetWhichEigenpairs()`.
303
304 The eigenvector is normalized to have unit norm.
305
306 Level: beginner
307
308 .seealso: [](ch:pep), `PEPSolve()`, `PEPGetConverged()`, `PEPSetWhichEigenpairs()`
309 @*/
310 16386 PetscErrorCode PEPGetEigenpair(PEP pep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
311 {
312 16386 PetscInt k;
313
314
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
16386 PetscFunctionBegin;
315
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
16386 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
316
27/62
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
16386 PetscValidLogicalCollectiveInt(pep,i,2);
317
17/46
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✓ Branch 26 taken 2 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✗ Branch 34 not taken.
✓ Branch 35 taken 2 times.
✗ Branch 36 not taken.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✓ Branch 39 taken 2 times.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 42 not taken.
✓ Branch 43 taken 2 times.
✗ Branch 44 not taken.
✗ Branch 45 not taken.
16386 if (Vr) { PetscValidHeaderSpecific(Vr,VEC_CLASSID,5); PetscCheckSameComm(pep,1,Vr,5); }
318
17/46
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✓ Branch 26 taken 2 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✗ Branch 34 not taken.
✓ Branch 35 taken 2 times.
✗ Branch 36 not taken.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✓ Branch 39 taken 2 times.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 42 not taken.
✓ Branch 43 taken 2 times.
✗ Branch 44 not taken.
✗ Branch 45 not taken.
16386 if (Vi) { PetscValidHeaderSpecific(Vi,VEC_CLASSID,6); PetscCheckSameComm(pep,1,Vi,6); }
319
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
16386 PEPCheckSolved(pep,1);
320
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
16386 PetscCheck(i>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
321
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
16386 PetscCheck(i<pep->nconv,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see PEPGetConverged()");
322
323
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.
16386 PetscCall(PEPComputeVectors(pep));
324 16386 k = pep->perm[i];
325
326 /* eigenvalue */
327 #if defined(PETSC_USE_COMPLEX)
328
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
8320 if (eigr) *eigr = pep->eigr[k];
329
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
8320 if (eigi) *eigi = 0;
330 #else
331
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
8066 if (eigr) *eigr = pep->eigr[k];
332
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
8066 if (eigi) *eigi = pep->eigi[k];
333 #endif
334
335 /* eigenvector */
336
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.
16386 PetscCall(BV_GetEigenvector(pep->V,k,pep->eigi[k],Vr,Vi));
337
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
3198 PetscFunctionReturn(PETSC_SUCCESS);
338 }
339
340 /*@
341 PEPGetErrorEstimate - Returns the error estimate associated to the `i`-th
342 computed eigenpair.
343
344 Not Collective
345
346 Input Parameters:
347 + pep - the polynomial eigensolver context
348 - i - index of eigenpair
349
350 Output Parameter:
351 . errest - the error estimate
352
353 Note:
354 This is the error estimate used internally by the eigensolver. The actual
355 error bound can be computed with `PEPComputeError()`.
356
357 Level: advanced
358
359 .seealso: [](ch:pep), `PEPComputeError()`
360 @*/
361 43 PetscErrorCode PEPGetErrorEstimate(PEP pep,PetscInt i,PetscReal *errest)
362 {
363
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
43 PetscFunctionBegin;
364
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
43 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
365
2/8
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
43 PetscAssertPointer(errest,3);
366
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
43 PEPCheckSolved(pep,1);
367
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
43 PetscCheck(i>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
368
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
43 PetscCheck(i<pep->nconv,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see PEPGetConverged()");
369 43 *errest = pep->errest[pep->perm[i]];
370
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
43 PetscFunctionReturn(PETSC_SUCCESS);
371 }
372
373 /*
374 PEPComputeResidualNorm_Private - Computes the norm of the residual vector
375 associated with an eigenpair.
376
377 Input Parameters:
378 kr,ki - eigenvalue
379 xr,xi - eigenvector
380 z - array of 4 work vectors (z[2],z[3] not referenced in complex scalars)
381 */
382 11837 PetscErrorCode PEPComputeResidualNorm_Private(PEP pep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm)
383 {
384 11837 Mat *A=pep->A;
385 11837 PetscInt i,nmat=pep->nmat;
386 11837 PetscScalar t[20],*vals=t,*ivals=NULL;
387 11837 Vec u,w;
388 #if !defined(PETSC_USE_COMPLEX)
389 5437 Vec ui,wi;
390 5437 PetscReal ni;
391 5437 PetscBool imag;
392 5437 PetscScalar it[20];
393 #endif
394
395
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
11837 PetscFunctionBegin;
396 11837 u = z[0]; w = z[1];
397
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.
11837 PetscCall(VecSet(u,0.0));
398 #if !defined(PETSC_USE_COMPLEX)
399 5437 ui = z[2]; wi = z[3];
400 5437 ivals = it;
401 #endif
402
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
11837 if (nmat>20) {
403 PetscCall(PetscMalloc1(nmat,&vals));
404 #if !defined(PETSC_USE_COMPLEX)
405 PetscCall(PetscMalloc1(nmat,&ivals));
406 #endif
407 }
408
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.
11837 PetscCall(PEPEvaluateBasis(pep,kr,ki,vals,ivals));
409 #if !defined(PETSC_USE_COMPLEX)
410
3/4
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✗ Branch 3 not taken.
5437 if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON))
411 imag = PETSC_FALSE;
412 else {
413 1139 imag = PETSC_TRUE;
414
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.
1139 PetscCall(VecSet(ui,0.0));
415 }
416 #endif
417
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
49718 for (i=0;i<nmat;i++) {
418
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
37881 if (vals[i]!=0.0) {
419
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.
37881 PetscCall(MatMult(A[i],xr,w));
420
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.
37881 PetscCall(VecAXPY(u,vals[i],w));
421 }
422 #if !defined(PETSC_USE_COMPLEX)
423
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
16636 if (imag) {
424
3/4
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✗ Branch 3 not taken.
3697 if (ivals[i]!=0 || vals[i]!=0) {
425
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.
3697 PetscCall(MatMult(A[i],xi,wi));
426
1/8
✗ Branch 0 not taken.
✓ Branch 1 taken 5 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.
3697 if (vals[i]==0) PetscCall(MatMult(A[i],xr,w));
427 }
428
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
3697 if (ivals[i]!=0) {
429
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.
2558 PetscCall(VecAXPY(u,-ivals[i],wi));
430
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.
2558 PetscCall(VecAXPY(ui,ivals[i],w));
431 }
432
5/8
✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
16636 if (vals[i]!=0) PetscCall(VecAXPY(ui,vals[i],wi));
433 }
434 #endif
435 }
436
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.
11837 PetscCall(VecNorm(u,NORM_2,norm));
437 #if !defined(PETSC_USE_COMPLEX)
438
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
5437 if (imag) {
439
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.
1139 PetscCall(VecNorm(ui,NORM_2,&ni));
440 1139 *norm = SlepcAbsEigenvalue(*norm,ni);
441 }
442 #endif
443
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
11837 if (nmat>20) {
444 PetscCall(PetscFree(vals));
445 #if !defined(PETSC_USE_COMPLEX)
446 PetscCall(PetscFree(ivals));
447 #endif
448 }
449
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
2326 PetscFunctionReturn(PETSC_SUCCESS);
450 }
451
452 /*@
453 PEPComputeError - Computes the error (based on the residual norm) associated
454 with the `i`-th computed eigenpair.
455
456 Collective
457
458 Input Parameters:
459 + pep - the polynomial eigensolver context
460 . i - the solution index
461 - type - the type of error to compute, see `PEPErrorType`
462
463 Output Parameter:
464 . error - the error
465
466 Note:
467 The error can be computed in various ways, all of them based on the residual
468 norm $\|P(\lambda)x\|_2$ where $(\lambda,x)$ is the approximate eigenpair.
469
470 Level: beginner
471
472 .seealso: [](ch:pep), `PEPErrorType`, `PEPSolve()`, `PEPGetErrorEstimate()`
473 @*/
474 10572 PetscErrorCode PEPComputeError(PEP pep,PetscInt i,PEPErrorType type,PetscReal *error)
475 {
476 10572 Vec xr,xi,w[4];
477 10572 PetscScalar kr,ki;
478 10572 PetscReal t,z=0.0;
479 10572 PetscInt j;
480 10572 PetscBool flg;
481
482
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
10572 PetscFunctionBegin;
483
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
10572 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
484
27/62
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
10572 PetscValidLogicalCollectiveInt(pep,i,2);
485
27/62
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
10572 PetscValidLogicalCollectiveEnum(pep,type,3);
486
2/8
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
10572 PetscAssertPointer(error,4);
487
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
10572 PEPCheckSolved(pep,1);
488
489 /* allocate work vectors */
490 #if defined(PETSC_USE_COMPLEX)
491
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.
5510 PetscCall(PEPSetWorkVecs(pep,3));
492 5510 xi = NULL;
493 5510 w[2] = NULL;
494 5510 w[3] = NULL;
495 #else
496
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.
5062 PetscCall(PEPSetWorkVecs(pep,6));
497 5062 xi = pep->work[3];
498 5062 w[2] = pep->work[4];
499 5062 w[3] = pep->work[5];
500 #endif
501 10572 xr = pep->work[0];
502 10572 w[0] = pep->work[1];
503 10572 w[1] = pep->work[2];
504
505 /* compute residual norms */
506
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.
10572 PetscCall(PEPGetEigenpair(pep,i,&kr,&ki,xr,xi));
507
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.
10572 PetscCall(PEPComputeResidualNorm_Private(pep,kr,ki,xr,xi,w,error));
508
509 /* compute error */
510
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
10572 switch (type) {
511 case PEP_ERROR_ABSOLUTE:
512 break;
513 327 case PEP_ERROR_RELATIVE:
514 327 *error /= SlepcAbsEigenvalue(kr,ki);
515 327 break;
516 10125 case PEP_ERROR_BACKWARD:
517 /* initialization of matrix norms */
518
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
10125 if (!pep->nrma[pep->nmat-1]) {
519
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
6432 for (j=0;j<pep->nmat;j++) {
520
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.
4879 PetscCall(MatHasOperation(pep->A[j],MATOP_NORM,&flg));
521
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
4879 PetscCheck(flg,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
522
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.
4879 PetscCall(MatNorm(pep->A[j],NORM_INFINITY,&pep->nrma[j]));
523 }
524 }
525 10125 t = SlepcAbsEigenvalue(kr,ki);
526
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
41920 for (j=pep->nmat-1;j>=0;j--) {
527 31795 z = z*t+pep->nrma[j];
528 }
529 10125 *error /= z;
530 10125 break;
531 default:
532 SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
533 }
534
6/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 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
2073 PetscFunctionReturn(PETSC_SUCCESS);
535 }
536