GCC Code Coverage Report


Directory: ./
File: src/eps/impls/power/power.c
Date: 2025-10-04 04:19:13
Exec Total Coverage
Lines: 634 681 93.1%
Functions: 37 40 92.5%
Branches: 1756 3314 53.0%

Line Branch Exec Source
1 /*
2 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 SLEPc - Scalable Library for Eigenvalue Problem Computations
4 Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5
6 This file is part of SLEPc.
7 SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 */
10 /*
11 SLEPc eigensolver: "power"
12
13 Method: Power Iteration
14
15 Algorithm:
16
17 This solver implements the power iteration for finding dominant
18 eigenpairs. It also includes the following well-known methods:
19 - Inverse Iteration: when used in combination with shift-and-invert
20 spectral transformation.
21 - Rayleigh Quotient Iteration (RQI): also with shift-and-invert plus
22 a variable shift.
23
24 It can also be used for nonlinear inverse iteration on the problem
25 A(x)*x=lambda*B(x)*x, where A and B are not constant but depend on x.
26
27 References:
28
29 [1] "Single Vector Iteration Methods in SLEPc", SLEPc Technical Report
30 STR-2, available at https://slepc.upv.es.
31 */
32
33 #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
34 #include <slepcblaslapack.h>
35 /* PETSc headers */
36 #include <petscdm.h>
37 #include <petscsnes.h>
38
39 static PetscErrorCode EPSPowerFormFunction_Update(SNES,Vec,Vec,void*);
40 static PetscErrorCode EPSSolve_Power(EPS);
41 static PetscErrorCode EPSSolve_TS_Power(EPS);
42
43 typedef struct {
44 EPSPowerShiftType shift_type;
45 PetscBool nonlinear;
46 PetscBool update;
47 PetscBool sign_normalization;
48 SNES snes;
49 PetscErrorCode (*formFunctionB)(SNES,Vec,Vec,void*);
50 void *formFunctionBctx;
51 PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
52 void *formFunctionActx;
53 PetscErrorCode (*formFunctionAB)(SNES,Vec,Vec,Vec,void*);
54 PetscErrorCode (*formNorm)(SNES,Vec,PetscReal*,void*);
55 void *formNormCtx;
56 PetscInt idx; /* index of the first nonzero entry in the iteration vector */
57 PetscMPIInt p; /* process id of the owner of idx */
58 PetscReal norm0; /* norm of initial vector */
59 } EPS_POWER;
60
61 930 static PetscErrorCode SNESMonitor_PowerUpdate(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
62 {
63 930 EPS eps;
64
65
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
930 PetscFunctionBegin;
66
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
930 PetscCall(PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps));
67
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
930 PetscCheck(eps,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_NULL,"No composed EPS");
68 /* Call EPS monitor on each SNES iteration */
69
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
930 PetscCall(EPSMonitor(eps,its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev)));
70
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.
186 PetscFunctionReturn(PETSC_SUCCESS);
71 }
72
73 506 static PetscErrorCode EPSSetUp_Power(EPS eps)
74 {
75 506 EPS_POWER *power = (EPS_POWER*)eps->data;
76 506 STMatMode mode;
77 506 Mat A,B,P;
78 506 Vec res;
79 506 PetscContainer container;
80 506 PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
81 506 PetscErrorCode (*formJacobianA)(SNES,Vec,Mat,Mat,void*);
82 506 PetscErrorCode (*formNorm)(SNES,Vec,PetscReal*,void*);
83 506 void *ctx;
84
85
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
506 PetscFunctionBegin;
86
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
506 EPSCheckNotStructured(eps);
87
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
506 if (eps->nev==0) eps->nev = 1;
88
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
506 if (eps->ncv!=PETSC_DETERMINE) {
89
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
50 PetscCheck(eps->ncv>=eps->nev,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev");
90 456 } else eps->ncv = eps->nev;
91
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.
506 if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
92
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
506 if (eps->max_it==PETSC_DETERMINE) {
93 /* SNES will directly return the solution for us, and we need to do only one iteration */
94
4/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
306 if (power->nonlinear && power->update) eps->max_it = 1;
95 176 else eps->max_it = PetscMax(1000*eps->nev,100*eps->n);
96 }
97
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.
506 if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
98
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
506 PetscCheck(eps->which==EPS_LARGEST_MAGNITUDE || eps->which==EPS_TARGET_MAGNITUDE,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only largest magnitude or target magnitude eigenvalues");
99
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
506 if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
100
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
50 PetscCheck(!power->nonlinear,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts not allowed in nonlinear problems");
101
6/10
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
50 EPSCheckSinvertCayleyCondition(eps,PETSC_TRUE," (with variable shifts)");
102
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
50 PetscCall(STGetMatMode(eps->st,&mode));
103
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
50 PetscCheck(mode!=ST_MATMODE_INPLACE,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
104 }
105
10/26
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 8 times.
✓ Branch 10 taken 2 times.
✓ Branch 11 taken 8 times.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 10 times.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 10 times.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
506 EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_CONVERGENCE | EPS_FEATURE_THRESHOLD);
106
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.
506 EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION);
107
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(EPSAllocateSolution(eps,0));
108
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
506 PetscCall(EPS_SetInnerProduct(eps));
109
110
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
506 if (power->nonlinear) {
111
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
330 PetscCheck(eps->nev==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration cannot compute more than one eigenvalue");
112
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(EPSSetWorkVecs(eps,3));
113
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.
330 PetscCheck(!power->update || eps->max_it==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"More than one iteration is not allowed for Newton eigensolver (SNES)");
114
115 /* Set up SNES solver */
116 /* If SNES was setup, we need to reset it so that it will be consistent with the current EPS */
117
5/8
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
330 if (power->snes) PetscCall(SNESReset(power->snes));
118 else PetscCall(EPSPowerGetSNES(eps,&power->snes));
119
120 /* For nonlinear solver (Newton), we should scale the initial vector back.
121 The initial vector will be scaled in EPSSetUp. */
122
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.
330 if (eps->IS) PetscCall(VecNorm(eps->IS[0],NORM_2,&power->norm0));
123
124
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(EPSGetOperators(eps,&A,&B));
125
126 /* form function */
127
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(PetscObjectQueryFunction((PetscObject)A,"formFunction",&formFunctionA));
128
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
330 PetscCheck(formFunctionA,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER,"For nonlinear inverse iteration you must compose a callback function 'formFunction' in matrix A");
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.
330 PetscCall(PetscObjectQuery((PetscObject)A,"formFunctionCtx",(PetscObject*)&container));
130
5/8
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
330 if (container) PetscCall(PetscContainerGetPointer(container,&ctx));
131 else ctx = NULL;
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.
330 PetscCall(MatCreateVecs(A,&res,NULL));
133 330 power->formFunctionA = formFunctionA;
134 330 power->formFunctionActx = ctx;
135
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
330 if (power->update) {
136
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(SNESSetFunction(power->snes,res,EPSPowerFormFunction_Update,ctx));
137
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(PetscObjectQueryFunction((PetscObject)A,"formFunctionAB",&power->formFunctionAB));
138
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(SNESMonitorSet(power->snes,SNESMonitor_PowerUpdate,NULL,NULL));
139 }
140
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 else PetscCall(SNESSetFunction(power->snes,res,formFunctionA,ctx));
141
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(VecDestroy(&res));
142
143 /* form Jacobian */
144
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(PetscObjectQueryFunction((PetscObject)A,"formJacobian",&formJacobianA));
145
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
330 PetscCheck(formJacobianA,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER,"For nonlinear inverse iteration you must compose a callback function 'formJacobian' in matrix A");
146
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(PetscObjectQuery((PetscObject)A,"formJacobianCtx",(PetscObject*)&container));
147
5/8
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
330 if (container) PetscCall(PetscContainerGetPointer(container,&ctx));
148 else ctx = NULL;
149 /* This allows users to compute a different preconditioning matrix than A.
150 * It is useful when A and B are shell matrices.
151 */
152
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(STGetPreconditionerMat(eps->st,&P));
153 /* If users do not provide a matrix, we simply use A */
154
7/8
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
330 PetscCall(SNESSetJacobian(power->snes,A,P? P:A,formJacobianA,ctx));
155
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(SNESSetFromOptions(power->snes));
156
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(SNESSetUp(power->snes));
157
158 330 ctx = NULL;
159 330 formNorm = NULL;
160
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
330 if (B) {
161
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(PetscObjectQueryFunction((PetscObject)B,"formFunction",&power->formFunctionB));
162
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(PetscObjectQuery((PetscObject)B,"formFunctionCtx",(PetscObject*)&container));
163
6/10
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ 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.
330 if (power->formFunctionB && container) PetscCall(PetscContainerGetPointer(container,&power->formFunctionBctx));
164 else power->formFunctionBctx = NULL;
165
166 /* form norm */
167
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(PetscObjectQueryFunction((PetscObject)B,"formNorm",&formNorm));
168
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
330 if (formNorm) {
169
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
90 PetscCall(PetscObjectQuery((PetscObject)B,"formNormCtx",(PetscObject*)&container));
170
5/8
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
90 if (container) PetscCall(PetscContainerGetPointer(container,&ctx));
171 }
172 }
173 330 power->formNorm = formNorm;
174 330 power->formNormCtx = ctx;
175 } else {
176
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.
176 if (eps->twosided) PetscCall(EPSSetWorkVecs(eps,3));
177
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 else PetscCall(EPSSetWorkVecs(eps,2));
178
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
176 PetscCall(DSSetType(eps->ds,DSNHEP));
179
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
176 PetscCall(DSAllocate(eps->ds,eps->nev));
180 }
181 /* dispatch solve method */
182
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
506 if (eps->twosided) {
183
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
40 PetscCheck(!power->nonlinear,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration does not have two-sided variant");
184
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
40 PetscCheck(power->shift_type!=EPS_POWER_SHIFT_WILKINSON,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Two-sided variant does not support Wilkinson shifts");
185 40 eps->ops->solve = EPSSolve_TS_Power;
186 466 } else eps->ops->solve = EPSSolve_Power;
187
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.
102 PetscFunctionReturn(PETSC_SUCCESS);
188 }
189
190 /*
191 Find the index of the first nonzero entry of x, and the process that owns it.
192 */
193 330 static PetscErrorCode FirstNonzeroIdx(Vec x,PetscInt *idx,PetscMPIInt *p)
194 {
195 330 PetscInt i,first,last,N;
196 330 PetscLayout map;
197 330 const PetscScalar *xx;
198
199
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
330 PetscFunctionBegin;
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.
330 PetscCall(VecGetSize(x,&N));
201
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(VecGetOwnershipRange(x,&first,&last));
202
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(VecGetArrayRead(x,&xx));
203
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
3630 for (i=first;i<last;i++) {
204
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
3630 if (PetscAbsScalar(xx[i-first])>10*PETSC_MACHINE_EPSILON) break;
205 }
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.
330 PetscCall(VecRestoreArrayRead(x,&xx));
207
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
330 if (i==last) i=N;
208
28/58
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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 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.
330 PetscCallMPI(MPIU_Allreduce(&i,idx,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)x)));
209
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
330 PetscCheck(*idx!=N,PetscObjectComm((PetscObject)x),PETSC_ERR_PLIB,"Zero vector found");
210
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(VecGetLayout(x,&map));
211
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(PetscLayoutFindOwner(map,*idx,p));
212
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.
66 PetscFunctionReturn(PETSC_SUCCESS);
213 }
214
215 /*
216 Normalize a vector x with respect to a given norm as well as, optionally, the
217 sign of the first nonzero entry (assumed to be idx in process p).
218 */
219 271919 static PetscErrorCode Normalize(Vec x,PetscReal norm,PetscInt idx,PetscMPIInt p,PetscBool sign_normalization,PetscScalar *sign)
220 {
221 271919 PetscScalar alpha=1.0;
222 271919 PetscInt first,last;
223 271919 PetscReal absv;
224 271919 const PetscScalar *xx;
225
226
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
271919 PetscFunctionBegin;
227
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
271919 if (sign_normalization) {
228
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
269449 PetscCall(VecGetOwnershipRange(x,&first,&last));
229
2/4
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
269449 if (idx>=first && idx<last) {
230
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
269449 PetscCall(VecGetArrayRead(x,&xx));
231 269449 absv = PetscAbsScalar(xx[idx-first]);
232
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
269449 if (absv>10*PETSC_MACHINE_EPSILON) alpha = xx[idx-first]/absv;
233
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
269449 PetscCall(VecRestoreArrayRead(x,&xx));
234 }
235
15/30
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ 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 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
538898 PetscCallMPI(MPI_Bcast(&alpha,1,MPIU_SCALAR,p,PetscObjectComm((PetscObject)x)));
236 }
237
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
271919 if (sign) *sign = alpha;
238 271919 alpha *= norm;
239
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
271919 PetscCall(VecScale(x,1.0/alpha));
240
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.
54501 PetscFunctionReturn(PETSC_SUCCESS);
241 }
242
243 5470 static PetscErrorCode EPSPowerUpdateFunctionB(EPS eps,Vec x,Vec Bx)
244 {
245 5470 EPS_POWER *power = (EPS_POWER*)eps->data;
246 5470 Mat B;
247
248
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
5470 PetscFunctionBegin;
249
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5470 PetscCall(STResetMatrixState(eps->st));
250
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5470 PetscCall(EPSGetOperators(eps,NULL,&B));
251
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
5470 if (B) {
252
5/8
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
5470 if (power->formFunctionB) PetscCall((*power->formFunctionB)(power->snes,x,Bx,power->formFunctionBctx));
253 else PetscCall(MatMult(B,x,Bx));
254 } else PetscCall(VecCopy(x,Bx));
255
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.
1094 PetscFunctionReturn(PETSC_SUCCESS);
256 }
257
258 3335 static PetscErrorCode EPSPowerUpdateFunctionA(EPS eps,Vec x,Vec Ax)
259 {
260 3335 EPS_POWER *power = (EPS_POWER*)eps->data;
261 3335 Mat A;
262
263
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
3335 PetscFunctionBegin;
264
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3335 PetscCall(STResetMatrixState(eps->st));
265
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3335 PetscCall(EPSGetOperators(eps,&A,NULL));
266
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
3335 PetscCheck(A,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_NULL,"Matrix A is required for an eigenvalue problem");
267
5/8
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
3335 if (power->formFunctionA) PetscCall((*power->formFunctionA)(power->snes,x,Ax,power->formFunctionActx));
268 else PetscCall(MatMult(A,x,Ax));
269
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.
667 PetscFunctionReturn(PETSC_SUCCESS);
270 }
271
272 5630 static PetscErrorCode EPSPowerFormFunction_Update(SNES snes,Vec x,Vec y,void *ctx)
273 {
274 5630 EPS eps;
275 5630 PetscReal bx;
276 5630 Vec Bx;
277 5630 PetscScalar sign;
278 5630 EPS_POWER *power;
279
280
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
5630 PetscFunctionBegin;
281
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5630 PetscCall(PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps));
282
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
5630 PetscCheck(eps,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_NULL,"No composed EPS");
283
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5630 PetscCall(STResetMatrixState(eps->st));
284 5630 power = (EPS_POWER*)eps->data;
285 5630 Bx = eps->work[2];
286
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.
5630 if (power->formFunctionAB) PetscCall((*power->formFunctionAB)(snes,x,y,Bx,ctx));
287 else {
288
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3335 PetscCall(EPSPowerUpdateFunctionA(eps,x,y));
289
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3335 PetscCall(EPSPowerUpdateFunctionB(eps,x,Bx));
290 }
291
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.
5630 if (power->formNorm) PetscCall((*power->formNorm)(snes,Bx,&bx,power->formNormCtx));
292
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3920 else PetscCall(VecNorm(Bx,NORM_2,&bx));
293
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5630 PetscCall(Normalize(Bx,bx,power->idx,power->p,power->sign_normalization,&sign));
294
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5630 PetscCall(VecAXPY(y,-1.0,Bx));
295 /* Keep tracking eigenvalue update. It would be useful when we want to monitor solver progress via snes monitor. */
296
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
5630 eps->eigr[(eps->nconv < eps->nev)? eps->nconv:(eps->nconv-1)] = 1.0/(bx*sign);
297
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.
5630 PetscFunctionReturn(PETSC_SUCCESS);
298 }
299
300 /*
301 Use SNES to compute y = A^{-1}*B*x for the nonlinear problem
302 */
303 1955 static PetscErrorCode EPSPowerApply_SNES(EPS eps,Vec x,Vec y)
304 {
305 1955 EPS_POWER *power = (EPS_POWER*)eps->data;
306 1955 Vec Bx;
307
308
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1955 PetscFunctionBegin;
309
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1955 PetscCall(VecCopy(x,y));
310
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.
1955 if (power->update) PetscCall(SNESSolve(power->snes,NULL,y));
311 else {
312 1805 Bx = eps->work[2];
313
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1805 PetscCall(SNESSolve(power->snes,Bx,y));
314 }
315
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.
391 PetscFunctionReturn(PETSC_SUCCESS);
316 }
317
318 /*
319 Use nonlinear inverse power to compute an initial guess
320 */
321 130 static PetscErrorCode EPSPowerComputeInitialGuess_Update(EPS eps)
322 {
323 130 EPS powereps;
324 130 Mat A,B,P;
325 130 Vec v1,v2;
326 130 SNES snes;
327 130 DM dm,newdm;
328 130 PetscBool sign_normalization;
329
330
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
130 PetscFunctionBegin;
331
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSCreate(PetscObjectComm((PetscObject)eps),&powereps));
332
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSGetOperators(eps,&A,&B));
333
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSSetType(powereps,EPSPOWER));
334
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSSetOperators(powereps,A,B));
335
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSSetTolerances(powereps,1e-6,4));
336
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSSetOptionsPrefix(powereps,((PetscObject)eps)->prefix));
337
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSAppendOptionsPrefix(powereps,"init_"));
338
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSSetProblemType(powereps,EPS_GNHEP));
339
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSSetWhichEigenpairs(powereps,EPS_TARGET_MAGNITUDE));
340
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSPowerSetNonlinear(powereps,PETSC_TRUE));
341
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(STGetPreconditionerMat(eps->st,&P));
342 /* attach dm to initial solve */
343
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSPowerGetSNES(eps,&snes));
344
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(SNESGetDM(snes,&dm));
345 /* use dmshell to temporarily store snes context */
346
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(DMCreate(PetscObjectComm((PetscObject)eps),&newdm));
347
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(DMSetType(newdm,DMSHELL));
348
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(DMSetUp(newdm));
349
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(DMCopyDMSNES(dm,newdm));
350
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSPowerGetSNES(powereps,&snes));
351
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(SNESSetDM(snes,dm));
352
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSPowerGetSignNormalization(eps,&sign_normalization));
353
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSPowerSetSignNormalization(powereps,sign_normalization));
354
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSSetFromOptions(powereps));
355
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.
130 if (P) PetscCall(STSetPreconditionerMat(powereps->st,P));
356
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSSolve(powereps));
357
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(BVGetColumn(eps->V,0,&v2));
358
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(BVGetColumn(powereps->V,0,&v1));
359
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(VecCopy(v1,v2));
360
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(BVRestoreColumn(powereps->V,0,&v1));
361
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(BVRestoreColumn(eps->V,0,&v2));
362
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
130 PetscCall(EPSDestroy(&powereps));
363 /* restore context back to the old nonlinear solver */
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.
130 PetscCall(DMCopyDMSNES(newdm,dm));
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.
130 PetscCall(DMDestroy(&newdm));
366
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.
26 PetscFunctionReturn(PETSC_SUCCESS);
367 }
368
369 466 static PetscErrorCode EPSSolve_Power(EPS eps)
370 {
371 466 EPS_POWER *power = (EPS_POWER*)eps->data;
372 466 PetscInt k,ld;
373 466 Vec v,y,e,Bx;
374 466 Mat A;
375 466 KSP ksp;
376 466 PetscReal relerr,norm,norm1,rt1,rt2,cs1;
377 466 PetscScalar theta,rho,delta,sigma,alpha2,beta1,sn1,*T,sign;
378 466 PetscBool breakdown;
379 466 KSPConvergedReason reason;
380 466 SNESConvergedReason snesreason;
381
382
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
466 PetscFunctionBegin;
383 466 e = eps->work[0];
384 466 y = eps->work[1];
385
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
466 if (power->nonlinear) Bx = eps->work[2];
386 else Bx = NULL;
387
388
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.
466 if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) PetscCall(STGetKSP(eps->st,&ksp));
389
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
466 if (power->nonlinear) {
390
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(PetscObjectCompose((PetscObject)power->snes,"eps",(PetscObject)eps));
391 /* Compute an initial guess only when users do not provide one */
392
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.
330 if (power->update && !eps->nini) PetscCall(EPSPowerComputeInitialGuess_Update(eps));
393
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 } else PetscCall(DSGetLeadingDimension(eps->ds,&ld));
394
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.
466 if (!power->update) PetscCall(EPSGetStartVector(eps,0,NULL));
395
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
466 if (power->nonlinear) {
396
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(BVGetColumn(eps->V,0,&v));
397
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
330 if (eps->nini) {
398 /* We scale the initial vector back if the initial vector was provided by users */
399
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
20 PetscCall(VecScale(v,power->norm0));
400 }
401 /* Do a couple things:
402 * 1) Determine the first non-zero index for Bx and the proc that owns it. This will be
403 * used if performing normalization by the sign of the first nonzero element of Bx.
404 * 2) Scale Bx by its norm. For non-update power iterations, Bx (in code terms) is used
405 * as the RHS argument to SNESSolve. And recall that the formula for generalized
406 * inverse power iterations in this case is: (Ax)_n = (Bx)_{n-1}/|(Bx)_{n-1}| (in
407 * math terms)
408 */
409
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(EPSPowerUpdateFunctionB(eps,v,Bx));
410
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(BVRestoreColumn(eps->V,0,&v));
411
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.
330 if (power->formNorm) PetscCall((*power->formNorm)(power->snes,Bx,&norm,power->formNormCtx));
412
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
240 else PetscCall(VecNorm(Bx,NORM_2,&norm));
413
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(FirstNonzeroIdx(Bx,&power->idx,&power->p));
414
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
330 PetscCall(Normalize(Bx,norm,power->idx,power->p,power->sign_normalization,NULL));
415 }
416
417
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(STGetShift(eps->st,&sigma)); /* original shift */
418 466 rho = sigma;
419
420
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
266575 while (eps->reason == EPS_CONVERGED_ITERATING) {
421 266109 eps->its++;
422 266109 k = eps->nconv;
423
424 /* y = OP v */
425
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
266109 PetscCall(BVGetColumn(eps->V,k,&v));
426
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.
266109 if (power->nonlinear) PetscCall(EPSPowerApply_SNES(eps,v,y));
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.
264154 else PetscCall(STApply(eps->st,v,y));
428
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
266109 PetscCall(BVRestoreColumn(eps->V,k,&v));
429
430 /* purge previously converged eigenvectors */
431
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
266109 if (PetscUnlikely(power->nonlinear)) {
432 /* We do not need to call this for Newton eigensolver since eigenvalue is
433 * updated in function evaluations.
434 */
435
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1955 if (!power->update) {
436
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1805 PetscCall(EPSPowerUpdateFunctionB(eps,y,Bx));
437
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.
1805 if (power->formNorm) PetscCall((*power->formNorm)(power->snes,Bx,&norm,power->formNormCtx));
438
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1125 else PetscCall(VecNorm(Bx,NORM_2,&norm));
439
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1805 PetscCall(Normalize(Bx,norm,power->idx,power->p,power->sign_normalization,&sign));
440 }
441 } else {
442
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
264154 PetscCall(DSGetArray(eps->ds,DS_MAT_A,&T));
443
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
264154 PetscCall(BVSetActiveColumns(eps->V,0,k));
444
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
264154 PetscCall(BVOrthogonalizeVec(eps->V,y,T+k*ld,&norm,NULL));
445 }
446
447 /* theta = (v,y)_B */
448
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
266109 PetscCall(BVSetActiveColumns(eps->V,k,k+1));
449
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
266109 PetscCall(BVDotVec(eps->V,y,&theta));
450
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
266109 if (!power->nonlinear) {
451 264154 T[k+k*ld] = theta;
452
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
264154 PetscCall(DSRestoreArray(eps->ds,DS_MAT_A,&T));
453 }
454
455 /* Eigenvalue is already stored in function evaluations.
456 * Assign eigenvalue to theta to make the rest of the code consistent
457 */
458
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
266109 if (power->update) theta = eps->eigr[eps->nconv];
459
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
265959 else if (power->nonlinear) theta = 1.0/(norm*sign); /* Eigenvalue: 1/|Bx| */
460
461
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
266109 if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
462
463 /* approximate eigenvalue is the Rayleigh quotient */
464 265764 eps->eigr[eps->nconv] = theta;
465
466 /**
467 * If the Newton method (update, SNES) is used, we do not compute "relerr"
468 * since SNES determines the convergence.
469 */
470
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
265764 if (PetscUnlikely(power->update)) relerr = 0.;
471 else {
472 /* compute relative error as ||y-theta v||_2/|theta| */
473
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
265614 PetscCall(VecCopy(y,e));
474
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
265614 PetscCall(BVGetColumn(eps->V,k,&v));
475
7/8
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
265614 PetscCall(VecAXPY(e,power->nonlinear?-1.0:-theta,v));
476
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
265614 PetscCall(BVRestoreColumn(eps->V,k,&v));
477
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
265614 PetscCall(VecNorm(e,NORM_2,&relerr));
478
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
265614 if (PetscUnlikely(power->nonlinear)) relerr *= PetscAbsScalar(theta);
479 263809 else relerr /= PetscAbsScalar(theta);
480 }
481
482 } else { /* RQI */
483
484 /* delta = ||y||_B */
485 345 delta = norm;
486
487 /* compute relative error */
488
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
345 if (rho == 0.0) relerr = PETSC_MAX_REAL;
489 345 else relerr = 1.0 / (norm*PetscAbsScalar(rho));
490
491 /* approximate eigenvalue is the shift */
492 345 eps->eigr[eps->nconv] = rho;
493
494 /* compute new shift */
495
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
345 if (relerr<eps->tol) {
496 80 rho = sigma; /* if converged, restore original shift */
497
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
80 PetscCall(STSetShift(eps->st,rho));
498 } else {
499 265 rho = rho + PetscConj(theta)/(delta*delta); /* Rayleigh quotient R(v) */
500
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
265 if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) {
501 /* beta1 is the norm of the residual associated with R(v) */
502
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(BVGetColumn(eps->V,k,&v));
503
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(VecAXPY(v,-PetscConj(theta)/(delta*delta),y));
504
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(BVRestoreColumn(eps->V,k,&v));
505
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(BVScaleColumn(eps->V,k,1.0/delta));
506
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(BVNormColumn(eps->V,k,NORM_2,&norm1));
507 95 beta1 = norm1;
508
509 /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
510
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(STGetMatrix(eps->st,0,&A));
511
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(BVGetColumn(eps->V,k,&v));
512
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(MatMult(A,v,e));
513
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(VecDot(v,e,&alpha2));
514
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(BVRestoreColumn(eps->V,k,&v));
515 95 alpha2 = alpha2 / (beta1 * beta1);
516
517 /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
518
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
519
10/20
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
95 PetscCallBLAS("LAPACKlaev2",LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1));
520
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
95 PetscCall(PetscFPTrapPop());
521
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
95 if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
522 40 else rho = rt2;
523 }
524 /* update operator according to new shift */
525
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
265 PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_FALSE));
526
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
265 PetscCall(STSetShift(eps->st,rho));
527
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
265 PetscCall(KSPGetConvergedReason(ksp,&reason));
528
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 10 times.
265 if (reason) {
529
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
2 PetscCall(PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n"));
530 2 rho *= 1+10*PETSC_MACHINE_EPSILON;
531
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
2 PetscCall(STSetShift(eps->st,rho));
532
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
2 PetscCall(KSPGetConvergedReason(ksp,&reason));
533
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
2 PetscCheck(!reason,PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Second factorization failed");
534 }
535
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
265 PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_TRUE));
536 }
537 }
538 266109 eps->errest[eps->nconv] = relerr;
539
540 /* normalize */
541
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.
266109 if (!power->nonlinear) PetscCall(Normalize(y,norm,power->idx,power->p,power->sign_normalization,NULL));
542
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
266109 PetscCall(BVInsertVec(eps->V,k,y));
543
544
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
266109 if (PetscUnlikely(power->update)) {
545
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(SNESGetConvergedReason(power->snes,&snesreason));
546 /* For Newton eigensolver, we are ready to return once SNES converged. */
547
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
150 if (snesreason>0) eps->nconv = 1;
548
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
265959 } else if (PetscUnlikely(relerr<eps->tol)) { /* accept eigenpair */
549 488 eps->nconv = eps->nconv + 1;
550
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
488 if (eps->nconv<eps->nev) {
551
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
302 PetscCall(EPSGetStartVector(eps,eps->nconv,&breakdown));
552
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
302 if (breakdown) {
553 eps->reason = EPS_DIVERGED_BREAKDOWN;
554 PetscCall(PetscInfo(eps,"Unable to generate more start vectors\n"));
555 break;
556 }
557 }
558 }
559 /* For Newton eigensolver, monitor will be called from SNES monitor */
560
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.
266109 if (!power->update) PetscCall(EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev)));
561
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
266109 PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx));
562
563 /**
564 * When a customized stopping test is used, and reason can be set to be converged (EPS_CONVERGED_USER).
565 * In this case, we need to increase eps->nconv to "1" so users can retrieve the solution.
566 */
567
4/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
266109 if (PetscUnlikely(power->nonlinear && eps->reason>0)) eps->nconv = 1;
568 }
569
570
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.
466 if (power->nonlinear) PetscCall(PetscObjectCompose((PetscObject)power->snes,"eps",NULL));
571 else {
572
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSSetDimensions(eps->ds,eps->nconv,0,0));
573
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
574 }
575
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.
94 PetscFunctionReturn(PETSC_SUCCESS);
576 }
577
578 40 static PetscErrorCode EPSSolve_TS_Power(EPS eps)
579 {
580 40 EPS_POWER *power = (EPS_POWER*)eps->data;
581 40 PetscInt k,ld;
582 40 Vec v,w,y,e,z;
583 40 KSP ksp;
584 40 PetscReal relerr=1.0,relerrl,delta;
585 40 PetscScalar theta,rho,alpha,sigma;
586 40 PetscBool breakdown,breakdownl;
587 40 KSPConvergedReason reason;
588
589
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
40 PetscFunctionBegin;
590 40 e = eps->work[0];
591 40 v = eps->work[1];
592 40 w = eps->work[2];
593
594
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.
40 if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) PetscCall(STGetKSP(eps->st,&ksp));
595
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(DSGetLeadingDimension(eps->ds,&ld));
596
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(EPSGetStartVector(eps,0,NULL));
597
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(EPSGetLeftStartVector(eps,0,NULL));
598
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(BVBiorthonormalizeColumn(eps->V,eps->W,0,NULL));
599
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(BVCopyVec(eps->V,0,v));
600
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(BVCopyVec(eps->W,0,w));
601
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(STGetShift(eps->st,&sigma)); /* original shift */
602 40 rho = sigma;
603
604 40 while (eps->reason == EPS_CONVERGED_ITERATING) {
605 4157 eps->its++;
606 4157 k = eps->nconv;
607
608 /* y = OP v, z = OP' w */
609
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVGetColumn(eps->V,k,&y));
610
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(STApply(eps->st,v,y));
611
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVRestoreColumn(eps->V,k,&y));
612
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVGetColumn(eps->W,k,&z));
613
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(STApplyHermitianTranspose(eps->st,w,z));
614
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVRestoreColumn(eps->W,k,&z));
615
616 /* purge previously converged eigenvectors */
617
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVBiorthogonalizeColumn(eps->V,eps->W,k));
618
619 /* theta = (w,y)_B */
620
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVSetActiveColumns(eps->V,k,k+1));
621
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVDotVec(eps->V,w,&theta));
622 4157 theta = PetscConj(theta);
623
624
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
4157 if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
625
626 /* approximate eigenvalue is the Rayleigh quotient */
627 3982 eps->eigr[eps->nconv] = theta;
628
629 /* compute relative errors as ||y-theta v||_2/|theta| and ||z-conj(theta) w||_2/|theta|*/
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.
3982 PetscCall(BVCopyVec(eps->V,k,e));
631
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3982 PetscCall(VecAXPY(e,-theta,v));
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.
3982 PetscCall(VecNorm(e,NORM_2,&relerr));
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.
3982 PetscCall(BVCopyVec(eps->W,k,e));
634
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3982 PetscCall(VecAXPY(e,-PetscConj(theta),w));
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.
3982 PetscCall(VecNorm(e,NORM_2,&relerrl));
636
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
7697 relerr = PetscMax(relerr,relerrl)/PetscAbsScalar(theta);
637 }
638
639 /* normalize */
640
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVSetActiveColumns(eps->V,k,k+1));
641
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVGetColumn(eps->W,k,&z));
642
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVDotVec(eps->V,z,&alpha));
643
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVRestoreColumn(eps->W,k,&z));
644 4157 delta = PetscSqrtReal(PetscAbsScalar(alpha));
645
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
4157 PetscCheck(delta!=0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Breakdown in two-sided Power/RQI");
646
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVScaleColumn(eps->V,k,1.0/PetscConj(alpha/delta)));
647
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVScaleColumn(eps->W,k,1.0/delta));
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.
4157 PetscCall(BVCopyVec(eps->V,k,v));
649
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(BVCopyVec(eps->W,k,w));
650
651
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
4157 if (power->shift_type == EPS_POWER_SHIFT_RAYLEIGH) { /* RQI */
652
653 /* compute relative error */
654
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
175 if (rho == 0.0) relerr = PETSC_MAX_REAL;
655 175 else relerr = 1.0 / PetscAbsScalar(delta*rho);
656
657 /* approximate eigenvalue is the shift */
658 175 eps->eigr[eps->nconv] = rho;
659
660 /* compute new shift */
661
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
175 if (relerr<eps->tol) {
662 35 rho = sigma; /* if converged, restore original shift */
663
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
35 PetscCall(STSetShift(eps->st,rho));
664 } else {
665 140 rho = rho + PetscConj(theta)/(delta*delta); /* Rayleigh quotient R(v) */
666 /* update operator according to new shift */
667
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
140 PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_FALSE));
668
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
140 PetscCall(STSetShift(eps->st,rho));
669
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
140 PetscCall(KSPGetConvergedReason(ksp,&reason));
670
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
140 if (reason) {
671
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
5 PetscCall(PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n"));
672 5 rho *= 1+10*PETSC_MACHINE_EPSILON;
673
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
5 PetscCall(STSetShift(eps->st,rho));
674
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
5 PetscCall(KSPGetConvergedReason(ksp,&reason));
675
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
5 PetscCheck(!reason,PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Second factorization failed");
676 }
677
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
140 PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_TRUE));
678 }
679 }
680 4157 eps->errest[eps->nconv] = relerr;
681
682 /* if relerr<tol, accept eigenpair */
683
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
4157 if (relerr<eps->tol) {
684 95 eps->nconv = eps->nconv + 1;
685
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
95 if (eps->nconv<eps->nev) {
686
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
55 PetscCall(EPSGetStartVector(eps,eps->nconv,&breakdown));
687
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
55 PetscCall(EPSGetLeftStartVector(eps,eps->nconv,&breakdownl));
688
2/4
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
55 if (breakdown || breakdownl) {
689 eps->reason = EPS_DIVERGED_BREAKDOWN;
690 PetscCall(PetscInfo(eps,"Unable to generate more start vectors\n"));
691 break;
692 }
693
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
55 PetscCall(BVBiorthonormalizeColumn(eps->V,eps->W,eps->nconv,NULL));
694 }
695 }
696
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4157 PetscCall(EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev)));
697
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.
4197 PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx));
698 }
699
700
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(DSSetDimensions(eps->ds,eps->nconv,0,0));
701
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
702
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.
8 PetscFunctionReturn(PETSC_SUCCESS);
703 }
704
705 270266 static PetscErrorCode EPSStopping_Power(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
706 {
707 270266 EPS_POWER *power = (EPS_POWER*)eps->data;
708 270266 SNESConvergedReason snesreason;
709
710
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
270266 PetscFunctionBegin;
711
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
270266 if (PetscUnlikely(power->update)) {
712
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(SNESGetConvergedReason(power->snes,&snesreason));
713
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
150 if (snesreason < 0) {
714 *reason = EPS_DIVERGED_BREAKDOWN;
715 PetscFunctionReturn(PETSC_SUCCESS);
716 }
717 }
718
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
270266 PetscCall(EPSStoppingBasic(eps,its,max_it,nconv,nev,reason,ctx));
719
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.
54166 PetscFunctionReturn(PETSC_SUCCESS);
720 }
721
722 176 static PetscErrorCode EPSBackTransform_Power(EPS eps)
723 {
724 176 EPS_POWER *power = (EPS_POWER*)eps->data;
725
726
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
176 PetscFunctionBegin;
727
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
176 if (power->nonlinear) eps->eigr[0] = 1.0/eps->eigr[0];
728
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.
176 else if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) PetscCall(EPSBackTransform_Default(eps));
729
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.
36 PetscFunctionReturn(PETSC_SUCCESS);
730 }
731
732 466 static PetscErrorCode EPSSetFromOptions_Power(EPS eps,PetscOptionItems PetscOptionsObject)
733 {
734 466 EPS_POWER *power = (EPS_POWER*)eps->data;
735 466 PetscBool flg,val;
736 466 EPSPowerShiftType shift;
737
738
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
466 PetscFunctionBegin;
739
1/12
✗ 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.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
466 PetscOptionsHeadBegin(PetscOptionsObject,"EPS Power Options");
740
741
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscOptionsEnum("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",EPSPowerShiftTypes,(PetscEnum)power->shift_type,(PetscEnum*)&shift,&flg));
742
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.
466 if (flg) PetscCall(EPSPowerSetShiftType(eps,shift));
743
744
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscOptionsBool("-eps_power_nonlinear","Use nonlinear inverse iteration","EPSPowerSetNonlinear",power->nonlinear,&val,&flg));
745
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.
466 if (flg) PetscCall(EPSPowerSetNonlinear(eps,val));
746
747
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscOptionsBool("-eps_power_update","Update residual monolithically","EPSPowerSetUpdate",power->update,&val,&flg));
748
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.
466 if (flg) PetscCall(EPSPowerSetUpdate(eps,val));
749
750
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscOptionsBool("-eps_power_sign_normalization","Normalize Bx with sign of first nonzero entry","EPSPowerSetSignNormalization",power->sign_normalization,&power->sign_normalization,&flg));
751
752
1/14
✗ Branch 0 not taken.
✓ 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.
466 PetscOptionsHeadEnd();
753
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.
94 PetscFunctionReturn(PETSC_SUCCESS);
754 }
755
756 80 static PetscErrorCode EPSPowerSetShiftType_Power(EPS eps,EPSPowerShiftType shift)
757 {
758 80 EPS_POWER *power = (EPS_POWER*)eps->data;
759
760
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
80 PetscFunctionBegin;
761
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
80 switch (shift) {
762 80 case EPS_POWER_SHIFT_CONSTANT:
763 case EPS_POWER_SHIFT_RAYLEIGH:
764 case EPS_POWER_SHIFT_WILKINSON:
765
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
80 if (power->shift_type != shift) {
766 50 power->shift_type = shift;
767 50 eps->state = EPS_STATE_INITIAL;
768 }
769 16 break;
770 default:
771 SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
772 }
773
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);
774 }
775
776 /*@
777 EPSPowerSetShiftType - Sets the type of shifts used during the power
778 iteration. This can be used to emulate the Rayleigh Quotient Iteration
779 (RQI) method.
780
781 Logically Collective
782
783 Input Parameters:
784 + eps - the eigenproblem solver context
785 - shift - the type of shift
786
787 Options Database Key:
788 . -eps_power_shift_type - Sets the shift type (either 'constant' or
789 'rayleigh' or 'wilkinson')
790
791 Notes:
792 By default, shifts are constant (EPS_POWER_SHIFT_CONSTANT) and the iteration
793 is the simple power method (or inverse iteration if a shift-and-invert
794 transformation is being used).
795
796 A variable shift can be specified (EPS_POWER_SHIFT_RAYLEIGH or
797 EPS_POWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
798 a cubic converging method such as RQI.
799
800 Level: advanced
801
802 .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
803 @*/
804 80 PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
805 {
806
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
80 PetscFunctionBegin;
807
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);
808
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.
80 PetscValidLogicalCollectiveEnum(eps,shift,2);
809
8/14
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
✓ 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.
80 PetscTryMethod(eps,"EPSPowerSetShiftType_C",(EPS,EPSPowerShiftType),(eps,shift));
810
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.
80 PetscFunctionReturn(PETSC_SUCCESS);
811 }
812
813 10 static PetscErrorCode EPSPowerGetShiftType_Power(EPS eps,EPSPowerShiftType *shift)
814 {
815 10 EPS_POWER *power = (EPS_POWER*)eps->data;
816
817
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
10 PetscFunctionBegin;
818 10 *shift = power->shift_type;
819
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.
10 PetscFunctionReturn(PETSC_SUCCESS);
820 }
821
822 /*@
823 EPSPowerGetShiftType - Gets the type of shifts used during the power
824 iteration.
825
826 Not Collective
827
828 Input Parameter:
829 . eps - the eigenproblem solver context
830
831 Output Parameter:
832 . shift - the type of shift
833
834 Level: advanced
835
836 .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
837 @*/
838 10 PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
839 {
840
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
10 PetscFunctionBegin;
841
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.
10 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
842
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.
10 PetscAssertPointer(shift,2);
843
9/16
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
10 PetscUseMethod(eps,"EPSPowerGetShiftType_C",(EPS,EPSPowerShiftType*),(eps,shift));
844
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.
10 PetscFunctionReturn(PETSC_SUCCESS);
845 }
846
847 310 static PetscErrorCode EPSPowerSetNonlinear_Power(EPS eps,PetscBool nonlinear)
848 {
849 310 EPS_POWER *power = (EPS_POWER*)eps->data;
850
851
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
310 PetscFunctionBegin;
852
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
310 if (power->nonlinear != nonlinear) {
853 310 power->nonlinear = nonlinear;
854 310 eps->useds = PetscNot(nonlinear);
855
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
310 eps->ops->setupsort = nonlinear? NULL: EPSSetUpSort_Default;
856 310 eps->state = EPS_STATE_INITIAL;
857 }
858
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
310 PetscFunctionReturn(PETSC_SUCCESS);
859 }
860
861 /*@
862 EPSPowerSetNonlinear - Sets a flag to indicate that the problem is nonlinear.
863
864 Logically Collective
865
866 Input Parameters:
867 + eps - the eigenproblem solver context
868 - nonlinear - whether the problem is nonlinear or not
869
870 Options Database Key:
871 . -eps_power_nonlinear - Sets the nonlinear flag
872
873 Notes:
874 If this flag is set, the solver assumes that the problem is nonlinear,
875 that is, the operators that define the eigenproblem are not constant
876 matrices, but depend on the eigenvector, A(x)*x=lambda*B(x)*x. This is
877 different from the case of nonlinearity with respect to the eigenvalue
878 (use the NEP solver class for this kind of problems).
879
880 The way in which nonlinear operators are specified is very similar to
881 the case of PETSc's SNES solver. The difference is that the callback
882 functions are provided via composed functions "formFunction" and
883 "formJacobian" in each of the matrix objects passed as arguments of
884 EPSSetOperators(). The application context required for these functions
885 can be attached via a composed PetscContainer.
886
887 Level: advanced
888
889 .seealso: EPSPowerGetNonlinear(), EPSSetOperators()
890 @*/
891 310 PetscErrorCode EPSPowerSetNonlinear(EPS eps,PetscBool nonlinear)
892 {
893
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
310 PetscFunctionBegin;
894
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
310 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
895
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.
310 PetscValidLogicalCollectiveBool(eps,nonlinear,2);
896
8/14
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
✓ 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.
310 PetscTryMethod(eps,"EPSPowerSetNonlinear_C",(EPS,PetscBool),(eps,nonlinear));
897
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
310 PetscFunctionReturn(PETSC_SUCCESS);
898 }
899
900 180 static PetscErrorCode EPSPowerGetNonlinear_Power(EPS eps,PetscBool *nonlinear)
901 {
902 180 EPS_POWER *power = (EPS_POWER*)eps->data;
903
904
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
180 PetscFunctionBegin;
905 180 *nonlinear = power->nonlinear;
906
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.
180 PetscFunctionReturn(PETSC_SUCCESS);
907 }
908
909 /*@
910 EPSPowerGetNonlinear - Returns a flag indicating if the problem is nonlinear.
911
912 Not Collective
913
914 Input Parameter:
915 . eps - the eigenproblem solver context
916
917 Output Parameter:
918 . nonlinear - the nonlinear flag
919
920 Level: advanced
921
922 .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
923 @*/
924 180 PetscErrorCode EPSPowerGetNonlinear(EPS eps,PetscBool *nonlinear)
925 {
926
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
180 PetscFunctionBegin;
927
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.
180 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
928
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.
180 PetscAssertPointer(nonlinear,2);
929
9/16
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
180 PetscUseMethod(eps,"EPSPowerGetNonlinear_C",(EPS,PetscBool*),(eps,nonlinear));
930
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.
180 PetscFunctionReturn(PETSC_SUCCESS);
931 }
932
933 130 static PetscErrorCode EPSPowerSetUpdate_Power(EPS eps,PetscBool update)
934 {
935 130 EPS_POWER *power = (EPS_POWER*)eps->data;
936
937
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
130 PetscFunctionBegin;
938
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
130 PetscCheck(power->nonlinear,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"This option does not make sense for linear problems");
939 130 power->update = update;
940 130 eps->state = EPS_STATE_INITIAL;
941
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.
130 PetscFunctionReturn(PETSC_SUCCESS);
942 }
943
944 /*@
945 EPSPowerSetUpdate - Sets a flag to indicate that the residual is updated monolithically
946 for nonlinear problems. This potentially has a better convergence.
947
948 Logically Collective
949
950 Input Parameters:
951 + eps - the eigenproblem solver context
952 - update - whether the residual is updated monolithically or not
953
954 Options Database Key:
955 . -eps_power_update - Sets the update flag
956
957 Level: advanced
958
959 .seealso: EPSPowerGetUpdate(), EPSPowerGetNonlinear(), EPSSetOperators()
960 @*/
961 130 PetscErrorCode EPSPowerSetUpdate(EPS eps,PetscBool update)
962 {
963
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
130 PetscFunctionBegin;
964
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.
130 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
965
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.
130 PetscValidLogicalCollectiveBool(eps,update,2);
966
8/14
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
✓ 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.
130 PetscTryMethod(eps,"EPSPowerSetUpdate_C",(EPS,PetscBool),(eps,update));
967
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.
130 PetscFunctionReturn(PETSC_SUCCESS);
968 }
969
970 180 static PetscErrorCode EPSPowerGetUpdate_Power(EPS eps,PetscBool *update)
971 {
972 180 EPS_POWER *power = (EPS_POWER*)eps->data;
973
974
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
180 PetscFunctionBegin;
975 180 *update = power->update;
976
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.
180 PetscFunctionReturn(PETSC_SUCCESS);
977 }
978
979 /*@
980 EPSPowerGetUpdate - Returns a flag indicating if the residual is updated monolithically
981 for nonlinear problems.
982
983 Not Collective
984
985 Input Parameter:
986 . eps - the eigenproblem solver context
987
988 Output Parameter:
989 . update - the update flag
990
991 Level: advanced
992
993 .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
994 @*/
995 180 PetscErrorCode EPSPowerGetUpdate(EPS eps,PetscBool *update)
996 {
997
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
180 PetscFunctionBegin;
998
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.
180 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
999
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.
180 PetscAssertPointer(update,2);
1000
9/16
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
180 PetscUseMethod(eps,"EPSPowerGetUpdate_C",(EPS,PetscBool*),(eps,update));
1001
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.
180 PetscFunctionReturn(PETSC_SUCCESS);
1002 }
1003
1004 310 static PetscErrorCode EPSPowerSetSignNormalization_Power(EPS eps,PetscBool sign_normalization)
1005 {
1006 310 EPS_POWER *power = (EPS_POWER*)eps->data;
1007
1008
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
310 PetscFunctionBegin;
1009 310 power->sign_normalization = sign_normalization;
1010
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
310 PetscFunctionReturn(PETSC_SUCCESS);
1011 }
1012
1013 /*@
1014 EPSPowerSetSignNormalization - Sets a flag to indicate whether the Bx vector
1015 should be normalized by the sign of the first non-zero element in the vector.
1016 E.g., if this is true, the post-normalization value of the first non-zero element
1017 in the vector is guaranteed to be positive.
1018
1019 Logically Collective
1020
1021 Input Parameters:
1022 + eps - the eigenproblem solver context
1023 - sign_normalization - whether Bx should be multiplied by the sign of the first non-zero
1024 element when performing normalization steps
1025
1026 Options Database Key:
1027 . -eps_power_sign_normalization - Sets the sign normalization flag
1028
1029 Level: advanced
1030
1031 .seealso: EPSPowerGetSignNormalization()
1032 @*/
1033 310 PetscErrorCode EPSPowerSetSignNormalization(EPS eps,PetscBool sign_normalization)
1034 {
1035
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
310 PetscFunctionBegin;
1036
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
310 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1037
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.
310 PetscValidLogicalCollectiveBool(eps,sign_normalization,2);
1038
8/14
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
✓ 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.
310 PetscTryMethod(eps,"EPSPowerSetSignNormalization_C",(EPS,PetscBool),(eps,sign_normalization));
1039
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
310 PetscFunctionReturn(PETSC_SUCCESS);
1040 }
1041
1042 130 static PetscErrorCode EPSPowerGetSignNormalization_Power(EPS eps,PetscBool *sign_normalization)
1043 {
1044 130 EPS_POWER *power = (EPS_POWER*)eps->data;
1045
1046
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
130 PetscFunctionBegin;
1047 130 *sign_normalization = power->sign_normalization;
1048
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.
130 PetscFunctionReturn(PETSC_SUCCESS);
1049 }
1050
1051 /*@
1052 EPSPowerGetSignNormalization - Returns a flag indicating whether the Bx vector
1053 is normalized by the sign of the first non-zero element in the vector. E.g.,
1054 if this is true, the post-normalization value of the first non-zero element in
1055 the vector is guaranteed to be positive.
1056
1057 Not Collective
1058
1059 Input Parameter:
1060 . eps - the eigenproblem solver context
1061
1062 Output Parameter:
1063 . sign_normalization - the sign normalization flag
1064
1065 Level: advanced
1066
1067 .seealso: EPSPowerSetSignNormalization()
1068 @*/
1069 130 PetscErrorCode EPSPowerGetSignNormalization(EPS eps,PetscBool *sign_normalization)
1070 {
1071
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
130 PetscFunctionBegin;
1072
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.
130 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1073
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.
130 PetscAssertPointer(sign_normalization,2);
1074
9/16
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
130 PetscUseMethod(eps,"EPSPowerGetSignNormalization_C",(EPS,PetscBool*),(eps,sign_normalization));
1075
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.
130 PetscFunctionReturn(PETSC_SUCCESS);
1076 }
1077
1078 static PetscErrorCode EPSPowerSetSNES_Power(EPS eps,SNES snes)
1079 {
1080 EPS_POWER *power = (EPS_POWER*)eps->data;
1081
1082 PetscFunctionBegin;
1083 PetscCall(PetscObjectReference((PetscObject)snes));
1084 PetscCall(SNESDestroy(&power->snes));
1085 power->snes = snes;
1086 eps->state = EPS_STATE_INITIAL;
1087 PetscFunctionReturn(PETSC_SUCCESS);
1088 }
1089
1090 /*@
1091 EPSPowerSetSNES - Associate a nonlinear solver object (SNES) to the
1092 eigenvalue solver (to be used in nonlinear inverse iteration).
1093
1094 Collective
1095
1096 Input Parameters:
1097 + eps - the eigenvalue solver
1098 - snes - the nonlinear solver object
1099
1100 Level: advanced
1101
1102 .seealso: EPSPowerGetSNES()
1103 @*/
1104 PetscErrorCode EPSPowerSetSNES(EPS eps,SNES snes)
1105 {
1106 PetscFunctionBegin;
1107 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1108 PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
1109 PetscCheckSameComm(eps,1,snes,2);
1110 PetscTryMethod(eps,"EPSPowerSetSNES_C",(EPS,SNES),(eps,snes));
1111 PetscFunctionReturn(PETSC_SUCCESS);
1112 }
1113
1114 440 static PetscErrorCode EPSPowerGetSNES_Power(EPS eps,SNES *snes)
1115 {
1116 440 EPS_POWER *power = (EPS_POWER*)eps->data;
1117
1118
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
440 PetscFunctionBegin;
1119
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
440 if (!power->snes) {
1120
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
310 PetscCall(SNESCreate(PetscObjectComm((PetscObject)eps),&power->snes));
1121
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
310 PetscCall(PetscObjectIncrementTabLevel((PetscObject)power->snes,(PetscObject)eps,1));
1122
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
310 PetscCall(SNESSetOptionsPrefix(power->snes,((PetscObject)eps)->prefix));
1123
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
310 PetscCall(SNESAppendOptionsPrefix(power->snes,"eps_power_"));
1124
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
310 PetscCall(PetscObjectSetOptions((PetscObject)power->snes,((PetscObject)eps)->options));
1125 }
1126 440 *snes = power->snes;
1127
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.
440 PetscFunctionReturn(PETSC_SUCCESS);
1128 }
1129
1130 /*@
1131 EPSPowerGetSNES - Retrieve the nonlinear solver object (SNES) associated
1132 with the eigenvalue solver.
1133
1134 Not Collective
1135
1136 Input Parameter:
1137 . eps - the eigenvalue solver
1138
1139 Output Parameter:
1140 . snes - the nonlinear solver object
1141
1142 Level: advanced
1143
1144 .seealso: EPSPowerSetSNES()
1145 @*/
1146 440 PetscErrorCode EPSPowerGetSNES(EPS eps,SNES *snes)
1147 {
1148
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
440 PetscFunctionBegin;
1149
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.
440 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1150
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.
440 PetscAssertPointer(snes,2);
1151
9/16
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
440 PetscUseMethod(eps,"EPSPowerGetSNES_C",(EPS,SNES*),(eps,snes));
1152
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.
440 PetscFunctionReturn(PETSC_SUCCESS);
1153 }
1154
1155 476 static PetscErrorCode EPSReset_Power(EPS eps)
1156 {
1157 476 EPS_POWER *power = (EPS_POWER*)eps->data;
1158
1159
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
476 PetscFunctionBegin;
1160
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.
476 if (power->snes) PetscCall(SNESReset(power->snes));
1161
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.
96 PetscFunctionReturn(PETSC_SUCCESS);
1162 }
1163
1164 466 static PetscErrorCode EPSDestroy_Power(EPS eps)
1165 {
1166 466 EPS_POWER *power = (EPS_POWER*)eps->data;
1167
1168
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
466 PetscFunctionBegin;
1169
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.
466 if (power->nonlinear) PetscCall(SNESDestroy(&power->snes));
1170
5/8
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
466 PetscCall(PetscFree(eps->data));
1171
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",NULL));
1172
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",NULL));
1173
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",NULL));
1174
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",NULL));
1175
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",NULL));
1176
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",NULL));
1177
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSignNormalization_C",NULL));
1178
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSignNormalization_C",NULL));
1179
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",NULL));
1180
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",NULL));
1181
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.
94 PetscFunctionReturn(PETSC_SUCCESS);
1182 }
1183
1184 static PetscErrorCode EPSView_Power(EPS eps,PetscViewer viewer)
1185 {
1186 EPS_POWER *power = (EPS_POWER*)eps->data;
1187 PetscBool isascii;
1188
1189 PetscFunctionBegin;
1190 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
1191 if (isascii) {
1192 if (power->nonlinear) {
1193 if (power->sign_normalization) PetscCall(PetscViewerASCIIPrintf(viewer," normalizing Bx by the sign of the first nonzero element\n"));
1194 else PetscCall(PetscViewerASCIIPrintf(viewer," not normalizing Bx by the sign of the first nonzero element\n"));
1195 PetscCall(PetscViewerASCIIPrintf(viewer," using nonlinear inverse iteration\n"));
1196 if (power->update) PetscCall(PetscViewerASCIIPrintf(viewer," updating the residual monolithically\n"));
1197 if (!power->snes) PetscCall(EPSPowerGetSNES(eps,&power->snes));
1198 PetscCall(PetscViewerASCIIPushTab(viewer));
1199 PetscCall(SNESView(power->snes,viewer));
1200 PetscCall(PetscViewerASCIIPopTab(viewer));
1201 } else PetscCall(PetscViewerASCIIPrintf(viewer," %s shifts\n",EPSPowerShiftTypes[power->shift_type]));
1202 }
1203 PetscFunctionReturn(PETSC_SUCCESS);
1204 }
1205
1206 376 static PetscErrorCode EPSComputeVectors_Power(EPS eps)
1207 {
1208 376 EPS_POWER *power = (EPS_POWER*)eps->data;
1209
1210
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
376 PetscFunctionBegin;
1211
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
376 if (eps->twosided) {
1212
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(EPSComputeVectors_Twosided(eps));
1213
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(BVNormalize(eps->V,NULL));
1214
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
40 PetscCall(BVNormalize(eps->W,NULL));
1215
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.
336 } else if (!power->nonlinear) PetscCall(EPSComputeVectors_Schur(eps));
1216
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.
76 PetscFunctionReturn(PETSC_SUCCESS);
1217 }
1218
1219 972 static PetscErrorCode EPSSetDefaultST_Power(EPS eps)
1220 {
1221 972 EPS_POWER *power = (EPS_POWER*)eps->data;
1222 972 KSP ksp;
1223 972 PC pc;
1224
1225
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
972 PetscFunctionBegin;
1226
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
972 if (power->nonlinear) {
1227 640 eps->categ=EPS_CATEGORY_PRECOND;
1228
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
640 PetscCall(STGetKSP(eps->st,&ksp));
1229 /* Set ST as STPRECOND so it can carry one preconditioning matrix
1230 * It is useful when A and B are shell matrices
1231 */
1232
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
640 PetscCall(STSetType(eps->st,STPRECOND));
1233
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
640 PetscCall(KSPGetPC(ksp,&pc));
1234
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
640 PetscCall(PCSetType(pc,PCNONE));
1235 }
1236
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.
196 PetscFunctionReturn(PETSC_SUCCESS);
1237 }
1238
1239 466 SLEPC_EXTERN PetscErrorCode EPSCreate_Power(EPS eps)
1240 {
1241 466 EPS_POWER *ctx;
1242
1243
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
466 PetscFunctionBegin;
1244
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscNew(&ctx));
1245 466 eps->data = (void*)ctx;
1246
1247 466 eps->useds = PETSC_TRUE;
1248 466 eps->categ = EPS_CATEGORY_OTHER;
1249
1250 466 eps->ops->setup = EPSSetUp_Power;
1251 466 eps->ops->setupsort = EPSSetUpSort_Default;
1252 466 eps->ops->setfromoptions = EPSSetFromOptions_Power;
1253 466 eps->ops->reset = EPSReset_Power;
1254 466 eps->ops->destroy = EPSDestroy_Power;
1255 466 eps->ops->view = EPSView_Power;
1256 466 eps->ops->backtransform = EPSBackTransform_Power;
1257 466 eps->ops->computevectors = EPSComputeVectors_Power;
1258 466 eps->ops->setdefaultst = EPSSetDefaultST_Power;
1259 466 eps->stopping = EPSStopping_Power;
1260 466 ctx->sign_normalization = PETSC_TRUE;
1261
1262
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",EPSPowerSetShiftType_Power));
1263
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",EPSPowerGetShiftType_Power));
1264
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",EPSPowerSetNonlinear_Power));
1265
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",EPSPowerGetNonlinear_Power));
1266
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",EPSPowerSetUpdate_Power));
1267
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",EPSPowerGetUpdate_Power));
1268
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSignNormalization_C",EPSPowerSetSignNormalization_Power));
1269
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSignNormalization_C",EPSPowerGetSignNormalization_Power));
1270
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",EPSPowerSetSNES_Power));
1271
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
466 PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",EPSPowerGetSNES_Power));
1272
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.
94 PetscFunctionReturn(PETSC_SUCCESS);
1273 }
1274