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 | EPS routines related to the solution process | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/ | ||
15 | #include <slepc/private/bvimpl.h> | ||
16 | #include <petscdraw.h> | ||
17 | |||
18 | 85048 | PetscErrorCode EPSComputeVectors(EPS eps) | |
19 | { | ||
20 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
85048 | PetscFunctionBegin; |
21 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
85048 | EPSCheckSolved(eps,1); |
22 |
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.
|
85048 | if (eps->state==EPS_STATE_SOLVED) PetscTryTypeMethod(eps,computevectors); |
23 | 85048 | eps->state = EPS_STATE_EIGENVECTORS; | |
24 |
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.
|
85048 | PetscFunctionReturn(PETSC_SUCCESS); |
25 | } | ||
26 | |||
27 | 9545 | static PetscErrorCode EPSComputeValues(EPS eps) | |
28 | { | ||
29 | 9545 | PetscBool injective,iscomp,isfilter; | |
30 | 9545 | PetscInt i,n,aux,nconv0; | |
31 | 9545 | Mat A,B=NULL,G,Z; | |
32 | |||
33 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
9545 | PetscFunctionBegin; |
34 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9545 | switch (eps->categ) { |
35 | 7425 | case EPS_CATEGORY_KRYLOV: | |
36 | case EPS_CATEGORY_OTHER: | ||
37 |
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.
|
7425 | PetscCall(STIsInjective(eps->st,&injective)); |
38 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
7425 | if (injective) { |
39 | /* one-to-one mapping: backtransform eigenvalues */ | ||
40 |
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.
|
7209 | PetscUseTypeMethod(eps,backtransform); |
41 | } else { | ||
42 | /* compute eigenvalues from Rayleigh quotient */ | ||
43 |
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.
|
216 | PetscCall(DSGetDimensions(eps->ds,&n,NULL,NULL,NULL)); |
44 |
2/2✓ Branch 0 taken 9 times.
✓ Branch 1 taken 1 times.
|
216 | if (!n) break; |
45 |
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.
|
216 | PetscCall(EPSGetOperators(eps,&A,&B)); |
46 |
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.
|
216 | PetscCall(BVSetActiveColumns(eps->V,0,n)); |
47 |
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.
|
216 | PetscCall(DSGetCompact(eps->ds,&iscomp)); |
48 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(DSSetCompact(eps->ds,PETSC_FALSE)); |
49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(DSGetMat(eps->ds,DS_MAT_A,&G)); |
50 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(BVMatProject(eps->V,A,eps->V,G)); |
51 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&G)); |
52 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
216 | if (B) { |
53 | ✗ | PetscCall(DSGetMat(eps->ds,DS_MAT_B,&G)); | |
54 | ✗ | PetscCall(BVMatProject(eps->V,B,eps->V,G)); | |
55 | ✗ | PetscCall(DSRestoreMat(eps->ds,DS_MAT_B,&G)); | |
56 | } | ||
57 |
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.
|
216 | PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi)); |
58 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL)); |
59 |
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.
|
216 | PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi)); |
60 |
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.
|
216 | PetscCall(DSSetCompact(eps->ds,iscomp)); |
61 |
2/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
216 | if (eps->ishermitian && (!eps->isgeneralized || eps->ispositive)) { /* V = V * Z */ |
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.
|
216 | PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL)); |
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.
|
216 | PetscCall(DSGetMat(eps->ds,DS_MAT_X,&Z)); |
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.
|
216 | PetscCall(BVMultInPlace(eps->V,Z,0,n)); |
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.
|
216 | PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&Z)); |
66 | } | ||
67 | /* in case of STFILTER discard computed eigenvalues that lie outside the wanted interval */ | ||
68 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
216 | PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilter)); |
69 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
216 | if (isfilter) { |
70 | 206 | nconv0 = eps->nconv; | |
71 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1896 | for (i=0;i<eps->nconv;i++) { |
72 |
2/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
1690 | if (PetscRealPart(eps->eigr[eps->perm[i]])<eps->inta || PetscRealPart(eps->eigr[eps->perm[i]])>eps->intb) { |
73 | ✗ | eps->nconv--; | |
74 | ✗ | if (i<eps->nconv) { SlepcSwap(eps->perm[i],eps->perm[eps->nconv],aux); i--; } | |
75 | } | ||
76 | } | ||
77 |
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.
|
206 | if (nconv0>eps->nconv) PetscCall(PetscInfo(eps,"Discarded %" PetscInt_FMT " computed eigenvalues lying outside the interval\n",nconv0-eps->nconv)); |
78 | } | ||
79 | } | ||
80 | break; | ||
81 | case EPS_CATEGORY_PRECOND: | ||
82 | case EPS_CATEGORY_CONTOUR: | ||
83 | /* eigenvalues already available as an output of the solver */ | ||
84 | break; | ||
85 | } | ||
86 |
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.
|
1840 | PetscFunctionReturn(PETSC_SUCCESS); |
87 | } | ||
88 | |||
89 | /*@ | ||
90 | EPSSolve - Solves the eigensystem. | ||
91 | |||
92 | Collective | ||
93 | |||
94 | Input Parameter: | ||
95 | . eps - eigensolver context obtained from EPSCreate() | ||
96 | |||
97 | Options Database Keys: | ||
98 | + -eps_view - print information about the solver used | ||
99 | . -eps_view_mat0 - view the first matrix (A) | ||
100 | . -eps_view_mat1 - view the second matrix (B) | ||
101 | . -eps_view_vectors - view the computed eigenvectors | ||
102 | . -eps_view_values - view the computed eigenvalues | ||
103 | . -eps_converged_reason - print reason for convergence, and number of iterations | ||
104 | . -eps_error_absolute - print absolute errors of each eigenpair | ||
105 | . -eps_error_relative - print relative errors of each eigenpair | ||
106 | - -eps_error_backward - print backward errors of each eigenpair | ||
107 | |||
108 | Notes: | ||
109 | All the command-line options listed above admit an optional argument specifying | ||
110 | the viewer type and options. For instance, use '-eps_view_mat0 binary:amatrix.bin' | ||
111 | to save the A matrix to a binary file, '-eps_view_values draw' to draw the computed | ||
112 | eigenvalues graphically, or '-eps_error_relative :myerr.m:ascii_matlab' to save | ||
113 | the errors in a file that can be executed in Matlab. | ||
114 | |||
115 | Level: beginner | ||
116 | |||
117 | .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances() | ||
118 | @*/ | ||
119 | 9545 | PetscErrorCode EPSSolve(EPS eps) | |
120 | { | ||
121 | 9545 | PetscInt i; | |
122 | 9545 | PetscBool hasname; | |
123 | 9545 | STMatMode matmode; | |
124 | 9545 | Mat A,B; | |
125 | |||
126 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
9545 | PetscFunctionBegin; |
127 |
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.
|
9545 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
128 |
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.
|
9545 | if (eps->state>=EPS_STATE_SOLVED) PetscFunctionReturn(PETSC_SUCCESS); |
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.
|
9545 | PetscCall(PetscLogEventBegin(EPS_Solve,eps,0,0,0)); |
130 | |||
131 | /* Call setup */ | ||
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.
|
9545 | PetscCall(EPSSetUp(eps)); |
133 | |||
134 | /* Safeguard for matrices of size 0 */ | ||
135 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
9545 | if (eps->n == 0) { |
136 | ✗ | eps->nconv = 0; | |
137 | ✗ | eps->reason = EPS_CONVERGED_TOL; | |
138 | ✗ | eps->state = EPS_STATE_SOLVED; | |
139 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
140 | } | ||
141 | |||
142 | 9545 | eps->nconv = 0; | |
143 | 9545 | eps->its = 0; | |
144 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
302334 | for (i=0;i<eps->ncv;i++) { |
145 | 292789 | eps->eigr[i] = 0.0; | |
146 | 292789 | eps->eigi[i] = 0.0; | |
147 | 292789 | eps->errest[i] = 0.0; | |
148 | 292789 | eps->perm[i] = i; | |
149 | } | ||
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.
|
9545 | PetscCall(EPSViewFromOptions(eps,NULL,"-eps_view_pre")); |
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.
|
9545 | PetscCall(RGViewFromOptions(eps->rg,NULL,"-rg_view")); |
152 | |||
153 | /* Call solver */ | ||
154 |
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.
|
9545 | PetscUseTypeMethod(eps,solve); |
155 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
9545 | PetscCheck(eps->reason,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason"); |
156 | 9545 | eps->state = EPS_STATE_SOLVED; | |
157 | |||
158 | /* Only the first nconv columns contain useful information (except in CISS) */ | ||
159 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
9545 | PetscCall(BVSetActiveColumns(eps->V,0,eps->nconv)); |
160 |
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.
|
9545 | if (eps->twosided) PetscCall(BVSetActiveColumns(eps->W,0,eps->nconv)); |
161 | |||
162 | /* If inplace, purify eigenvectors before reverting operator */ | ||
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.
|
9545 | PetscCall(STGetMatMode(eps->st,&matmode)); |
164 |
1/10✗ 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.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
9545 | if (matmode == ST_MATMODE_INPLACE && eps->ispositive) PetscCall(EPSComputeVectors(eps)); |
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.
|
9545 | PetscCall(STPostSolve(eps->st)); |
166 | |||
167 | /* Map eigenvalues back to the original problem if appropriate */ | ||
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.
|
9545 | PetscCall(EPSComputeValues(eps)); |
169 | |||
170 | #if !defined(PETSC_USE_COMPLEX) | ||
171 | /* Reorder conjugate eigenvalues (positive imaginary first) */ | ||
172 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
55318 | for (i=0;i<eps->nconv-1;i++) { |
173 |
6/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 5 times.
✓ Branch 5 taken 5 times.
|
50350 | if (eps->eigi[i] != 0 && (eps->problem_type!=EPS_HAMILT || eps->eigr[i]!=0)) { |
174 | /* conjugate eigenvalues */ | ||
175 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
926 | if (eps->eigi[i] < 0) { |
176 | 116 | eps->eigi[i] = -eps->eigi[i]; | |
177 | 116 | eps->eigi[i+1] = -eps->eigi[i+1]; | |
178 | /* the next correction only works with eigenvectors */ | ||
179 |
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.
|
116 | PetscCall(EPSComputeVectors(eps)); |
180 |
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.
|
116 | PetscCall(BVScaleColumn(eps->V,i+1,-1.0)); |
181 |
6/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ 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.
|
116 | if (eps->W) PetscCall(BVScaleColumn(eps->W,i+1,-1.0)); |
182 | } | ||
183 | 926 | i++; | |
184 | } | ||
185 | } | ||
186 | #endif | ||
187 | |||
188 | /* Sort eigenvalues according to eps->which parameter */ | ||
189 |
6/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
✓ 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.
|
9545 | if (eps->problem_type==EPS_HAMILT) PetscCall(SlepcSortEigenvaluesSpecial(eps->sc,eps->nconv,eps->eigr,eps->eigi,eps->perm)); |
190 |
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.
|
9510 | else PetscCall(SlepcSortEigenvalues(eps->sc,eps->nconv,eps->eigr,eps->eigi,eps->perm)); |
191 |
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.
|
9545 | PetscCall(PetscLogEventEnd(EPS_Solve,eps,0,0,0)); |
192 | |||
193 | /* Various viewers */ | ||
194 |
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.
|
9545 | PetscCall(EPSViewFromOptions(eps,NULL,"-eps_view")); |
195 |
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.
|
9545 | PetscCall(EPSConvergedReasonViewFromOptions(eps)); |
196 |
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.
|
9545 | PetscCall(EPSErrorViewFromOptions(eps)); |
197 |
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.
|
9545 | PetscCall(EPSValuesViewFromOptions(eps)); |
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.
|
9545 | PetscCall(EPSVectorsViewFromOptions(eps)); |
199 | |||
200 |
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.
|
9545 | PetscCall(PetscOptionsHasName(((PetscObject)eps)->options,((PetscObject)eps)->prefix,"-eps_view_mat0",&hasname)); |
201 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
9545 | if (hasname) { |
202 | ✗ | PetscCall(EPSGetOperators(eps,&A,NULL)); | |
203 | ✗ | PetscCall(MatViewFromOptions(A,(PetscObject)eps,"-eps_view_mat0")); | |
204 | } | ||
205 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9545 | if (eps->isgeneralized) { |
206 |
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.
|
2080 | PetscCall(PetscOptionsHasName(((PetscObject)eps)->options,((PetscObject)eps)->prefix,"-eps_view_mat1",&hasname)); |
207 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
2080 | if (hasname) { |
208 | ✗ | PetscCall(EPSGetOperators(eps,NULL,&B)); | |
209 | ✗ | PetscCall(MatViewFromOptions(B,(PetscObject)eps,"-eps_view_mat1")); | |
210 | } | ||
211 | } | ||
212 | |||
213 | /* Remove deflation and initial subspaces */ | ||
214 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9545 | if (eps->nds) { |
215 |
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.
|
148 | PetscCall(BVSetNumConstraints(eps->V,0)); |
216 | 148 | eps->nds = 0; | |
217 | } | ||
218 | 9545 | eps->nini = 0; | |
219 |
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.
|
9545 | PetscFunctionReturn(PETSC_SUCCESS); |
220 | } | ||
221 | |||
222 | /*@ | ||
223 | EPSGetIterationNumber - Gets the current iteration number. If the | ||
224 | call to EPSSolve() is complete, then it returns the number of iterations | ||
225 | carried out by the solution method. | ||
226 | |||
227 | Not Collective | ||
228 | |||
229 | Input Parameter: | ||
230 | . eps - the eigensolver context | ||
231 | |||
232 | Output Parameter: | ||
233 | . its - number of iterations | ||
234 | |||
235 | Note: | ||
236 | During the i-th iteration this call returns i-1. If EPSSolve() is | ||
237 | complete, then parameter "its" contains either the iteration number at | ||
238 | which convergence was successfully reached, or failure was detected. | ||
239 | Call EPSGetConvergedReason() to determine if the solver converged or | ||
240 | failed and why. | ||
241 | |||
242 | Level: intermediate | ||
243 | |||
244 | .seealso: EPSGetConvergedReason(), EPSSetTolerances() | ||
245 | @*/ | ||
246 | 1550 | PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its) | |
247 | { | ||
248 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1550 | PetscFunctionBegin; |
249 |
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.
|
1550 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
250 |
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.
|
1550 | PetscAssertPointer(its,2); |
251 | 1550 | *its = eps->its; | |
252 |
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.
|
1550 | PetscFunctionReturn(PETSC_SUCCESS); |
253 | } | ||
254 | |||
255 | /*@ | ||
256 | EPSGetConverged - Gets the number of converged eigenpairs. | ||
257 | |||
258 | Not Collective | ||
259 | |||
260 | Input Parameter: | ||
261 | . eps - the eigensolver context | ||
262 | |||
263 | Output Parameter: | ||
264 | . nconv - number of converged eigenpairs | ||
265 | |||
266 | Note: | ||
267 | This function should be called after EPSSolve() has finished. | ||
268 | |||
269 | Level: beginner | ||
270 | |||
271 | .seealso: EPSSetDimensions(), EPSSolve(), EPSGetEigenpair() | ||
272 | @*/ | ||
273 | 5635 | PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv) | |
274 | { | ||
275 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
5635 | PetscFunctionBegin; |
276 |
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.
|
5635 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
277 |
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.
|
5635 | PetscAssertPointer(nconv,2); |
278 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
5635 | EPSCheckSolved(eps,1); |
279 |
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.
|
5635 | PetscCall(EPS_GetActualConverged(eps,nconv)); |
280 |
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.
|
1083 | PetscFunctionReturn(PETSC_SUCCESS); |
281 | } | ||
282 | |||
283 | /*@ | ||
284 | EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was | ||
285 | stopped. | ||
286 | |||
287 | Not Collective | ||
288 | |||
289 | Input Parameter: | ||
290 | . eps - the eigensolver context | ||
291 | |||
292 | Output Parameter: | ||
293 | . reason - negative value indicates diverged, positive value converged | ||
294 | |||
295 | Options Database Key: | ||
296 | . -eps_converged_reason - print the reason to a viewer | ||
297 | |||
298 | Notes: | ||
299 | Possible values for reason are | ||
300 | + EPS_CONVERGED_TOL - converged up to tolerance | ||
301 | . EPS_CONVERGED_USER - converged due to a user-defined condition | ||
302 | . EPS_DIVERGED_ITS - required more than max_it iterations to reach convergence | ||
303 | . EPS_DIVERGED_BREAKDOWN - generic breakdown in method | ||
304 | - EPS_DIVERGED_SYMMETRY_LOST - pseudo-Lanczos was not able to keep symmetry | ||
305 | |||
306 | Can only be called after the call to EPSSolve() is complete. | ||
307 | |||
308 | Level: intermediate | ||
309 | |||
310 | .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason | ||
311 | @*/ | ||
312 | 1509 | PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason) | |
313 | { | ||
314 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1509 | 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.
|
1509 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
316 |
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.
|
1509 | PetscAssertPointer(reason,2); |
317 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
1509 | EPSCheckSolved(eps,1); |
318 | 1509 | *reason = eps->reason; | |
319 |
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.
|
1509 | PetscFunctionReturn(PETSC_SUCCESS); |
320 | } | ||
321 | |||
322 | /*@ | ||
323 | EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant | ||
324 | subspace. | ||
325 | |||
326 | Collective | ||
327 | |||
328 | Input Parameter: | ||
329 | . eps - the eigensolver context | ||
330 | |||
331 | Output Parameter: | ||
332 | . v - an array of vectors | ||
333 | |||
334 | Notes: | ||
335 | This function should be called after EPSSolve() has finished. | ||
336 | |||
337 | The user should provide in v an array of nconv vectors, where nconv is | ||
338 | the value returned by EPSGetConverged(). | ||
339 | |||
340 | The first k vectors returned in v span an invariant subspace associated | ||
341 | with the first k computed eigenvalues (note that this is not true if the | ||
342 | k-th eigenvalue is complex and matrix A is real; in this case the first | ||
343 | k+1 vectors should be used). An invariant subspace X of A satisfies Ax | ||
344 | in X for all x in X (a similar definition applies for generalized | ||
345 | eigenproblems). | ||
346 | |||
347 | Level: intermediate | ||
348 | |||
349 | .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve() | ||
350 | @*/ | ||
351 | 80 | PetscErrorCode EPSGetInvariantSubspace(EPS eps,Vec v[]) | |
352 | { | ||
353 | 80 | PetscInt i; | |
354 | 80 | BV V=eps->V; | |
355 | 80 | Vec w; | |
356 | |||
357 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
80 | PetscFunctionBegin; |
358 |
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.
|
80 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
359 |
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.
|
80 | PetscAssertPointer(v,2); |
360 |
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.
|
80 | PetscValidHeaderSpecific(*v,VEC_CLASSID,2); |
361 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
80 | EPSCheckSolved(eps,1); |
362 |
2/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
80 | PetscCheck(eps->ishermitian || eps->state!=EPS_STATE_EIGENVECTORS,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSGetInvariantSubspace must be called before EPSGetEigenpair,EPSGetEigenvector or EPSComputeError"); |
363 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
80 | if (eps->balance!=EPS_BALANCE_NONE && eps->D) { |
364 |
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(BVDuplicateResize(eps->V,eps->nconv,&V)); |
365 |
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(BVSetActiveColumns(eps->V,0,eps->nconv)); |
366 |
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(BVCopy(eps->V,V)); |
367 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
55 | for (i=0;i<eps->nconv;i++) { |
368 |
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.
|
45 | PetscCall(BVGetColumn(V,i,&w)); |
369 |
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.
|
45 | PetscCall(VecPointwiseDivide(w,w,eps->D)); |
370 |
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.
|
45 | PetscCall(BVRestoreColumn(V,i,&w)); |
371 | } | ||
372 |
3/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
|
10 | PetscCall(BVOrthogonalize(V,NULL)); |
373 | } | ||
374 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
625 | for (i=0;i<eps->nconv;i++) PetscCall(BVCopyVec(V,i,v[i])); |
375 |
7/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ 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.
|
80 | if (eps->balance!=EPS_BALANCE_NONE && eps->D) PetscCall(BVDestroy(&V)); |
376 |
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.
|
16 | PetscFunctionReturn(PETSC_SUCCESS); |
377 | } | ||
378 | |||
379 | /*@ | ||
380 | EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by | ||
381 | EPSSolve(). The solution consists in both the eigenvalue and the eigenvector. | ||
382 | |||
383 | Collective | ||
384 | |||
385 | Input Parameters: | ||
386 | + eps - eigensolver context | ||
387 | - i - index of the solution | ||
388 | |||
389 | Output Parameters: | ||
390 | + eigr - real part of eigenvalue | ||
391 | . eigi - imaginary part of eigenvalue | ||
392 | . Vr - real part of eigenvector | ||
393 | - Vi - imaginary part of eigenvector | ||
394 | |||
395 | Notes: | ||
396 | It is allowed to pass NULL for Vr and Vi, if the eigenvector is not | ||
397 | required. Otherwise, the caller must provide valid Vec objects, i.e., | ||
398 | they must be created by the calling program with e.g. MatCreateVecs(). | ||
399 | |||
400 | If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is | ||
401 | configured with complex scalars the eigenvalue is stored | ||
402 | directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is | ||
403 | set to zero). In both cases, the user can pass NULL in eigi and Vi. | ||
404 | |||
405 | The index i should be a value between 0 and nconv-1 (see EPSGetConverged()). | ||
406 | Eigenpairs are indexed according to the ordering criterion established | ||
407 | with EPSSetWhichEigenpairs(). | ||
408 | |||
409 | The 2-norm of the eigenvector is one unless the problem is generalized | ||
410 | Hermitian. In this case the eigenvector is normalized with respect to the | ||
411 | norm defined by the B matrix. | ||
412 | |||
413 | Level: beginner | ||
414 | |||
415 | .seealso: EPSGetEigenvalue(), EPSGetEigenvector(), EPSGetLeftEigenvector(), EPSSolve(), | ||
416 | EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace() | ||
417 | @*/ | ||
418 | 55066 | PetscErrorCode EPSGetEigenpair(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi) | |
419 | { | ||
420 | 55066 | PetscInt nconv; | |
421 | |||
422 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
55066 | PetscFunctionBegin; |
423 |
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.
|
55066 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
424 |
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.
|
55066 | PetscValidLogicalCollectiveInt(eps,i,2); |
425 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
55066 | EPSCheckSolved(eps,1); |
426 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
55066 | PetscCheck(i>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative"); |
427 |
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.
|
55066 | PetscCall(EPS_GetActualConverged(eps,&nconv)); |
428 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
55066 | PetscCheck(i<nconv,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see EPSGetConverged()"); |
429 |
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.
|
55066 | PetscCall(EPSGetEigenvalue(eps,i,eigr,eigi)); |
430 |
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.
|
55066 | if (Vr || Vi) PetscCall(EPSGetEigenvector(eps,i,Vr,Vi)); |
431 |
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.
|
11090 | PetscFunctionReturn(PETSC_SUCCESS); |
432 | } | ||
433 | |||
434 | /*@ | ||
435 | EPSGetEigenvalue - Gets the i-th eigenvalue as computed by EPSSolve(). | ||
436 | |||
437 | Not Collective | ||
438 | |||
439 | Input Parameters: | ||
440 | + eps - eigensolver context | ||
441 | - i - index of the solution | ||
442 | |||
443 | Output Parameters: | ||
444 | + eigr - real part of eigenvalue | ||
445 | - eigi - imaginary part of eigenvalue | ||
446 | |||
447 | Notes: | ||
448 | If the eigenvalue is real, then eigi is set to zero. If PETSc is | ||
449 | configured with complex scalars the eigenvalue is stored | ||
450 | directly in eigr (eigi is set to zero). | ||
451 | |||
452 | The index i should be a value between 0 and nconv-1 (see EPSGetConverged()). | ||
453 | Eigenpairs are indexed according to the ordering criterion established | ||
454 | with EPSSetWhichEigenpairs(). | ||
455 | |||
456 | Level: beginner | ||
457 | |||
458 | .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair() | ||
459 | @*/ | ||
460 | 107963 | PetscErrorCode EPSGetEigenvalue(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi) | |
461 | { | ||
462 | 107963 | PetscInt k,nconv; | |
463 | #if !defined(PETSC_USE_COMPLEX) | ||
464 | 55729 | PetscInt k2, iquad; | |
465 | #endif | ||
466 | |||
467 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
107963 | PetscFunctionBegin; |
468 |
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.
|
107963 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
469 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
107963 | EPSCheckSolved(eps,1); |
470 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
107963 | PetscCheck(i>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative"); |
471 |
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.
|
107963 | PetscCall(EPS_GetActualConverged(eps,&nconv)); |
472 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
107963 | PetscCheck(i<nconv,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see EPSGetConverged()"); |
473 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
107963 | if (nconv==eps->nconv) { |
474 | 101771 | k = eps->perm[i]; | |
475 | #if defined(PETSC_USE_COMPLEX) | ||
476 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
49274 | if (eigr) *eigr = eps->eigr[k]; |
477 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
49274 | if (eigi) *eigi = 0; |
478 | #else | ||
479 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
52497 | if (eigr) *eigr = eps->eigr[k]; |
480 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
52497 | if (eigi) *eigi = eps->eigi[k]; |
481 | #endif | ||
482 | } else { | ||
483 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6192 | PetscCheck(eps->problem_type==EPS_BSE || eps->problem_type==EPS_HAMILT,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Problem type should be BSE or Hamiltonian"); |
484 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
|
6192 | if (eps->problem_type==EPS_BSE) { |
485 | /* BSE problem, even index is +lambda, odd index is -lambda */ | ||
486 | 5672 | k = eps->perm[i/2]; | |
487 | #if defined(PETSC_USE_COMPLEX) | ||
488 |
3/4✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
2960 | if (eigr) *eigr = (i%2)? -eps->eigr[k]: eps->eigr[k]; |
489 |
1/2✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
|
2960 | if (eigi) *eigi = 0; |
490 | #else | ||
491 |
3/4✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
2712 | if (eigr) *eigr = (i%2)? -eps->eigr[k]: eps->eigr[k]; |
492 |
1/2✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
|
2712 | if (eigi) *eigi = eps->eigi[k]; |
493 | #endif | ||
494 |
1/2✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
|
520 | } else if (eps->problem_type==EPS_HAMILT) { |
495 | /* Hamiltonian eigenproblem */ | ||
496 | 520 | k = eps->perm[i/2]; | |
497 | #if defined(PETSC_USE_COMPLEX) | ||
498 | ✗ | if (eigr) *eigr = (i%2)? -eps->eigr[k]: eps->eigr[k]; | |
499 | ✗ | if (eigi) *eigi = 0; | |
500 | #else | ||
501 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
520 | if (eps->eigi[k]==0.0) { /* real eigenvalue */ |
502 | ✗ | if (eigr) *eigr = (i%2)? -eps->eigr[k]: eps->eigr[k]; | |
503 | ✗ | if (eigi) *eigi = 0.0; | |
504 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
520 | } else if (eps->eigr[k]==0.0) { /* purely imaginary eigenvalue */ |
505 |
1/2✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
|
120 | if (eigr) *eigr = 0.0; |
506 |
3/4✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
120 | if (eigi) *eigi = (i%2)? -eps->eigi[k]: eps->eigi[k]; |
507 | } else { /* quadruple eigenvalue (-conj(lambda),-lambda,lambda,conj(lambda)) */ | ||
508 | 400 | iquad = i%2; /* index within the 4 values */ | |
509 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
400 | if (i>1) { |
510 | 320 | k2 = eps->perm[(i-2)/2]; | |
511 |
3/4✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✗ Branch 3 not taken.
|
320 | if (eps->eigr[k]==eps->eigr[k2] && eps->eigi[k]==-eps->eigi[k2]) iquad += 2; |
512 | } | ||
513 |
3/4✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
400 | if (eigr) *eigr = (iquad<2)? -eps->eigr[k]: eps->eigr[k]; |
514 |
3/4✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
400 | if (eigi) *eigi = (iquad%3)? -eps->eigi[k]: eps->eigi[k]; |
515 | } | ||
516 | #endif | ||
517 | } | ||
518 | } | ||
519 |
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.
|
22063 | PetscFunctionReturn(PETSC_SUCCESS); |
520 | } | ||
521 | |||
522 | /*@ | ||
523 | EPSGetEigenvector - Gets the i-th right eigenvector as computed by EPSSolve(). | ||
524 | |||
525 | Collective | ||
526 | |||
527 | Input Parameters: | ||
528 | + eps - eigensolver context | ||
529 | - i - index of the solution | ||
530 | |||
531 | Output Parameters: | ||
532 | + Vr - real part of eigenvector | ||
533 | - Vi - imaginary part of eigenvector | ||
534 | |||
535 | Notes: | ||
536 | The caller must provide valid Vec objects, i.e., they must be created | ||
537 | by the calling program with e.g. MatCreateVecs(). | ||
538 | |||
539 | If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is | ||
540 | configured with complex scalars the eigenvector is stored | ||
541 | directly in Vr (Vi is set to zero). In any case, the user can pass NULL in Vr | ||
542 | or Vi if one of them is not required. | ||
543 | |||
544 | The index i should be a value between 0 and nconv-1 (see EPSGetConverged()). | ||
545 | Eigenpairs are indexed according to the ordering criterion established | ||
546 | with EPSSetWhichEigenpairs(). | ||
547 | |||
548 | The 2-norm of the eigenvector is one unless the problem is generalized | ||
549 | Hermitian. In this case the eigenvector is normalized with respect to the | ||
550 | norm defined by the B matrix. | ||
551 | |||
552 | Level: beginner | ||
553 | |||
554 | .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair(), EPSGetLeftEigenvector() | ||
555 | @*/ | ||
556 | 79915 | PetscErrorCode EPSGetEigenvector(EPS eps,PetscInt i,Vec Vr,Vec Vi) | |
557 | { | ||
558 | 79915 | PetscInt nconv; | |
559 | |||
560 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
79915 | PetscFunctionBegin; |
561 |
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.
|
79915 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
562 |
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.
|
79915 | PetscValidLogicalCollectiveInt(eps,i,2); |
563 |
16/46✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ 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.
|
79915 | if (Vr) { PetscValidHeaderSpecific(Vr,VEC_CLASSID,3); PetscCheckSameComm(eps,1,Vr,3); } |
564 |
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.
|
79915 | if (Vi) { PetscValidHeaderSpecific(Vi,VEC_CLASSID,4); PetscCheckSameComm(eps,1,Vi,4); } |
565 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
79915 | EPSCheckSolved(eps,1); |
566 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
79915 | PetscCheck(i>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative"); |
567 |
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.
|
79915 | PetscCall(EPS_GetActualConverged(eps,&nconv)); |
568 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
79915 | PetscCheck(i<nconv,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see EPSGetConverged()"); |
569 |
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.
|
79915 | PetscCall(EPSComputeVectors(eps)); |
570 |
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.
|
79915 | PetscCall(EPS_GetEigenvector(eps,eps->V,i,Vr,Vi)); |
571 |
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.
|
14325 | PetscFunctionReturn(PETSC_SUCCESS); |
572 | } | ||
573 | |||
574 | /*@ | ||
575 | EPSGetLeftEigenvector - Gets the i-th left eigenvector as computed by EPSSolve(). | ||
576 | |||
577 | Collective | ||
578 | |||
579 | Input Parameters: | ||
580 | + eps - eigensolver context | ||
581 | - i - index of the solution | ||
582 | |||
583 | Output Parameters: | ||
584 | + Wr - real part of left eigenvector | ||
585 | - Wi - imaginary part of left eigenvector | ||
586 | |||
587 | Notes: | ||
588 | The caller must provide valid Vec objects, i.e., they must be created | ||
589 | by the calling program with e.g. MatCreateVecs(). | ||
590 | |||
591 | If the corresponding eigenvalue is real, then Wi is set to zero. If PETSc is | ||
592 | configured with complex scalars the eigenvector is stored directly in Wr | ||
593 | (Wi is set to zero). In any case, the user can pass NULL in Wr or Wi if | ||
594 | one of them is not required. | ||
595 | |||
596 | The index i should be a value between 0 and nconv-1 (see EPSGetConverged()). | ||
597 | Eigensolutions are indexed according to the ordering criterion established | ||
598 | with EPSSetWhichEigenpairs(). | ||
599 | |||
600 | Left eigenvectors are available only if the twosided flag was set, see | ||
601 | EPSSetTwoSided(). | ||
602 | |||
603 | Level: intermediate | ||
604 | |||
605 | .seealso: EPSGetEigenvector(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSSetTwoSided() | ||
606 | @*/ | ||
607 | 4505 | PetscErrorCode EPSGetLeftEigenvector(EPS eps,PetscInt i,Vec Wr,Vec Wi) | |
608 | { | ||
609 | 4505 | PetscInt nconv; | |
610 | 4505 | PetscBool trivial; | |
611 | 4505 | Mat H; | |
612 | 4505 | IS is[2]; | |
613 | 4505 | Vec v; | |
614 | |||
615 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4505 | PetscFunctionBegin; |
616 |
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.
|
4505 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
617 |
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.
|
4505 | PetscValidLogicalCollectiveInt(eps,i,2); |
618 |
16/46✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ 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.
|
4505 | if (Wr) { PetscValidHeaderSpecific(Wr,VEC_CLASSID,3); PetscCheckSameComm(eps,1,Wr,3); } |
619 |
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.
|
4505 | if (Wi) { PetscValidHeaderSpecific(Wi,VEC_CLASSID,4); PetscCheckSameComm(eps,1,Wi,4); } |
620 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
4505 | EPSCheckSolved(eps,1); |
621 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
4505 | PetscCheck(i>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative"); |
622 |
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.
|
4505 | PetscCall(EPS_GetActualConverged(eps,&nconv)); |
623 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
4505 | PetscCheck(i<nconv,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see EPSGetConverged()"); |
624 | |||
625 | 4505 | trivial = (eps->problem_type==EPS_HEP || eps->problem_type==EPS_GHEP || eps->problem_type==EPS_BSE)? PETSC_TRUE: PETSC_FALSE; | |
626 |
3/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
4505 | if (!trivial) PetscCheck(eps->twosided,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must request left vectors with EPSSetTwoSided"); |
627 | |||
628 |
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.
|
4505 | PetscCall(EPSComputeVectors(eps)); |
629 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
4505 | if (trivial) { |
630 |
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.
|
3612 | PetscCall(EPS_GetEigenvector(eps,eps->V,i,Wr,Wi)); |
631 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
3612 | if (eps->problem_type==EPS_BSE) { /* change sign of bottom part of the vector */ |
632 |
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.
|
3612 | PetscCall(STGetMatrix(eps->st,0,&H)); |
633 |
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.
|
3612 | PetscCall(MatNestGetISs(H,is,NULL)); |
634 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
3612 | if (Wr) { |
635 |
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.
|
3612 | PetscCall(VecGetSubVector(Wr,is[1],&v)); |
636 |
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.
|
3612 | PetscCall(VecScale(v,-1.0)); |
637 |
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.
|
3612 | PetscCall(VecRestoreSubVector(Wr,is[1],&v)); |
638 | } | ||
639 | #if !defined(PETSC_USE_COMPLEX) | ||
640 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
1768 | if (Wi) { |
641 | ✗ | PetscCall(VecGetSubVector(Wi,is[1],&v)); | |
642 | ✗ | PetscCall(VecScale(v,-1.0)); | |
643 | ✗ | PetscCall(VecRestoreSubVector(Wi,is[1],&v)); | |
644 | } | ||
645 | #endif | ||
646 | } | ||
647 | } else { | ||
648 |
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.
|
893 | PetscCall(EPS_GetEigenvector(eps,eps->W,i,Wr,Wi)); |
649 | } | ||
650 |
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.
|
817 | PetscFunctionReturn(PETSC_SUCCESS); |
651 | } | ||
652 | |||
653 | /*@ | ||
654 | EPSGetErrorEstimate - Returns the error estimate associated to the i-th | ||
655 | computed eigenpair. | ||
656 | |||
657 | Not Collective | ||
658 | |||
659 | Input Parameters: | ||
660 | + eps - eigensolver context | ||
661 | - i - index of eigenpair | ||
662 | |||
663 | Output Parameter: | ||
664 | . errest - the error estimate | ||
665 | |||
666 | Notes: | ||
667 | This is the error estimate used internally by the eigensolver. The actual | ||
668 | error bound can be computed with EPSComputeError(). See also the users | ||
669 | manual for details. | ||
670 | |||
671 | Level: advanced | ||
672 | |||
673 | .seealso: EPSComputeError() | ||
674 | @*/ | ||
675 | 612 | PetscErrorCode EPSGetErrorEstimate(EPS eps,PetscInt i,PetscReal *errest) | |
676 | { | ||
677 | 612 | PetscInt nconv; | |
678 | |||
679 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
612 | PetscFunctionBegin; |
680 |
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.
|
612 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
681 |
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.
|
612 | PetscAssertPointer(errest,3); |
682 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
612 | EPSCheckSolved(eps,1); |
683 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
612 | PetscCheck(i>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative"); |
684 |
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.
|
612 | PetscCall(EPS_GetActualConverged(eps,&nconv)); |
685 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
612 | PetscCheck(i<nconv,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see EPSGetConverged()"); |
686 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
612 | if (nconv==eps->nconv) { |
687 | 612 | *errest = eps->errest[eps->perm[i]]; | |
688 | } else { | ||
689 | ✗ | PetscCheck(eps->problem_type==EPS_BSE,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Problem type should be BSE"); | |
690 | /* BSE problem, even index is +lambda, odd index is -lambda, assume both have same error */ | ||
691 | ✗ | *errest = eps->errest[eps->perm[i/2]]; | |
692 | } | ||
693 |
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.
|
116 | PetscFunctionReturn(PETSC_SUCCESS); |
694 | } | ||
695 | |||
696 | /* | ||
697 | EPSComputeResidualNorm_Private - Computes the norm of the residual vector | ||
698 | associated with an eigenpair. | ||
699 | |||
700 | Input Parameters: | ||
701 | trans - whether A' must be used instead of A | ||
702 | kr,ki - eigenvalue | ||
703 | xr,xi - eigenvector | ||
704 | z - three work vectors (the second one not referenced in complex scalars) | ||
705 | */ | ||
706 | 44419 | PetscErrorCode EPSComputeResidualNorm_Private(EPS eps,PetscBool trans,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm) | |
707 | { | ||
708 | 44419 | PetscInt nmat; | |
709 | 44419 | Mat A,B; | |
710 | 44419 | Vec u,w; | |
711 | 44419 | PetscScalar alpha; | |
712 | #if !defined(PETSC_USE_COMPLEX) | ||
713 | 22558 | Vec v; | |
714 | 22558 | PetscReal ni,nr; | |
715 | #endif | ||
716 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
44419 | PetscErrorCode (*matmult)(Mat,Vec,Vec) = trans? MatMultHermitianTranspose: MatMult; |
717 | |||
718 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
44419 | PetscFunctionBegin; |
719 | 44419 | u = z[0]; w = z[2]; | |
720 |
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.
|
44419 | PetscCall(STGetNumMatrices(eps->st,&nmat)); |
721 |
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.
|
44419 | PetscCall(STGetMatrix(eps->st,0,&A)); |
722 |
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.
|
44419 | if (nmat>1) PetscCall(STGetMatrix(eps->st,1,&B)); |
723 | |||
724 | #if !defined(PETSC_USE_COMPLEX) | ||
725 | 22558 | v = z[1]; | |
726 |
3/4✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 5 times.
|
22558 | if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) { |
727 | #endif | ||
728 |
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.
|
43247 | PetscCall((*matmult)(A,xr,u)); /* u=A*x */ |
729 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
43247 | if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) { |
730 |
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.
|
43247 | if (nmat>1) PetscCall((*matmult)(B,xr,w)); |
731 |
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.
|
21478 | else PetscCall(VecCopy(xr,w)); /* w=B*x */ |
732 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
43247 | alpha = trans? -PetscConj(kr): -kr; |
733 |
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.
|
43247 | PetscCall(VecAXPY(u,alpha,w)); /* u=A*x-k*B*x */ |
734 | } | ||
735 |
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.
|
43247 | PetscCall(VecNorm(u,NORM_2,norm)); |
736 | #if !defined(PETSC_USE_COMPLEX) | ||
737 | } else { | ||
738 |
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.
|
1172 | PetscCall((*matmult)(A,xr,u)); /* u=A*xr */ |
739 |
1/2✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
|
1172 | if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) { |
740 |
6/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ 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.
|
1172 | if (nmat>1) PetscCall((*matmult)(B,xr,v)); |
741 |
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.
|
1156 | else PetscCall(VecCopy(xr,v)); /* v=B*xr */ |
742 |
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.
|
1172 | PetscCall(VecAXPY(u,-kr,v)); /* u=A*xr-kr*B*xr */ |
743 |
6/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ 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.
|
1172 | if (nmat>1) PetscCall((*matmult)(B,xi,w)); |
744 |
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.
|
1156 | else PetscCall(VecCopy(xi,w)); /* w=B*xi */ |
745 |
7/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
1172 | PetscCall(VecAXPY(u,trans?-ki:ki,w)); /* u=A*xr-kr*B*xr+ki*B*xi */ |
746 | } | ||
747 |
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.
|
1172 | PetscCall(VecNorm(u,NORM_2,&nr)); |
748 |
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.
|
1172 | PetscCall((*matmult)(A,xi,u)); /* u=A*xi */ |
749 |
1/2✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
|
1172 | if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) { |
750 |
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.
|
1172 | PetscCall(VecAXPY(u,-kr,w)); /* u=A*xi-kr*B*xi */ |
751 |
7/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
1172 | PetscCall(VecAXPY(u,trans?ki:-ki,v)); /* u=A*xi-kr*B*xi-ki*B*xr */ |
752 | } | ||
753 |
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.
|
1172 | PetscCall(VecNorm(u,NORM_2,&ni)); |
754 | 1172 | *norm = SlepcAbsEigenvalue(nr,ni); | |
755 | } | ||
756 | #endif | ||
757 |
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.
|
8816 | PetscFunctionReturn(PETSC_SUCCESS); |
758 | } | ||
759 | |||
760 | /*@ | ||
761 | EPSComputeError - Computes the error (based on the residual norm) associated | ||
762 | with the i-th computed eigenpair. | ||
763 | |||
764 | Collective | ||
765 | |||
766 | Input Parameters: | ||
767 | + eps - the eigensolver context | ||
768 | . i - the solution index | ||
769 | - type - the type of error to compute | ||
770 | |||
771 | Output Parameter: | ||
772 | . error - the error | ||
773 | |||
774 | Notes: | ||
775 | The error can be computed in various ways, all of them based on the residual | ||
776 | norm ||Ax-kBx||_2 where k is the eigenvalue and x is the eigenvector. | ||
777 | |||
778 | Level: beginner | ||
779 | |||
780 | .seealso: EPSErrorType, EPSSolve(), EPSGetErrorEstimate() | ||
781 | @*/ | ||
782 | 41034 | PetscErrorCode EPSComputeError(EPS eps,PetscInt i,EPSErrorType type,PetscReal *error) | |
783 | { | ||
784 | 41034 | Mat A,B; | |
785 | 41034 | Vec xr,xi,w[3]; | |
786 | 41034 | PetscReal t,vecnorm=1.0,errorl; | |
787 | 41034 | PetscScalar kr,ki; | |
788 | 41034 | PetscBool flg; | |
789 | |||
790 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
41034 | PetscFunctionBegin; |
791 |
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.
|
41034 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
792 |
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.
|
41034 | PetscValidLogicalCollectiveInt(eps,i,2); |
793 |
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.
|
41034 | PetscValidLogicalCollectiveEnum(eps,type,3); |
794 |
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.
|
41034 | PetscAssertPointer(error,4); |
795 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
41034 | EPSCheckSolved(eps,1); |
796 | |||
797 | /* allocate work vectors */ | ||
798 | #if defined(PETSC_USE_COMPLEX) | ||
799 |
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.
|
20106 | PetscCall(EPSSetWorkVecs(eps,3)); |
800 | 20106 | xi = NULL; | |
801 | 20106 | w[1] = NULL; | |
802 | #else | ||
803 |
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.
|
20928 | PetscCall(EPSSetWorkVecs(eps,5)); |
804 | 20928 | xi = eps->work[3]; | |
805 | 20928 | w[1] = eps->work[4]; | |
806 | #endif | ||
807 | 41034 | xr = eps->work[0]; | |
808 | 41034 | w[0] = eps->work[1]; | |
809 | 41034 | w[2] = eps->work[2]; | |
810 | |||
811 | /* compute residual norm */ | ||
812 |
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.
|
41034 | PetscCall(EPSGetEigenpair(eps,i,&kr,&ki,xr,xi)); |
813 |
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.
|
41034 | PetscCall(EPSComputeResidualNorm_Private(eps,PETSC_FALSE,kr,ki,xr,xi,w,error)); |
814 | |||
815 | /* compute 2-norm of eigenvector */ | ||
816 |
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.
|
41034 | if (eps->problem_type==EPS_GHEP) PetscCall(VecNorm(xr,NORM_2,&vecnorm)); |
817 | |||
818 | /* if two-sided, compute left residual norm and take the maximum */ | ||
819 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
41034 | if (eps->twosided) { |
820 |
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.
|
450 | PetscCall(EPSGetLeftEigenvector(eps,i,xr,xi)); |
821 |
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.
|
450 | PetscCall(EPSComputeResidualNorm_Private(eps,PETSC_TRUE,kr,ki,xr,xi,w,&errorl)); |
822 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
673 | *error = PetscMax(*error,errorl); |
823 | } | ||
824 | |||
825 | /* compute error */ | ||
826 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
41034 | switch (type) { |
827 | case EPS_ERROR_ABSOLUTE: | ||
828 | break; | ||
829 | 23584 | case EPS_ERROR_RELATIVE: | |
830 | 23584 | *error /= SlepcAbsEigenvalue(kr,ki)*vecnorm; | |
831 | 23584 | break; | |
832 | 17198 | case EPS_ERROR_BACKWARD: | |
833 | /* initialization of matrix norms */ | ||
834 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
17198 | if (!eps->nrma) { |
835 |
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.
|
395 | PetscCall(STGetMatrix(eps->st,0,&A)); |
836 |
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.
|
395 | PetscCall(MatHasOperation(A,MATOP_NORM,&flg)); |
837 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
395 | PetscCheck(flg,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation"); |
838 |
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.
|
395 | PetscCall(MatNorm(A,NORM_INFINITY,&eps->nrma)); |
839 | } | ||
840 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
17198 | if (eps->isgeneralized) { |
841 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
16945 | if (!eps->nrmb) { |
842 |
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.
|
336 | PetscCall(STGetMatrix(eps->st,1,&B)); |
843 |
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.
|
336 | PetscCall(MatHasOperation(B,MATOP_NORM,&flg)); |
844 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
336 | PetscCheck(flg,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation"); |
845 |
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.
|
336 | PetscCall(MatNorm(B,NORM_INFINITY,&eps->nrmb)); |
846 | } | ||
847 | 253 | } else eps->nrmb = 1.0; | |
848 | 17198 | t = SlepcAbsEigenvalue(kr,ki); | |
849 | 17198 | *error /= (eps->nrma+t*eps->nrmb)*vecnorm; | |
850 | 17198 | break; | |
851 | ✗ | default: | |
852 | ✗ | SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type"); | |
853 | } | ||
854 |
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.
|
8100 | PetscFunctionReturn(PETSC_SUCCESS); |
855 | } | ||
856 | |||
857 | /* | ||
858 | EPSGetStartVector - Generate a suitable vector to be used as the starting vector | ||
859 | for the recurrence that builds the right subspace. | ||
860 | |||
861 | Collective | ||
862 | |||
863 | Input Parameters: | ||
864 | + eps - the eigensolver context | ||
865 | - i - iteration number | ||
866 | |||
867 | Output Parameters: | ||
868 | . breakdown - flag indicating that a breakdown has occurred | ||
869 | |||
870 | Notes: | ||
871 | The start vector is computed from another vector: for the first step (i=0), | ||
872 | the first initial vector is used (see EPSSetInitialSpace()); otherwise a random | ||
873 | vector is created. Then this vector is forced to be in the range of OP (only | ||
874 | for generalized definite problems) and orthonormalized with respect to all | ||
875 | V-vectors up to i-1. The resulting vector is placed in V[i]. | ||
876 | |||
877 | The flag breakdown is set to true if either i=0 and the vector belongs to the | ||
878 | deflation space, or i>0 and the vector is linearly dependent with respect | ||
879 | to the V-vectors. | ||
880 | */ | ||
881 | 7002 | PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,PetscBool *breakdown) | |
882 | { | ||
883 | 7002 | PetscReal norm; | |
884 | 7002 | PetscBool lindep; | |
885 | 7002 | Vec w,z; | |
886 | |||
887 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
7002 | PetscFunctionBegin; |
888 |
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.
|
7002 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
889 |
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.
|
7002 | PetscValidLogicalCollectiveInt(eps,i,2); |
890 | |||
891 | /* For the first step, use the first initial vector, otherwise a random one */ | ||
892 |
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.
|
7002 | if (i>0 || eps->nini==0) PetscCall(BVSetRandomColumn(eps->V,i)); |
893 | |||
894 | /* Force the vector to be in the range of OP for generalized problems with B-inner product */ | ||
895 |
6/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
|
7002 | if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) { |
896 |
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.
|
1033 | PetscCall(BVCreateVec(eps->V,&w)); |
897 |
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.
|
1033 | PetscCall(BVCopyVec(eps->V,i,w)); |
898 |
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.
|
1033 | PetscCall(BVGetColumn(eps->V,i,&z)); |
899 |
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.
|
1033 | PetscCall(STApply(eps->st,w,z)); |
900 |
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.
|
1033 | PetscCall(BVRestoreColumn(eps->V,i,&z)); |
901 |
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.
|
1033 | PetscCall(VecDestroy(&w)); |
902 | } | ||
903 | |||
904 | /* Orthonormalize the vector with respect to previous vectors */ | ||
905 |
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.
|
7002 | PetscCall(BVOrthogonalizeColumn(eps->V,i,NULL,&norm,&lindep)); |
906 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
7002 | if (breakdown) *breakdown = lindep; |
907 |
2/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
6645 | else if (lindep || norm == 0.0) { |
908 | ✗ | PetscCheck(i,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Initial vector is zero or belongs to the deflation space"); | |
909 | ✗ | PetscCheck(!i,PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Unable to generate more start vectors"); | |
910 | } | ||
911 |
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.
|
7002 | PetscCall(BVScaleColumn(eps->V,i,1.0/norm)); |
912 |
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.
|
1346 | PetscFunctionReturn(PETSC_SUCCESS); |
913 | } | ||
914 | |||
915 | /* | ||
916 | EPSGetLeftStartVector - Generate a suitable vector to be used as the left starting | ||
917 | vector for the recurrence that builds the left subspace. See EPSGetStartVector(). | ||
918 | */ | ||
919 | 254 | PetscErrorCode EPSGetLeftStartVector(EPS eps,PetscInt i,PetscBool *breakdown) | |
920 | { | ||
921 | 254 | PetscReal norm; | |
922 | 254 | PetscBool lindep; | |
923 | 254 | Vec w,z; | |
924 | |||
925 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
254 | PetscFunctionBegin; |
926 |
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.
|
254 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
927 |
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.
|
254 | PetscValidLogicalCollectiveInt(eps,i,2); |
928 | |||
929 | /* For the first step, use the first initial vector, otherwise a random one */ | ||
930 |
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.
|
254 | if (i>0 || eps->ninil==0) PetscCall(BVSetRandomColumn(eps->W,i)); |
931 | |||
932 | /* Force the vector to be in the range of OP' for generalized problems with B-inner product */ | ||
933 |
5/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 8 times.
|
254 | if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) { |
934 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13 | PetscCall(BVCreateVec(eps->W,&w)); |
935 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13 | PetscCall(BVCopyVec(eps->W,i,w)); |
936 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13 | PetscCall(BVGetColumn(eps->W,i,&z)); |
937 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13 | PetscCall(STApplyHermitianTranspose(eps->st,w,z)); |
938 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13 | PetscCall(BVRestoreColumn(eps->W,i,&z)); |
939 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13 | PetscCall(VecDestroy(&w)); |
940 | } | ||
941 | |||
942 | /* Orthonormalize the vector with respect to previous vectors */ | ||
943 |
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.
|
254 | PetscCall(BVOrthogonalizeColumn(eps->W,i,NULL,&norm,&lindep)); |
944 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
254 | if (breakdown) *breakdown = lindep; |
945 |
2/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
199 | else if (lindep || norm == 0.0) { |
946 | ✗ | PetscCheck(i,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Left initial vector is zero"); | |
947 | ✗ | PetscCheck(!i,PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Unable to generate more left start vectors"); | |
948 | } | ||
949 |
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.
|
254 | PetscCall(BVScaleColumn(eps->W,i,1.0/norm)); |
950 |
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.
|
50 | PetscFunctionReturn(PETSC_SUCCESS); |
951 | } | ||
952 |