GCC Code Coverage Report


Directory: ./
File: src/pep/interface/pepsolve.c
Date: 2026-05-12 09:16:45
Exec Total Coverage
Lines: 188 193 97.4%
Functions: 10 10 100.0%
Branches: 480 1022 47.0%

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 16952 PetscErrorCode PEPComputeVectors(PEP pep)
38 {
39
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
16952 PetscFunctionBegin;
40
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
16952 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.
16952 if (pep->state==PEP_STATE_SOLVED) PetscTryTypeMethod(pep,computevectors);
42 16952 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.
16952 PetscFunctionReturn(PETSC_SUCCESS);
44 }
45
46 1780 PetscErrorCode PEPExtractVectors(PEP pep)
47 {
48
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1780 PetscFunctionBegin;
49
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
1780 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.
1780 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.
348 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 See `PetscObjectViewFromOptions()` for more details.
87
88 Level: beginner
89
90 .seealso: [](ch:pep), `PEPCreate()`, `PEPSetUp()`, `PEPDestroy()`, `PEPSetTolerances()`, `PEPGetConverged()`, `PEPGetConvergedReason()`
91 @*/
92 1980 PetscErrorCode PEPSolve(PEP pep)
93 {
94 1980 PetscInt i,k;
95 1980 PetscBool flg,islinear;
96 1980 char str[16];
97
98
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1980 PetscFunctionBegin;
99
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.
1980 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
100
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.
1980 if (pep->state>=PEP_STATE_SOLVED) PetscFunctionReturn(PETSC_SUCCESS);
101
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.
1980 PetscCall(PetscLogEventBegin(PEP_Solve,pep,0,0,0));
102
103 /* call setup */
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.
1980 PetscCall(PEPSetUp(pep));
105 1980 pep->nconv = 0;
106 1980 pep->its = 0;
107
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1980 k = pep->lineariz? pep->ncv: pep->ncv*(pep->nmat-1);
108
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
134517 for (i=0;i<k;i++) {
109 132537 pep->eigr[i] = 0.0;
110 132537 pep->eigi[i] = 0.0;
111 132537 pep->errest[i] = 0.0;
112 132537 pep->perm[i] = i;
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.
1980 PetscCall(PEPViewFromOptions(pep,NULL,"-pep_view_pre"));
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.
1980 PetscCall(RGViewFromOptions(pep->rg,NULL,"-rg_view"));
116
117 /* Call solver */
118
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.
1980 PetscUseTypeMethod(pep,solve);
119
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
1980 PetscCheck(pep->reason,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
120 1980 pep->state = PEP_STATE_SOLVED;
121
122 /* Only the first nconv columns contain useful information */
123
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.
1980 PetscCall(BVSetActiveColumns(pep->V,0,pep->nconv));
124
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.
1980 PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&islinear));
126
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1980 if (!islinear) {
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.
1567 PetscCall(STPostSolve(pep->st));
128 /* Map eigenvalues back to the original problem */
129
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.
1567 PetscCall(STGetTransform(pep->st,&flg));
130
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.
1567 if (flg) PetscTryTypeMethod(pep,backtransform);
131 }
132
133 #if !defined(PETSC_USE_COMPLEX)
134 /* reorder conjugate eigenvalues (positive imaginary first) */
135
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
5727 for (i=0;i<pep->nconv-1;i++) {
136
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
4774 if (pep->eigi[i] != 0) {
137
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
715 if (pep->eigi[i] < 0) {
138 175 pep->eigi[i] = -pep->eigi[i];
139 175 pep->eigi[i+1] = -pep->eigi[i+1];
140 /* the next correction only works with eigenvectors */
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(PEPComputeVectors(pep));
142
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));
143 }
144 715 i++;
145 }
146 }
147 #endif
148
149
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.
1980 if (pep->refine!=PEP_REFINE_NONE) PetscCall(PetscCitationsRegister(citation,&cited));
150
151
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.
1980 if (pep->refine==PEP_REFINE_SIMPLE && pep->rits>0 && pep->nconv>0) {
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(PEPComputeVectors(pep));
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.
156 PetscCall(PEPNewtonRefinementSimple(pep,&pep->rits,pep->rtol,pep->nconv));
154 }
155
156 /* sort eigenvalues according to pep->which parameter */
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.
1980 PetscCall(SlepcSortEigenvalues(pep->sc,pep->nconv,pep->eigr,pep->eigi,pep->perm));
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.
1980 PetscCall(PetscLogEventEnd(PEP_Solve,pep,0,0,0));
159
160 /* various viewers */
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.
1980 PetscCall(PEPViewFromOptions(pep,NULL,"-pep_view"));
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.
1980 PetscCall(PEPConvergedReasonViewFromOptions(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.
1980 PetscCall(PEPErrorViewFromOptions(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.
1980 PetscCall(PEPValuesViewFromOptions(pep));
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.
1980 PetscCall(PEPVectorsViewFromOptions(pep));
166
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
9400 for (i=0;i<pep->nmat;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.
7420 PetscCall(PetscSNPrintf(str,sizeof(str),"-pep_view_mat%" PetscInt_FMT,i));
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.
7420 PetscCall(MatViewFromOptions(pep->A[i],(PetscObject)pep,str));
169 }
170
171 /* Remove the initial subspace */
172 1980 pep->nini = 0;
173
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.
1980 PetscFunctionReturn(PETSC_SUCCESS);
174 }
175
176 /*@
177 PEPGetIterationNumber - Gets the current iteration number. If the
178 call to `PEPSolve()` is complete, then it returns the number of iterations
179 carried out by the solution method.
180
181 Not Collective
182
183 Input Parameter:
184 . pep - the polynomial eigensolver context
185
186 Output Parameter:
187 . its - number of iterations
188
189 Note:
190 During the $i$-th iteration this call returns $i-1$. If `PEPSolve()` is
191 complete, then parameter `its` contains either the iteration number at
192 which convergence was successfully reached, or failure was detected.
193 Call `PEPGetConvergedReason()` to determine if the solver converged or
194 failed and why.
195
196 Level: intermediate
197
198 .seealso: [](ch:pep), `PEPGetConvergedReason()`, `PEPSetTolerances()`
199 @*/
200 316 PetscErrorCode PEPGetIterationNumber(PEP pep,PetscInt *its)
201 {
202
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
316 PetscFunctionBegin;
203
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.
316 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
204
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.
316 PetscAssertPointer(its,2);
205 316 *its = pep->its;
206
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.
316 PetscFunctionReturn(PETSC_SUCCESS);
207 }
208
209 /*@
210 PEPGetConverged - Gets the number of converged eigenpairs.
211
212 Not Collective
213
214 Input Parameter:
215 . pep - the polynomial eigensolver context
216
217 Output Parameter:
218 . nconv - number of converged eigenpairs
219
220 Notes:
221 This function should be called after `PEPSolve()` has finished.
222
223 The value `nconv` may be different from the number of requested solutions
224 `nev`, but not larger than `ncv`, see `PEPSetDimensions()`.
225
226 Level: beginner
227
228 .seealso: [](ch:pep), `PEPSetDimensions()`, `PEPSolve()`, `PEPGetEigenpair()`
229 @*/
230 942 PetscErrorCode PEPGetConverged(PEP pep,PetscInt *nconv)
231 {
232
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
942 PetscFunctionBegin;
233
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.
942 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
234
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.
942 PetscAssertPointer(nconv,2);
235
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
942 PEPCheckSolved(pep,1);
236 942 *nconv = pep->nconv;
237
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.
942 PetscFunctionReturn(PETSC_SUCCESS);
238 }
239
240 /*@
241 PEPGetConvergedReason - Gets the reason why the `PEPSolve()` iteration was
242 stopped.
243
244 Not Collective
245
246 Input Parameter:
247 . pep - the polynomial eigensolver context
248
249 Output Parameter:
250 . reason - negative value indicates diverged, positive value converged, see
251 `PEPConvergedReason` for the possible values
252
253 Options Database Key:
254 . -pep_converged_reason - print reason for convergence/divergence, and number of iterations
255
256 Note:
257 If this routine is called before or doing the `PEPSolve()` the value of
258 `PEP_CONVERGED_ITERATING` is returned.
259
260 Level: intermediate
261
262 .seealso: [](ch:pep), `PEPSetTolerances()`, `PEPSolve()`, `PEPConvergedReason`
263 @*/
264 270 PetscErrorCode PEPGetConvergedReason(PEP pep,PEPConvergedReason *reason)
265 {
266
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
270 PetscFunctionBegin;
267
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.
270 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
268
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.
270 PetscAssertPointer(reason,2);
269
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
270 PEPCheckSolved(pep,1);
270 270 *reason = pep->reason;
271
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.
270 PetscFunctionReturn(PETSC_SUCCESS);
272 }
273
274 /*@
275 PEPGetEigenpair - Gets the `i`-th solution of the eigenproblem as computed by
276 `PEPSolve()`. The solution consists in both the eigenvalue and the eigenvector.
277
278 Collective
279
280 Input Parameters:
281 + pep - the polynomial eigensolver context
282 - i - index of the solution
283
284 Output Parameters:
285 + eigr - real part of eigenvalue
286 . eigi - imaginary part of eigenvalue
287 . Vr - real part of eigenvector
288 - Vi - imaginary part of eigenvector
289
290 Notes:
291 It is allowed to pass `NULL` for `Vr` and `Vi`, if the eigenvector is not
292 required. Otherwise, the caller must provide valid `Vec` objects, i.e.,
293 they must be created by the calling program with e.g. `MatCreateVecs()`.
294
295 If the eigenvalue is real, then `eigi` and `Vi` are set to zero. If PETSc is
296 configured with complex scalars the eigenvalue is stored
297 directly in `eigr` (`eigi` is set to zero) and the eigenvector in `Vr` (`Vi` is
298 set to zero). In any case, the user can pass `NULL` in `Vr` or `Vi` if one of
299 them is not required.
300
301 The index `i` should be a value between 0 and `nconv`-1 (see `PEPGetConverged()`).
302 Eigenpairs are indexed according to the ordering criterion established
303 with `PEPSetWhichEigenpairs()`.
304
305 The eigenvector is normalized to have unit norm.
306
307 Level: beginner
308
309 .seealso: [](ch:pep), `PEPSolve()`, `PEPGetConverged()`, `PEPSetWhichEigenpairs()`
310 @*/
311 16611 PetscErrorCode PEPGetEigenpair(PEP pep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
312 {
313 16611 PetscInt k;
314
315
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
16611 PetscFunctionBegin;
316
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.
16611 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
317
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.
16611 PetscValidLogicalCollectiveInt(pep,i,2);
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.
16611 if (Vr) { PetscValidHeaderSpecific(Vr,VEC_CLASSID,5); PetscCheckSameComm(pep,1,Vr,5); }
319
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.
16611 if (Vi) { PetscValidHeaderSpecific(Vi,VEC_CLASSID,6); PetscCheckSameComm(pep,1,Vi,6); }
320
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
16611 PEPCheckSolved(pep,1);
321
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
16611 PetscCheck(i>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
322
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
16611 PetscCheck(i<pep->nconv,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see PEPGetConverged()");
323
324
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.
16611 PetscCall(PEPComputeVectors(pep));
325 16611 k = pep->perm[i];
326
327 /* eigenvalue */
328 #if defined(PETSC_USE_COMPLEX)
329
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
8413 if (eigr) *eigr = pep->eigr[k];
330
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
8413 if (eigi) *eigi = 0;
331 #else
332
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
8198 if (eigr) *eigr = pep->eigr[k];
333
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
8198 if (eigi) *eigi = pep->eigi[k];
334 #endif
335
336 /* eigenvector */
337
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.
16611 PetscCall(BV_GetEigenvector(pep->V,k,pep->eigi[k],Vr,Vi));
338
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.
3240 PetscFunctionReturn(PETSC_SUCCESS);
339 }
340
341 /*@
342 PEPGetErrorEstimate - Returns the error estimate associated to the `i`-th
343 computed eigenpair.
344
345 Not Collective
346
347 Input Parameters:
348 + pep - the polynomial eigensolver context
349 - i - index of eigenpair
350
351 Output Parameter:
352 . errest - the error estimate
353
354 Note:
355 This is the error estimate used internally by the eigensolver. The actual
356 error bound can be computed with `PEPComputeError()`.
357
358 Level: advanced
359
360 .seealso: [](ch:pep), `PEPComputeError()`
361 @*/
362 45 PetscErrorCode PEPGetErrorEstimate(PEP pep,PetscInt i,PetscReal *errest)
363 {
364
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
45 PetscFunctionBegin;
365
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.
45 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
366
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.
45 PetscAssertPointer(errest,3);
367
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
45 PEPCheckSolved(pep,1);
368
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
45 PetscCheck(i>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
369
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
45 PetscCheck(i<pep->nconv,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see PEPGetConverged()");
370 45 *errest = pep->errest[pep->perm[i]];
371
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.
45 PetscFunctionReturn(PETSC_SUCCESS);
372 }
373
374 /*
375 PEPComputeResidualNorm_Private - Computes the norm of the residual vector
376 associated with an eigenpair.
377
378 Input Parameters:
379 kr,ki - eigenvalue
380 xr,xi - eigenvector
381 z - array of 4 work vectors (z[2],z[3] not referenced in complex scalars)
382 */
383 12003 PetscErrorCode PEPComputeResidualNorm_Private(PEP pep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm)
384 {
385 12003 Mat *A=pep->A;
386 12003 PetscInt i,nmat=pep->nmat;
387 12003 PetscScalar t[20],*vals=t,*ivals=NULL;
388 12003 Vec u,w;
389 #if !defined(PETSC_USE_COMPLEX)
390 5559 Vec ui,wi;
391 5559 PetscReal ni;
392 5559 PetscBool imag;
393 5559 PetscScalar it[20];
394 #endif
395
396
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
12003 PetscFunctionBegin;
397 12003 u = z[0]; w = z[1];
398
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.
12003 PetscCall(VecSet(u,0.0));
399 #if !defined(PETSC_USE_COMPLEX)
400 5559 ui = z[2]; wi = z[3];
401 5559 ivals = it;
402 #endif
403
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
12003 if (nmat>20) {
404 PetscCall(PetscMalloc1(nmat,&vals));
405 #if !defined(PETSC_USE_COMPLEX)
406 PetscCall(PetscMalloc1(nmat,&ivals));
407 #endif
408 }
409
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.
12003 PetscCall(PEPEvaluateBasis(pep,kr,ki,vals,ivals));
410 #if !defined(PETSC_USE_COMPLEX)
411
3/4
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✗ Branch 3 not taken.
5559 if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON))
412 imag = PETSC_FALSE;
413 else {
414 1157 imag = PETSC_TRUE;
415
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.
1157 PetscCall(VecSet(ui,0.0));
416 }
417 #endif
418
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
50382 for (i=0;i<nmat;i++) {
419
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
38379 if (vals[i]!=0.0) {
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.
38379 PetscCall(MatMult(A[i],xr,w));
421
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.
38379 PetscCall(VecAXPY(u,vals[i],w));
422 }
423 #if !defined(PETSC_USE_COMPLEX)
424
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
17002 if (imag) {
425
3/4
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✗ Branch 3 not taken.
3751 if (ivals[i]!=0 || vals[i]!=0) {
426
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.
3751 PetscCall(MatMult(A[i],xi,wi));
427
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.
3751 if (vals[i]==0) PetscCall(MatMult(A[i],xr,w));
428 }
429
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
3751 if (ivals[i]!=0) {
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.
2594 PetscCall(VecAXPY(u,-ivals[i],wi));
431
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.
2594 PetscCall(VecAXPY(ui,ivals[i],w));
432 }
433
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.
17002 if (vals[i]!=0) PetscCall(VecAXPY(ui,vals[i],wi));
434 }
435 #endif
436 }
437
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.
12003 PetscCall(VecNorm(u,NORM_2,norm));
438 #if !defined(PETSC_USE_COMPLEX)
439
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
5559 if (imag) {
440
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.
1157 PetscCall(VecNorm(ui,NORM_2,&ni));
441 1157 *norm = SlepcAbsEigenvalue(*norm,ni);
442 }
443 #endif
444
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
12003 if (nmat>20) {
445 PetscCall(PetscFree(vals));
446 #if !defined(PETSC_USE_COMPLEX)
447 PetscCall(PetscFree(ivals));
448 #endif
449 }
450
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.
2344 PetscFunctionReturn(PETSC_SUCCESS);
451 }
452
453 /*@
454 PEPComputeError - Computes the error (based on the residual norm) associated
455 with the `i`-th computed eigenpair.
456
457 Collective
458
459 Input Parameters:
460 + pep - the polynomial eigensolver context
461 . i - the solution index
462 - type - the type of error to compute, see `PEPErrorType`
463
464 Output Parameter:
465 . error - the error
466
467 Note:
468 The error can be computed in various ways, all of them based on the residual
469 norm $\|P(\lambda)x\|_2$ where $(\lambda,x)$ is the approximate eigenpair.
470
471 Level: beginner
472
473 .seealso: [](ch:pep), `PEPErrorType`, `PEPSolve()`, `PEPGetErrorEstimate()`
474 @*/
475 10738 PetscErrorCode PEPComputeError(PEP pep,PetscInt i,PEPErrorType type,PetscReal *error)
476 {
477 10738 Vec xr,xi,w[4];
478 10738 PetscScalar kr,ki;
479 10738 PetscReal t,z=0.0;
480 10738 PetscInt j;
481 10738 PetscBool flg;
482
483
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
10738 PetscFunctionBegin;
484
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.
10738 PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
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.
10738 PetscValidLogicalCollectiveInt(pep,i,2);
486
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.
10738 PetscValidLogicalCollectiveEnum(pep,type,3);
487
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.
10738 PetscAssertPointer(error,4);
488
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
10738 PEPCheckSolved(pep,1);
489
490 /* allocate work vectors */
491 #if defined(PETSC_USE_COMPLEX)
492
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.
5554 PetscCall(PEPSetWorkVecs(pep,3));
493 5554 xi = NULL;
494 5554 w[2] = NULL;
495 5554 w[3] = NULL;
496 #else
497
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.
5184 PetscCall(PEPSetWorkVecs(pep,6));
498 5184 xi = pep->work[3];
499 5184 w[2] = pep->work[4];
500 5184 w[3] = pep->work[5];
501 #endif
502 10738 xr = pep->work[0];
503 10738 w[0] = pep->work[1];
504 10738 w[1] = pep->work[2];
505
506 /* compute residual norms */
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.
10738 PetscCall(PEPGetEigenpair(pep,i,&kr,&ki,xr,xi));
508
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.
10738 PetscCall(PEPComputeResidualNorm_Private(pep,kr,ki,xr,xi,w,error));
509
510 /* compute error */
511
3/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
10738 switch (type) {
512 case PEP_ERROR_ABSOLUTE:
513 break;
514 329 case PEP_ERROR_RELATIVE:
515 329 *error /= SlepcAbsEigenvalue(kr,ki);
516 329 break;
517 10289 case PEP_ERROR_BACKWARD:
518 /* initialization of matrix norms */
519
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
10289 if (!pep->nrma[pep->nmat-1]) {
520
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
6504 for (j=0;j<pep->nmat;j++) {
521
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.
4933 PetscCall(MatHasOperation(pep->A[j],MATOP_NORM,&flg));
522
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
4933 PetscCheck(flg,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
523
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.
4933 PetscCall(MatNorm(pep->A[j],NORM_INFINITY,&pep->nrma[j]));
524 }
525 }
526 10289 t = SlepcAbsEigenvalue(kr,ki);
527
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
42576 for (j=pep->nmat-1;j>=0;j--) {
528 32287 z = z*t+pep->nrma[j];
529 }
530 10289 *error /= z;
531 10289 break;
532 default:
533 SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
534 }
535
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.
2091 PetscFunctionReturn(PETSC_SUCCESS);
536 }
537