GCC Code Coverage Report


Directory: ./
File: src/eps/interface/epssetup.c
Date: 2026-02-22 03:58:10
Exec Total Coverage
Lines: 383 387 99.0%
Functions: 18 18 100.0%
Branches: 1007 1855 54.3%

Line Branch Exec Source
1 /*
2 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 SLEPc - Scalable Library for Eigenvalue Problem Computations
4 Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5
6 This file is part of SLEPc.
7 SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 */
10 /*
11 EPS routines related to problem setup
12 */
13
14 #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
15
16 /*
17 Let the solver choose the ST type that should be used by default,
18 otherwise set it to SHIFT.
19 This is called at EPSSetFromOptions (before STSetFromOptions)
20 and also at EPSSetUp (in case EPSSetFromOptions was not called).
21 */
22 13448 PetscErrorCode EPSSetDefaultST(EPS eps)
23 {
24
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
13448 PetscFunctionBegin;
25
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
13448 PetscTryTypeMethod(eps,setdefaultst);
26
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
13448 if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STSHIFT));
27
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.
3252 PetscFunctionReturn(PETSC_SUCCESS);
28 }
29
30 /*
31 This is done by preconditioned eigensolvers that use the PC only.
32 It sets STPRECOND with KSPPREONLY.
33 */
34 836 PetscErrorCode EPSSetDefaultST_Precond(EPS eps)
35 {
36 836 KSP ksp;
37
38
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
836 PetscFunctionBegin;
39
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
836 if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STPRECOND));
40
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
836 PetscCall(STGetKSP(eps->st,&ksp));
41
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
836 if (!((PetscObject)ksp)->type_name) PetscCall(KSPSetType(ksp,KSPPREONLY));
42
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.
208 PetscFunctionReturn(PETSC_SUCCESS);
43 }
44
45 /*
46 This is done by preconditioned eigensolvers that can also use the KSP.
47 It sets STPRECOND with the default KSP (GMRES) and maxit=5.
48 */
49 794 PetscErrorCode EPSSetDefaultST_GMRES(EPS eps)
50 {
51 794 KSP ksp;
52
53
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
794 PetscFunctionBegin;
54
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
794 if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STPRECOND));
55
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
794 PetscCall(STPrecondSetKSPHasMat(eps->st,PETSC_TRUE));
56
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
794 PetscCall(STGetKSP(eps->st,&ksp));
57
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
794 if (!((PetscObject)ksp)->type_name) {
58
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
322 PetscCall(KSPSetType(ksp,KSPGMRES));
59
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
322 PetscCall(KSPSetTolerances(ksp,PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT,5));
60 }
61
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.
216 PetscFunctionReturn(PETSC_SUCCESS);
62 }
63
64 #if defined(SLEPC_HAVE_SCALAPACK) || defined(SLEPC_HAVE_ELPA) || defined(SLEPC_HAVE_ELEMENTAL) || defined(SLEPC_HAVE_EVSL)
65 /*
66 This is for direct eigensolvers that work with A and B directly, so
67 no need to factorize B.
68 */
69 173 PetscErrorCode EPSSetDefaultST_NoFactor(EPS eps)
70 {
71 173 KSP ksp;
72 173 PC pc;
73
74
1/2
✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
173 PetscFunctionBegin;
75
6/8
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
173 if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STSHIFT));
76
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.
173 PetscCall(STGetKSP(eps->st,&ksp));
77
6/8
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
173 if (!((PetscObject)ksp)->type_name) PetscCall(KSPSetType(ksp,KSPPREONLY));
78
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.
173 PetscCall(KSPGetPC(ksp,&pc));
79
6/8
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
173 if (!((PetscObject)pc)->type_name) PetscCall(PCSetType(pc,PCNONE));
80
6/12
✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
24 PetscFunctionReturn(PETSC_SUCCESS);
81 }
82 #endif
83
84 /*
85 Check that the ST selected by the user is compatible with the EPS solver and options
86 */
87 7998 static PetscErrorCode EPSCheckCompatibleST(EPS eps)
88 {
89 7998 PetscBool precond,shift,sinvert,cayley,lyapii;
90 #if defined(PETSC_USE_COMPLEX)
91 3842 PetscScalar sigma;
92 #endif
93
94
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
7998 PetscFunctionBegin;
95
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&precond));
96
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&shift));
97
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&sinvert));
98
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&cayley));
99
100 /* preconditioned eigensolvers */
101
3/6
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
7998 PetscCheck(eps->categ!=EPS_CATEGORY_PRECOND || precond,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires ST=PRECOND");
102
3/6
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
7998 PetscCheck(eps->categ==EPS_CATEGORY_PRECOND || !precond,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"STPRECOND is intended for preconditioned eigensolvers only");
103
104 /* harmonic extraction */
105
5/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 8 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
7998 PetscCheck(precond || shift || !eps->extraction || eps->extraction==EPS_RITZ,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cannot use a spectral transformation combined with harmonic extraction");
106
107 /* real shifts in Hermitian problems */
108 #if defined(PETSC_USE_COMPLEX)
109
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 3 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
3842 PetscCall(STGetShift(eps->st,&sigma));
110
3/6
✓ Branch 0 taken 4 times.
✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
3842 PetscCheck(!eps->ishermitian || PetscImaginaryPart(sigma)==0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Hermitian problems are not compatible with complex shifts");
111 #endif
112
113 /* Cayley with PGNHEP */
114
3/6
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
7998 PetscCheck(!cayley || eps->problem_type!=EPS_PGNHEP,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");
115
116 /* make sure that the user does not specify smallest magnitude with shift-and-invert */
117
6/6
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 8 times.
✓ Branch 5 taken 8 times.
7998 if ((cayley || sinvert) && (eps->categ==EPS_CATEGORY_KRYLOV || eps->categ==EPS_CATEGORY_OTHER)) {
118
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1025 PetscCall(PetscObjectTypeCompare((PetscObject)eps,EPSLYAPII,&lyapii));
119
10/14
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 8 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 8 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 8 times.
✓ Branch 9 taken 8 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 8 times.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
1025 PetscCheck(lyapii || eps->which==EPS_TARGET_MAGNITUDE || eps->which==EPS_TARGET_REAL || eps->which==EPS_TARGET_IMAGINARY || eps->which==EPS_ALL || eps->which==EPS_WHICH_USER,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"Shift-and-invert requires a target 'which' (see EPSSetWhichEigenpairs), for instance -st_type sinvert -eps_target 0 -eps_target_magnitude");
120 }
121
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.
1926 PetscFunctionReturn(PETSC_SUCCESS);
122 }
123
124 /*
125 MatEstimateSpectralRange_EPS: estimate the spectral range [left,right] of a
126 symmetric/Hermitian matrix A using an auxiliary EPS object
127 */
128 112 PetscErrorCode MatEstimateSpectralRange_EPS(Mat A,PetscReal *left,PetscReal *right)
129 {
130 112 PetscInt nconv;
131 112 PetscScalar eig0;
132 112 PetscReal tol=1e-3,errest=tol;
133 112 EPS eps;
134
135
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
112 PetscFunctionBegin;
136 112 *left = 0.0; *right = 0.0;
137
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSCreate(PetscObjectComm((PetscObject)A),&eps));
138
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSSetOptionsPrefix(eps,"eps_filter_"));
139
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSSetOperators(eps,A,NULL));
140
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSSetProblemType(eps,EPS_HEP));
141
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSSetTolerances(eps,tol,50));
142
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_ABS));
143
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
144
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSSolve(eps));
145
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSGetConverged(eps,&nconv));
146
1/2
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
112 if (nconv>0) {
147
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSGetEigenvalue(eps,0,&eig0,NULL));
148
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSGetErrorEstimate(eps,0,&errest));
149 } else eig0 = eps->eigr[0];
150 112 *left = PetscRealPart(eig0)-errest;
151
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
152
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSSolve(eps));
153
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSGetConverged(eps,&nconv));
154
1/2
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
112 if (nconv>0) {
155
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSGetEigenvalue(eps,0,&eig0,NULL));
156
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSGetErrorEstimate(eps,0,&errest));
157 } else eig0 = eps->eigr[0];
158 112 *right = PetscRealPart(eig0)+errest;
159
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
112 PetscCall(EPSDestroy(&eps));
160
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);
161 }
162
163 /*
164 EPSSetUpSort_Basic: configure the EPS sorting criterion according to 'which'
165 */
166 7513 PetscErrorCode EPSSetUpSort_Basic(EPS eps)
167 {
168
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
7513 PetscFunctionBegin;
169
11/11
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 8 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 8 times.
✓ Branch 7 taken 8 times.
✓ Branch 8 taken 4 times.
✓ Branch 9 taken 8 times.
✓ Branch 10 taken 8 times.
7513 switch (eps->which) {
170 3098 case EPS_LARGEST_MAGNITUDE:
171 3098 eps->sc->comparison = SlepcCompareLargestMagnitude;
172 3098 eps->sc->comparisonctx = NULL;
173 3098 break;
174 374 case EPS_SMALLEST_MAGNITUDE:
175 374 eps->sc->comparison = SlepcCompareSmallestMagnitude;
176 374 eps->sc->comparisonctx = NULL;
177 374 break;
178 945 case EPS_LARGEST_REAL:
179 945 eps->sc->comparison = SlepcCompareLargestReal;
180 945 eps->sc->comparisonctx = NULL;
181 945 break;
182 1311 case EPS_SMALLEST_REAL:
183 1311 eps->sc->comparison = SlepcCompareSmallestReal;
184 1311 eps->sc->comparisonctx = NULL;
185 1311 break;
186 12 case EPS_LARGEST_IMAGINARY:
187 12 eps->sc->comparison = SlepcCompareLargestImaginary;
188 12 eps->sc->comparisonctx = NULL;
189 12 break;
190 8 case EPS_SMALLEST_IMAGINARY:
191 8 eps->sc->comparison = SlepcCompareSmallestImaginary;
192 8 eps->sc->comparisonctx = NULL;
193 8 break;
194 999 case EPS_TARGET_MAGNITUDE:
195 999 eps->sc->comparison = SlepcCompareTargetMagnitude;
196 999 eps->sc->comparisonctx = &eps->target;
197 999 break;
198 44 case EPS_TARGET_REAL:
199 44 eps->sc->comparison = SlepcCompareTargetReal;
200 44 eps->sc->comparisonctx = &eps->target;
201 44 break;
202 4 case EPS_TARGET_IMAGINARY:
203 4 eps->sc->comparison = SlepcCompareTargetImaginary;
204 4 eps->sc->comparisonctx = &eps->target;
205 4 break;
206 515 case EPS_ALL:
207 515 eps->sc->comparison = SlepcCompareSmallestReal;
208 515 eps->sc->comparisonctx = NULL;
209 515 break;
210 case EPS_WHICH_USER:
211 break;
212 }
213 7513 eps->sc->map = NULL;
214 7513 eps->sc->mapobj = NULL;
215
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.
7513 PetscFunctionReturn(PETSC_SUCCESS);
216 }
217
218 /*
219 EPSSetUpSort_Default: configure both EPS and DS sorting criterion
220 */
221 7270 PetscErrorCode EPSSetUpSort_Default(EPS eps)
222 {
223 7270 SlepcSC sc;
224 7270 PetscBool istrivial;
225
226
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
7270 PetscFunctionBegin;
227 /* fill sorting criterion context */
228
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7270 PetscCall(EPSSetUpSort_Basic(eps));
229 /* fill sorting criterion for DS */
230
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7270 PetscCall(DSGetSlepcSC(eps->ds,&sc));
231
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7270 PetscCall(RGIsTrivial(eps->rg,&istrivial));
232
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7270 sc->rg = istrivial? NULL: eps->rg;
233 7270 sc->comparison = eps->sc->comparison;
234 7270 sc->comparisonctx = eps->sc->comparisonctx;
235 7270 sc->map = SlepcMap_ST;
236 7270 sc->mapobj = (PetscObject)eps->st;
237
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
7270 PetscFunctionReturn(PETSC_SUCCESS);
238 }
239
240 /*@
241 EPSSetDSType - Sets the type of the internal `DS` object based on the current
242 settings of the eigenvalue solver.
243
244 Collective
245
246 Input Parameter:
247 . eps - the linear eigensolver context
248
249 Note:
250 This function need not be called explicitly, since it will be called at
251 both `EPSSetFromOptions()` and `EPSSetUp()`.
252
253 Level: developer
254
255 .seealso: [](ch:eps), `EPSSetFromOptions()`, `EPSSetUp()`
256 @*/
257 12296 PetscErrorCode EPSSetDSType(EPS eps)
258 {
259
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
12296 PetscFunctionBegin;
260
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.
12296 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
261
1/8
✗ Branch 0 not taken.
✓ Branch 1 taken 8 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.
12296 PetscTryTypeMethod(eps,setdstype);
262
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.
2964 PetscFunctionReturn(PETSC_SUCCESS);
263 }
264
265 /*@
266 EPSSetUp - Sets up all the internal data structures necessary for the
267 execution of the eigensolver. Then calls `STSetUp()` for any set-up
268 operations associated to the internal `ST` object.
269
270 Collective
271
272 Input Parameter:
273 . eps - the linear eigensolver context
274
275 Notes:
276 This function need not be called explicitly in most cases, since `EPSSolve()`
277 calls it. It can be useful when one wants to measure the set-up time
278 separately from the solve time.
279
280 Level: developer
281
282 .seealso: [](ch:eps), `EPSCreate()`, `EPSSolve()`, `EPSDestroy()`, `STSetUp()`, `EPSSetInitialSpace()`, `EPSSetDeflationSpace()`
283 @*/
284 9472 PetscErrorCode EPSSetUp(EPS eps)
285 {
286 9472 Mat A,B;
287 9472 PetscInt k,nmat;
288 9472 PetscBool flg;
289 9472 EPSStoppingCtx ctx;
290
291
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
9472 PetscFunctionBegin;
292
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.
9472 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
293
8/14
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 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.
9472 if (eps->state) PetscFunctionReturn(PETSC_SUCCESS);
294
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(PetscLogEventBegin(EPS_SetUp,eps,0,0,0));
295
296 /* reset the convergence flag from the previous solves */
297 7998 eps->reason = EPS_CONVERGED_ITERATING;
298
299 /* Set default solver type (EPSSetFromOptions was not called) */
300
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
7998 if (!((PetscObject)eps)->type_name) PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
301
1/8
✗ Branch 0 not taken.
✓ Branch 1 taken 8 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.
7998 if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
302
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(EPSSetDefaultST(eps));
303
304
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(STSetTransform(eps->st,PETSC_TRUE));
305
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(STSetStructured(eps->st,PETSC_FALSE));
306
8/10
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 6 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
7998 if (eps->useds && !eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
307
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
7998 if (eps->useds) PetscCall(EPSSetDSType(eps));
308
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (eps->twosided) {
309
1/8
✗ Branch 0 not taken.
✓ Branch 1 taken 8 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.
169 PetscCheck(!eps->ishermitian || (eps->isgeneralized && !eps->ispositive),PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Two-sided methods are not intended for %s problems",SLEPC_STRING_HERMITIAN);
310 }
311
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
7998 if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
312
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
7998 if (!((PetscObject)eps->rg)->type_name) PetscCall(RGSetType(eps->rg,RGINTERVAL));
313
314 /* Set problem dimensions */
315
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(STGetNumMatrices(eps->st,&nmat));
316
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
7998 PetscCheck(nmat,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSSetOperators must be called first");
317
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(STMatGetSize(eps->st,&eps->n,NULL));
318
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(STMatGetLocalSize(eps->st,&eps->nloc,NULL));
319
320 /* Set default problem type */
321
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (!eps->problem_type) {
322
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
365 if (nmat==1) PetscCall(EPSSetProblemType(eps,EPS_NHEP));
323
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
34 else PetscCall(EPSSetProblemType(eps,EPS_GNHEP));
324
4/4
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 3 times.
✓ Branch 3 taken 8 times.
7633 } else if (nmat==1 && eps->isgeneralized) {
325
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
3 PetscCall(PetscInfo(eps,"Eigenproblem set as generalized but no matrix B was provided; reverting to a standard eigenproblem\n"));
326 3 eps->isgeneralized = PETSC_FALSE;
327
1/2
✓ Branch 0 taken 3 times.
✗ Branch 1 not taken.
6 eps->problem_type = eps->ishermitian? EPS_HEP: EPS_NHEP;
328
3/6
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
7630 } else PetscCheck(nmat==1 || eps->isgeneralized,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Inconsistent EPS state: the problem type does not match the number of matrices");
329
330
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (eps->isstructured) {
331 /* make sure the user has set the appropriate matrix */
332
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
408 PetscCall(STGetMatrix(eps->st,0,&A));
333
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
408 if (eps->problem_type==EPS_BSE) PetscCall(SlepcCheckMatStruct(A,SLEPC_MAT_STRUCT_BSE,NULL));
334
6/8
✓ Branch 0 taken 4 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 3 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
408 if (eps->problem_type==EPS_HAMILT) PetscCall(SlepcCheckMatStruct(A,SLEPC_MAT_STRUCT_HAMILT,NULL));
335 }
336
337 /* safeguard for small problems */
338
2/14
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
7998 if (eps->n == 0) PetscFunctionReturn(PETSC_SUCCESS);
339
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
7998 if (eps->nev > eps->n) eps->nev = eps->n;
340
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (eps->problem_type == EPS_BSE) {
341
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
380 if (2*eps->ncv > eps->n) eps->ncv = eps->n/2;
342 } else {
343
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7618 if (eps->ncv > eps->n) eps->ncv = eps->n;
344 }
345
346 /* check some combinations of eps->which */
347
7/10
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 8 times.
✓ Branch 5 taken 8 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 8 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
7998 PetscCheck(!eps->ishermitian || (eps->isgeneralized && !eps->ispositive) || (eps->which!=EPS_LARGEST_IMAGINARY && eps->which!=EPS_SMALLEST_IMAGINARY && eps->which!=EPS_TARGET_IMAGINARY),PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Sorting the eigenvalues along the imaginary axis is not allowed when all eigenvalues are real");
348
349 /* initialization of matrix norms */
350
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (eps->conv==EPS_CONV_NORM) {
351
1/2
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
382 if (!eps->nrma) {
352
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
382 PetscCall(STGetMatrix(eps->st,0,&A));
353
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
382 PetscCall(MatNorm(A,NORM_INFINITY,&eps->nrma));
354 }
355
3/4
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✗ Branch 3 not taken.
382 if (nmat>1 && !eps->nrmb) {
356
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
374 PetscCall(STGetMatrix(eps->st,1,&B));
357
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
374 PetscCall(MatNorm(B,NORM_INFINITY,&eps->nrmb));
358 }
359 }
360
361 /* call specific solver setup */
362
5/10
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 6 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
7998 PetscUseTypeMethod(eps,setup);
363
364 /* threshold stopping test */
365
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (eps->stop==EPS_STOP_THRESHOLD) {
366
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
120 PetscCheck(eps->thres!=PETSC_MIN_REAL,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"Must give a threshold value with EPSSetThreshold()");
367
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
120 PetscCheck(eps->which==EPS_LARGEST_MAGNITUDE || eps->which==EPS_SMALLEST_MAGNITUDE || eps->which==EPS_LARGEST_REAL || eps->which==EPS_SMALLEST_REAL || eps->which==EPS_TARGET_MAGNITUDE,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Threshold stopping test can only be used with largest/smallest/target magnitude or largest/smallest real selection of eigenvalues");
368
1/6
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
120 if (eps->which==EPS_LARGEST_REAL || eps->which==EPS_SMALLEST_REAL) PetscCheck(eps->problem_type==EPS_HEP || eps->problem_type==EPS_GHEP || eps->problem_type==EPS_BSE,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Threshold stopping test with largest/smallest real can only be used in problems that have all eigenvaues real");
369
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
120 else PetscCheck(eps->thres>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"In case of largest/smallest/target magnitude the threshold value must be positive");
370
3/6
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
120 PetscCheck(eps->which==EPS_LARGEST_MAGNITUDE || eps->which==EPS_TARGET_MAGNITUDE || !eps->threlative,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"Can only use a relative threshold with largest/target magnitude selection of eigenvalues");
371
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
120 PetscCall(PetscNew(&ctx));
372 120 ctx->thres = eps->thres;
373 120 ctx->threlative = eps->threlative;
374 120 ctx->which = eps->which;
375
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
120 PetscCall(EPSSetStoppingTestFunction(eps,EPSStoppingThreshold,ctx,PetscCtxDestroyDefault));
376 }
377
378 /* if purification is set, check that it really makes sense */
379
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (eps->purify) {
380
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
5687 if (eps->categ==EPS_CATEGORY_PRECOND || eps->categ==EPS_CATEGORY_CONTOUR) eps->purify = PETSC_FALSE;
381 else {
382
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
4342 if (!eps->isgeneralized) eps->purify = PETSC_FALSE;
383
4/4
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 7 times.
1183 else if (!eps->ishermitian && !eps->ispositive) eps->purify = PETSC_FALSE;
384 else {
385
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
924 PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&flg));
386
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
924 if (flg) eps->purify = PETSC_FALSE;
387 }
388 }
389 }
390
391 /* set tolerance if not yet set */
392
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (eps->tol==(PetscReal)PETSC_DETERMINE) eps->tol = SLEPC_DEFAULT_TOL;
393
394 /* set up sorting criterion */
395
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
7998 PetscTryTypeMethod(eps,setupsort);
396
397 /* Build balancing matrix if required */
398
1/2
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
7998 if (eps->balance!=EPS_BALANCE_USER) {
399
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(STSetBalanceMatrix(eps->st,NULL));
400
4/4
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
7998 if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
401
5/8
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
104 if (!eps->D) PetscCall(BVCreateVec(eps->V,&eps->D));
402
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
104 PetscCall(EPSBuildBalance_Krylov(eps));
403
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
104 PetscCall(STSetBalanceMatrix(eps->st,eps->D));
404 }
405 }
406
407 /* Setup ST */
408
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(STSetUp(eps->st));
409
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(EPSCheckCompatibleST(eps));
410
411 /* process deflation and initial vectors */
412
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (eps->nds<0) {
413
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
120 PetscCheck(!eps->isstructured,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Deflation space is not supported in structured eigenproblems");
414 120 k = -eps->nds;
415
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
120 PetscCall(BVInsertConstraints(eps->V,&k,eps->defl));
416
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
120 PetscCall(SlepcBasisDestroy_Private(&eps->nds,&eps->defl));
417 120 eps->nds = k;
418
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
120 PetscCall(STCheckNullSpace(eps->st,eps->V));
419 }
420
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (eps->nini<0) {
421 2112 k = -eps->nini;
422
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
2112 PetscCheck(k<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The number of initial vectors is larger than ncv");
423
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
2112 PetscCall(BVInsertVecs(eps->V,0,&k,eps->IS,PETSC_TRUE));
424
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
2112 PetscCall(SlepcBasisDestroy_Private(&eps->nini,&eps->IS));
425 2112 eps->nini = k;
426 }
427
4/4
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
7998 if (eps->twosided && eps->ninil<0) {
428 32 k = -eps->ninil;
429
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
32 PetscCheck(k<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The number of left initial vectors is larger than ncv");
430
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
32 PetscCall(BVInsertVecs(eps->W,0,&k,eps->ISL,PETSC_TRUE));
431
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
32 PetscCall(SlepcBasisDestroy_Private(&eps->ninil,&eps->ISL));
432 32 eps->ninil = k;
433 }
434
435
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(PetscLogEventEnd(EPS_SetUp,eps,0,0,0));
436 7998 eps->state = EPS_STATE_SETUP;
437
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.
7998 PetscFunctionReturn(PETSC_SUCCESS);
438 }
439
440 /*@
441 EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
442
443 Collective
444
445 Input Parameters:
446 + eps - the linear eigensolver context
447 . A - the matrix associated with the eigensystem
448 - B - the second matrix in the case of generalized eigenproblems
449
450 Notes:
451 To specify a standard eigenproblem, use `NULL` for parameter `B`.
452
453 It must be called before `EPSSetUp()`. If it is called again after `EPSSetUp()` and
454 the matrix sizes have changed then the `EPS` object is reset.
455
456 For structured eigenproblem types such as `EPS_BSE` (see `EPSSetProblemType()`), the
457 provided matrices must have been created with the corresponding helper function,
458 i.e., `MatCreateBSE()`.
459
460 Level: beginner
461
462 .seealso: [](ch:eps), `EPSSolve()`, `EPSSetUp()`, `EPSReset()`, `EPSSetProblemType()`
463 @*/
464 6723 PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
465 {
466 6723 PetscInt m,n,m0,mloc,nloc,mloc0,nmat;
467 6723 Mat mat[2];
468
469
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
6723 PetscFunctionBegin;
470
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.
6723 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
471
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.
6723 PetscValidHeaderSpecific(A,MAT_CLASSID,2);
472
4/14
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
6723 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
473
13/32
✓ 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 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
6723 PetscCheckSameComm(eps,1,A,2);
474
15/34
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ 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 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.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✗ Branch 32 not taken.
✗ Branch 33 not taken.
6723 if (B) PetscCheckSameComm(eps,1,B,3);
475
476 /* Check matrix sizes */
477
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
6723 PetscCall(MatGetSize(A,&m,&n));
478
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
6723 PetscCall(MatGetLocalSize(A,&mloc,&nloc));
479
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
6723 PetscCheck(m==n,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",m,n);
480
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
6723 PetscCheck(mloc==nloc,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A does not have equal row and column sizes (%" PetscInt_FMT ", %" PetscInt_FMT ")",mloc,nloc);
481
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
6723 if (B) {
482
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1837 PetscCall(MatGetSize(B,&m0,&n));
483
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1837 PetscCall(MatGetLocalSize(B,&mloc0,&nloc));
484
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
1837 PetscCheck(m0==n,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",m0,n);
485
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
1837 PetscCheck(mloc0==nloc,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B does not have equal row and column local sizes (%" PetscInt_FMT ", %" PetscInt_FMT ")",mloc0,nloc);
486
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
1837 PetscCheck(m==m0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Dimensions of A and B do not match (%" PetscInt_FMT ", %" PetscInt_FMT ")",m,m0);
487
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
1837 PetscCheck(mloc==mloc0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Local dimensions of A and B do not match (%" PetscInt_FMT ", %" PetscInt_FMT ")",mloc,mloc0);
488 }
489
10/12
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 8 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 6 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
6723 if (eps->state && (n!=eps->n || nloc!=eps->nloc)) PetscCall(EPSReset(eps));
490 6723 eps->nrma = 0.0;
491 6723 eps->nrmb = 0.0;
492
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
6723 if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
493 6723 mat[0] = A;
494
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
6723 if (B) {
495 1837 mat[1] = B;
496 1837 nmat = 2;
497 } else nmat = 1;
498
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
6723 PetscCall(STSetMatrices(eps->st,nmat,mat));
499 6723 eps->state = EPS_STATE_INITIAL;
500
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.
6723 PetscFunctionReturn(PETSC_SUCCESS);
501 }
502
503 /*@
504 EPSGetOperators - Gets the matrices associated with the eigensystem.
505
506 Collective
507
508 Input Parameter:
509 . eps - the linear eigensolver context
510
511 Output Parameters:
512 + A - the matrix associated with the eigensystem
513 - B - the second matrix in the case of generalized eigenproblems
514
515 Note:
516 Does not increase the reference count of the matrices, so you should not destroy them.
517
518 Level: intermediate
519
520 .seealso: [](ch:eps), `EPSSolve()`, `EPSGetST()`, `STGetMatrix()`, `STSetMatrices()`
521 @*/
522 9084 PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)
523 {
524 9084 ST st;
525 9084 PetscInt k;
526
527
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
9084 PetscFunctionBegin;
528
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.
9084 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
529
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
9084 PetscCall(EPSGetST(eps,&st));
530
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
9084 PetscCall(STGetNumMatrices(st,&k));
531
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
9084 if (A) {
532
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
4690 if (k<1) *A = NULL;
533
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4690 else PetscCall(STGetMatrix(st,0,A));
534 }
535
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
9084 if (B) {
536
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
5806 if (k<2) *B = NULL;
537
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5481 else PetscCall(STGetMatrix(st,1,B));
538 }
539
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.
2239 PetscFunctionReturn(PETSC_SUCCESS);
540 }
541
542 /*@
543 EPSSetDeflationSpace - Specify a basis of vectors that constitute the deflation
544 space.
545
546 Collective
547
548 Input Parameters:
549 + eps - the linear eigensolver context
550 . n - number of vectors
551 - v - set of basis vectors of the deflation space
552
553 Notes:
554 When a deflation space is given, the eigensolver seeks the eigensolution
555 in the restriction of the problem to the orthogonal complement of this
556 space. This can be used for instance in the case that an invariant
557 subspace is known beforehand (such as the nullspace of the matrix).
558
559 These vectors do not persist from one `EPSSolve()` call to the other, so the
560 deflation space should be set every time.
561
562 The vectors do not need to be mutually orthonormal, since they are explicitly
563 orthonormalized internally.
564
565 Level: intermediate
566
567 .seealso: [](ch:eps), `EPSSetInitialSpace()`
568 @*/
569 135 PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec v[])
570 {
571
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
135 PetscFunctionBegin;
572
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.
135 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
573
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.
135 PetscValidLogicalCollectiveInt(eps,n,2);
574
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
135 PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
575
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
135 if (n>0) {
576
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.
131 PetscAssertPointer(v,3);
577
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.
131 PetscValidHeaderSpecific(*v,VEC_CLASSID,3);
578 }
579
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
135 PetscCall(SlepcBasisReference_Private(n,v,&eps->nds,&eps->defl));
580
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 7 times.
135 if (n>0) eps->state = EPS_STATE_INITIAL;
581
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.
33 PetscFunctionReturn(PETSC_SUCCESS);
582 }
583
584 /*@
585 EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
586 space, that is, the subspace from which the solver starts to iterate.
587
588 Collective
589
590 Input Parameters:
591 + eps - the linear eigensolver context
592 . n - number of vectors
593 - is - set of basis vectors of the initial space
594
595 Notes:
596 Some solvers such as `EPSKRYLOVSCHUR` start to iterate on a single vector
597 (initial vector). In that case, only `is[0]` is taken into account and the
598 other vectors are ignored. But other solvers such as `EPSSUBSPACE` are
599 able to make use of the whole initial subspace as an initial guess.
600
601 These vectors do not persist from one `EPSSolve()` call to the other, so the
602 initial space should be set every time.
603
604 The vectors do not need to be mutually orthonormal, since they are explicitly
605 orthonormalized internally.
606
607 Common usage of this function is when the user can provide a rough approximation
608 of the wanted eigenspace. Then, convergence may be faster.
609
610 Level: intermediate
611
612 .seealso: [](ch:eps), `EPSSetLeftInitialSpace()`, `EPSSetDeflationSpace()`
613 @*/
614 2131 PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec is[])
615 {
616
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2131 PetscFunctionBegin;
617
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.
2131 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
618
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.
2131 PetscValidLogicalCollectiveInt(eps,n,2);
619
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
2131 PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
620
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
2131 if (n>0) {
621
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.
2123 PetscAssertPointer(is,3);
622
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.
2123 PetscValidHeaderSpecific(*is,VEC_CLASSID,3);
623 }
624
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
2131 PetscCall(SlepcBasisReference_Private(n,is,&eps->nini,&eps->IS));
625
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 7 times.
2131 if (n>0) eps->state = EPS_STATE_INITIAL;
626
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.
501 PetscFunctionReturn(PETSC_SUCCESS);
627 }
628
629 /*@
630 EPSSetLeftInitialSpace - Specify a basis of vectors that constitute the left
631 initial space, used by two-sided solvers to start the left subspace.
632
633 Collective
634
635 Input Parameters:
636 + eps - the linear eigensolver context
637 . n - number of vectors
638 - isl - set of basis vectors of the left initial space
639
640 Notes:
641 Left initial vectors are used to initiate the left search space in two-sided
642 eigensolvers. Users should pass here an approximation of the left eigenspace,
643 if available.
644
645 The same comments in `EPSSetInitialSpace()` are applicable here.
646
647 Level: intermediate
648
649 .seealso: [](ch:eps), `EPSSetInitialSpace()`, `EPSSetTwoSided()`
650 @*/
651 32 PetscErrorCode EPSSetLeftInitialSpace(EPS eps,PetscInt n,Vec isl[])
652 {
653
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
32 PetscFunctionBegin;
654
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.
32 PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
655
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.
32 PetscValidLogicalCollectiveInt(eps,n,2);
656
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
32 PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
657
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
32 if (n>0) {
658
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.
32 PetscAssertPointer(isl,3);
659
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.
32 PetscValidHeaderSpecific(*isl,VEC_CLASSID,3);
660 }
661
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
32 PetscCall(SlepcBasisReference_Private(n,isl,&eps->ninil,&eps->ISL));
662
1/2
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
32 if (n>0) eps->state = EPS_STATE_INITIAL;
663
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);
664 }
665
666 /*
667 EPSSetDimensions_Default - Set reasonable values for ncv, mpd if not set
668 by the user. This is called at setup.
669 */
670 5741 PetscErrorCode EPSSetDimensions_Default(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
671 {
672 5741 PetscBool krylov;
673
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
5741 PetscInt nev2, n = eps->problem_type==EPS_BSE? eps->n/2: eps->n;
674
675
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
5741 PetscFunctionBegin;
676
4/4
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
5741 if (*nev==0 && eps->stop!=EPS_STOP_THRESHOLD) *nev = 1;
677
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
5741 nev2 = eps->problem_type==EPS_BSE? (*nev+1)/2: *nev;
678
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
5741 if (*ncv!=PETSC_DETERMINE) { /* ncv set */
679
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3133 PetscCall(PetscObjectTypeCompareAny((PetscObject)eps,&krylov,EPSKRYLOVSCHUR,EPSARNOLDI,EPSLANCZOS,""));
680
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
3133 if (krylov) {
681
4/8
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
2993 PetscCheck(*ncv>=nev2+1 || (*ncv==nev2 && *ncv==n),PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev+1");
682 } else {
683
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
140 PetscCheck(*ncv>=nev2,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev");
684 }
685
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
2608 } else if (*mpd!=PETSC_DETERMINE) { /* mpd set */
686 48 *ncv = PetscMin(n,nev2+(*mpd));
687 } else { /* neither set: defaults depend on nev being small or large */
688
3/4
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
2560 if (nev2<500) *ncv = PetscMin(n,PetscMax(2*(nev2),nev2+15));
689 else {
690 *mpd = 500;
691 *ncv = PetscMin(n,nev2+(*mpd));
692 }
693 }
694
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
5741 if (*mpd==PETSC_DETERMINE) *mpd = *ncv;
695
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.
1360 PetscFunctionReturn(PETSC_SUCCESS);
696 }
697
698 /*@
699 EPSAllocateSolution - Allocate memory storage for common variables such
700 as eigenvalues and eigenvectors.
701
702 Collective
703
704 Input Parameters:
705 + eps - the linear eigensolver context
706 - extra - number of additional positions, used for methods that require a
707 working basis slightly larger than `ncv`
708
709 Developer Note:
710 This is `SLEPC_EXTERN` because it may be required by user plugin `EPS`
711 implementations.
712
713 Level: developer
714
715 .seealso: [](ch:eps), `EPSSetUp()`, `EPSSetDimensions()`
716 @*/
717 7998 PetscErrorCode EPSAllocateSolution(EPS eps,PetscInt extra)
718 {
719 7998 PetscInt oldsize,requested;
720 7998 PetscRandom rand;
721 7998 Vec t;
722
723
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
7998 PetscFunctionBegin;
724 7998 requested = eps->ncv + extra;
725
726 /* oldsize is zero if this is the first time setup is called */
727
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7998 PetscCall(BVGetSizes(eps->V,NULL,NULL,&oldsize));
728
729 /* allocate space for eigenvalues and friends */
730
3/4
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 8 times.
7998 if (requested != oldsize || !eps->eigr) {
731
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5994 PetscCall(PetscFree4(eps->eigr,eps->eigi,eps->errest,eps->perm));
732
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5994 PetscCall(PetscMalloc4(requested,&eps->eigr,requested,&eps->eigi,requested,&eps->errest,requested,&eps->perm));
733 }
734
735 /* workspace for the case of arbitrary selection */
736
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (eps->arbitrary) {
737
1/8
✗ Branch 0 not taken.
✓ Branch 1 taken 8 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.
66 if (eps->rr) PetscCall(PetscFree2(eps->rr,eps->ri));
738
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
66 PetscCall(PetscMalloc2(requested,&eps->rr,requested,&eps->ri));
739 }
740
741 /* allocate V */
742
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
7998 if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
743
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (!oldsize) {
744
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
5883 if (!((PetscObject)eps->V)->type_name) PetscCall(BVSetType(eps->V,BVMAT));
745
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5883 PetscCall(STMatCreateVecsEmpty(eps->st,&t,NULL));
746
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5883 PetscCall(BVSetSizesFromVec(eps->V,t,requested));
747
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5883 PetscCall(VecDestroy(&t));
748
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
2115 } else PetscCall(BVResize(eps->V,requested,PETSC_FALSE));
749
750 /* allocate W */
751
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
7998 if (eps->twosided) {
752
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
169 PetscCall(BVGetRandomContext(eps->V,&rand)); /* make sure the random context is available when duplicating */
753
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
169 PetscCall(BVDestroy(&eps->W));
754
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
169 PetscCall(BVDuplicate(eps->V,&eps->W));
755 }
756
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.
1926 PetscFunctionReturn(PETSC_SUCCESS);
757 }
758
759 /*@
760 EPSReallocateSolution - Reallocate memory storage for common variables such
761 as the eigenvalues and the basis vectors.
762
763 Collective
764
765 Input Parameters:
766 + eps - the linear eigensolver context
767 - newsize - new size
768
769 Developer Notes:
770 This is `SLEPC_EXTERN` because it may be required by user plugin `EPS`
771 implementations.
772
773 This is called during the iteration in case the threshold stopping test has
774 been selected.
775
776 Level: developer
777
778 .seealso: [](ch:eps), `EPSAllocateSolution()`, `EPSSetThreshold()`
779 @*/
780 59 PetscErrorCode EPSReallocateSolution(EPS eps,PetscInt newsize)
781 {
782 59 PetscInt oldsize,*nperm;
783 59 PetscReal *nerrest;
784 59 PetscScalar *neigr,*neigi;
785
786
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
59 PetscFunctionBegin;
787
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
59 PetscCall(BVGetSizes(eps->V,NULL,NULL,&oldsize));
788
2/14
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
59 if (oldsize>=newsize) PetscFunctionReturn(PETSC_SUCCESS);
789
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
59 PetscCall(PetscInfo(eps,"Reallocating basis vectors to %" PetscInt_FMT " columns\n",newsize));
790
791 /* reallocate eigenvalues */
792
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
59 PetscCall(PetscMalloc4(newsize,&neigr,newsize,&neigi,newsize,&nerrest,newsize,&nperm));
793
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
59 PetscCall(PetscArraycpy(neigr,eps->eigr,oldsize));
794
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
59 PetscCall(PetscArraycpy(neigi,eps->eigi,oldsize));
795
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
59 PetscCall(PetscArraycpy(nerrest,eps->errest,oldsize));
796
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
59 PetscCall(PetscArraycpy(nperm,eps->perm,oldsize));
797
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
59 PetscCall(PetscFree4(eps->eigr,eps->eigi,eps->errest,eps->perm));
798 59 eps->eigr = neigr;
799 59 eps->eigi = neigi;
800 59 eps->errest = nerrest;
801 59 eps->perm = nperm;
802 /* reallocate V,W */
803
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
59 PetscCall(BVResize(eps->V,newsize,PETSC_TRUE));
804
1/8
✗ Branch 0 not taken.
✓ Branch 1 taken 8 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.
59 if (eps->twosided) PetscCall(BVResize(eps->W,newsize,PETSC_TRUE));
805
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.
15 PetscFunctionReturn(PETSC_SUCCESS);
806 }
807