GCC Code Coverage Report


Directory: ./
File: src/sys/classes/st/impls/sinvert/sinvert.c
Date: 2025-10-03 04:28:47
Exec Total Coverage
Lines: 146 146 100.0%
Functions: 7 7 100.0%
Branches: 373 544 68.6%

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 Implements the shift-and-invert technique for eigenvalue problems
12 */
13
14 #include <slepc/private/stimpl.h>
15
16 /*
17 Special STApply() for the BSE structured matrix
18
19 H = [ R C; -C^H -R^T ].
20
21 Assumes that H is a MATNEST and x,y are VECNEST.
22
23 The bottom part of the input vector x2 is computed as
24 either conj(x1) or -conj(x1), where the sign is given by
25 s in the context SlepcMatStruct.
26 */
27 3822 static PetscErrorCode STApply_Sinvert_BSE(ST st,Vec x,Vec y)
28 {
29 3822 Vec x1,x2;
30 3822 PetscContainer container;
31 3822 SlepcMatStruct mctx;
32
33
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
3822 PetscFunctionBegin;
34
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3822 PetscCall(PetscObjectQuery((PetscObject)st->A[0],"SlepcMatStruct",(PetscObject*)&container));
35
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3822 PetscCall(PetscContainerGetPointer(container,(void**)&mctx));
36
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3822 PetscCall(VecNestGetSubVec(x,0,&x1));
37
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3822 PetscCall(VecNestGetSubVec(x,1,&x2));
38
39 /* x2 = +/-conj(x1) */
40
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3822 PetscCall(VecCopy(x1,x2));
41
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3822 PetscCall(VecConjugate(x2));
42
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.
3822 if (mctx->s==-1.0) PetscCall(VecScale(x2,-1.0));
43
44
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
3822 PetscCall(STApply_Generic(st,x,y));
45
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.
294 PetscFunctionReturn(PETSC_SUCCESS);
46 }
47
48 1739292 static PetscErrorCode STBackTransform_Sinvert(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
49 {
50 1739292 PetscInt j;
51 #if !defined(PETSC_USE_COMPLEX)
52 682538 PetscScalar t;
53 #endif
54
55
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1739292 PetscFunctionBegin;
56 #if !defined(PETSC_USE_COMPLEX)
57
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
2093199 for (j=0;j<n;j++) {
58
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
1410661 if (eigi[j] == 0) eigr[j] = 1.0 / eigr[j] + st->sigma;
59 else {
60 262988 t = eigr[j] * eigr[j] + eigi[j] * eigi[j];
61 262988 eigr[j] = eigr[j] / t + st->sigma;
62 262988 eigi[j] = - eigi[j] / t;
63 }
64 }
65 #else
66
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
3225247 for (j=0;j<n;j++) {
67 2168493 eigr[j] = 1.0 / eigr[j] + st->sigma;
68 }
69 #endif
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.
1739292 PetscFunctionReturn(PETSC_SUCCESS);
71 }
72
73 2197 static PetscErrorCode STPostSolve_Sinvert(ST st)
74 {
75
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2197 PetscFunctionBegin;
76
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
2197 if (st->matmode == ST_MATMODE_INPLACE) {
77
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.
25 if (st->nmat>1) PetscCall(MatAXPY(st->A[0],st->sigma,st->A[1],st->str));
78
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 else PetscCall(MatShift(st->A[0],st->sigma));
79 25 st->Astate[0] = ((PetscObject)st->A[0])->state;
80 25 st->state = ST_STATE_INITIAL;
81 25 st->opready = PETSC_FALSE;
82 }
83
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.
431 PetscFunctionReturn(PETSC_SUCCESS);
84 }
85
86 /*
87 Operator (sinvert):
88 Op P M
89 if nmat=1: (A-sI)^-1 A-sI NULL
90 if nmat=2: (A-sB)^-1 B A-sB B
91 */
92 1586 static PetscErrorCode STComputeOperator_Sinvert(ST st)
93 {
94
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1586 PetscFunctionBegin;
95 /* if the user did not set the shift, use the target value */
96
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1586 if (!st->sigma_set) st->sigma = st->defsigma;
97
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1586 PetscCall(PetscObjectReference((PetscObject)st->A[1]));
98
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1586 PetscCall(MatDestroy(&st->T[0]));
99 1586 st->T[0] = st->A[1];
100
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1586 PetscCall(STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[1]));
101
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1586 PetscCall(PetscObjectReference((PetscObject)st->T[1]));
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.
1586 PetscCall(MatDestroy(&st->P));
103 1586 st->P = st->T[1];
104
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1586 st->M = (st->nmat>1)? st->T[0]: NULL;
105
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1586 if (st->Psplit) { /* build custom preconditioner from the split matrices */
106
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
70 PetscCall(STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PETSC_TRUE,PETSC_TRUE,&st->Pmat));
107 }
108
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.
305 PetscFunctionReturn(PETSC_SUCCESS);
109 }
110
111 2744 static PetscErrorCode STSetUp_Sinvert(ST st)
112 {
113 2744 PetscInt k,nc,nmat=st->nmat,nsub;
114 2744 PetscScalar *coeffs=NULL;
115 2744 PetscBool islu;
116 2744 KSP *subksp;
117 2744 PC pc,pcsub;
118 2744 Mat A00;
119 2744 MatType type;
120 2744 PetscBool flg;
121 2744 char str[64];
122
123
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2744 PetscFunctionBegin;
124
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.
2744 if (nmat>1) PetscCall(STSetWorkVecs(st,1));
125 /* if the user did not set the shift, use the target value */
126
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
2744 if (!st->sigma_set) st->sigma = st->defsigma;
127
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
2744 if (nmat>2) { /* set-up matrices for polynomial eigenproblems */
128
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1124 if (st->transform) {
129 225 nc = (nmat*(nmat+1))/2;
130
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
225 PetscCall(PetscMalloc1(nc,&coeffs));
131 /* Compute coeffs */
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.
225 PetscCall(STCoeffs_Monomial(st,coeffs));
133 /* T[0] = A_n */
134 225 k = nmat-1;
135
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
225 PetscCall(PetscObjectReference((PetscObject)st->A[k]));
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.
225 PetscCall(MatDestroy(&st->T[0]));
137 225 st->T[0] = st->A[k];
138
8/10
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 8 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 2 times.
705 for (k=1;k<nmat;k++) PetscCall(STMatMAXPY_Private(st,nmat>2?st->sigma:-st->sigma,0.0,nmat-k-1,coeffs?coeffs+(k*(k+1))/2:NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[k]));
139
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.
225 PetscCall(PetscFree(coeffs));
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.
225 PetscCall(PetscObjectReference((PetscObject)st->T[nmat-1]));
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.
225 PetscCall(MatDestroy(&st->P));
142 225 st->P = st->T[nmat-1];
143
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
225 if (st->Psplit) { /* build custom preconditioner from the split matrices */
144
6/8
✓ 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.
20 PetscCall(STMatMAXPY_Private(st,st->sigma,0.0,0,coeffs?coeffs+((nmat-1)*nmat)/2:NULL,PETSC_TRUE,PETSC_TRUE,&st->Pmat));
145 }
146
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.
225 PetscCall(ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P));
147 } else {
148
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
4934 for (k=0;k<nmat;k++) {
149
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4035 PetscCall(PetscObjectReference((PetscObject)st->A[k]));
150
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4035 PetscCall(MatDestroy(&st->T[k]));
151 4035 st->T[k] = st->A[k];
152 }
153 }
154 }
155
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
2744 if (st->structured) {
156 /* ./ex55 -st_type sinvert -eps_target 0 -eps_target_magnitude
157 -st_ksp_type preonly -st_pc_type fieldsplit
158 -st_fieldsplit_0_ksp_type preonly -st_fieldsplit_0_pc_type lu
159 -st_fieldsplit_1_ksp_type preonly -st_fieldsplit_1_pc_type lu
160 -st_pc_fieldsplit_type schur -st_pc_fieldsplit_schur_fact_type full
161 -st_pc_fieldsplit_schur_precondition full */
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.
78 PetscCall(KSPGetPC(st->ksp,&pc));
163
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCLU,&islu));
164
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
78 if (islu) { /* assume PCLU means the user has not set any options */
165
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(KSPSetType(st->ksp,KSPPREONLY));
166
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(PCSetType(pc,PCFIELDSPLIT));
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.
78 PetscCall(PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR));
168
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(PCFieldSplitSetSchurFactType(pc,PC_FIELDSPLIT_SCHUR_FACT_FULL));
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.
78 PetscCall(PCFieldSplitSetSchurPre(pc,PC_FIELDSPLIT_SCHUR_PRE_FULL,NULL));
170 /* hack to set Mat type of Schur complement equal to A00, assumes default prefixes */
171
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(PetscOptionsHasName(((PetscObject)st)->options,((PetscObject)st)->prefix,"-st_fieldsplit_1_explicit_operator_mat_type",&flg));
172
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
78 if (!flg) {
173
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(MatNestGetSubMat(st->A[0],0,0,&A00));
174
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(MatGetType(A00,&type));
175
6/8
✓ 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.
78 PetscCall(PetscSNPrintf(str,sizeof(str),"-%sst_fieldsplit_1_explicit_operator_mat_type %s",((PetscObject)st)->prefix?((PetscObject)st)->prefix:"",type));
176
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(PetscOptionsInsertString(((PetscObject)st)->options,str));
177 }
178 /* set preonly+lu on block solvers */
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.
78 PetscCall(KSPSetUp(st->ksp));
180
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(PCFieldSplitGetSubKSP(pc,&nsub,&subksp));
181
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(KSPSetType(subksp[0],KSPPREONLY));
182
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(KSPGetPC(subksp[0],&pcsub));
183
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(PCSetType(pcsub,PCLU));
184
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(KSPSetType(subksp[1],KSPPREONLY));
185
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(KSPGetPC(subksp[1],&pcsub));
186
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
78 PetscCall(PCSetType(pcsub,PCLU));
187
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.
78 PetscCall(PetscFree(subksp));
188 }
189 78 st->ops->apply = STApply_Sinvert_BSE;
190 } else {
191
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.
2666 if (st->P) PetscCall(KSPSetUp(st->ksp));
192 }
193
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.
535 PetscFunctionReturn(PETSC_SUCCESS);
194 }
195
196 5106 static PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
197 {
198 5106 PetscInt nmat=PetscMax(st->nmat,2),k,nc;
199 5106 PetscScalar *coeffs=NULL;
200
201
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
5106 PetscFunctionBegin;
202
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
5106 if (st->transform) {
203
4/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
4783 if (st->matmode == ST_MATMODE_COPY && nmat>2) {
204 20 nc = (nmat*(nmat+1))/2;
205
4/6
✓ Branch 0 taken 2 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(PetscMalloc1(nc,&coeffs));
206 /* Compute coeffs */
207
3/6
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
20 PetscCall(STCoeffs_Monomial(st,coeffs));
208 }
209
13/14
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 10 times.
✓ Branch 8 taken 10 times.
✓ Branch 9 taken 10 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✓ Branch 13 taken 2 times.
9626 for (k=1;k<nmat;k++) PetscCall(STMatMAXPY_Private(st,nmat>2?newshift:-newshift,nmat>2?st->sigma:-st->sigma,nmat-k-1,coeffs?coeffs+(k*(k+1))/2:NULL,PETSC_FALSE,PETSC_FALSE,&st->T[k]));
210
9/12
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 10 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 8 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
4783 if (st->matmode == ST_MATMODE_COPY && nmat>2) PetscCall(PetscFree(coeffs));
211
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
4783 if (st->P!=st->T[nmat-1]) {
212
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
148 PetscCall(PetscObjectReference((PetscObject)st->T[nmat-1]));
213
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
148 PetscCall(MatDestroy(&st->P));
214 148 st->P = st->T[nmat-1];
215 }
216
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
4783 if (st->Psplit) { /* build custom preconditioner from the split matrices */
217
11/12
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 8 times.
✓ Branch 5 taken 10 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 10 times.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
350 PetscCall(STMatMAXPY_Private(st,nmat>2?newshift:-newshift,nmat>2?st->sigma:-st->sigma,0,coeffs?coeffs+((nmat-1)*nmat)/2:NULL,PETSC_FALSE,PETSC_TRUE,&st->Pmat));
218 }
219 }
220
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
5106 if (st->P) {
221
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.
8867 PetscCall(ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P));
222
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5086 PetscCall(KSPSetUp(st->ksp));
223 }
224
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.
1029 PetscFunctionReturn(PETSC_SUCCESS);
225 }
226
227 2622 SLEPC_EXTERN PetscErrorCode STCreate_Sinvert(ST st)
228 {
229
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2622 PetscFunctionBegin;
230 2622 st->usesksp = PETSC_TRUE;
231
232 2622 st->ops->apply = STApply_Generic;
233 2622 st->ops->applytrans = STApplyTranspose_Generic;
234 2622 st->ops->applyhermtrans = STApplyHermitianTranspose_Generic;
235 2622 st->ops->backtransform = STBackTransform_Sinvert;
236 2622 st->ops->setshift = STSetShift_Sinvert;
237 2622 st->ops->getbilinearform = STGetBilinearForm_Default;
238 2622 st->ops->setup = STSetUp_Sinvert;
239 2622 st->ops->computeoperator = STComputeOperator_Sinvert;
240 2622 st->ops->postsolve = STPostSolve_Sinvert;
241 2622 st->ops->checknullspace = STCheckNullSpace_Default;
242 2622 st->ops->setdefaultksp = STSetDefaultKSP_Default;
243
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.
2622 PetscFunctionReturn(PETSC_SUCCESS);
244 }
245