Line | Branch | Exec | Source |
---|---|---|---|
1 | /* | ||
2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
5 | |||
6 | This file is part of SLEPc. | ||
7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
9 | */ | ||
10 | /* | ||
11 | SLEPc polynomial eigensolver: "toar" | ||
12 | |||
13 | Method: TOAR | ||
14 | |||
15 | Algorithm: | ||
16 | |||
17 | Two-Level Orthogonal Arnoldi. | ||
18 | |||
19 | References: | ||
20 | |||
21 | [1] Y. Su, J. Zhang and Z. Bai, "A compact Arnoldi algorithm for | ||
22 | polynomial eigenvalue problems", talk presented at RANMEP 2008. | ||
23 | |||
24 | [2] C. Campos and J.E. Roman, "Parallel Krylov solvers for the | ||
25 | polynomial eigenvalue problem in SLEPc", SIAM J. Sci. Comput. | ||
26 | 38(5):S385-S411, 2016. | ||
27 | |||
28 | [3] D. Lu, Y. Su and Z. Bai, "Stability analysis of the two-level | ||
29 | orthogonal Arnoldi procedure", SIAM J. Matrix Anal. App. | ||
30 | 37(1):195-214, 2016. | ||
31 | */ | ||
32 | |||
33 | #include <slepc/private/pepimpl.h> /*I "slepcpep.h" I*/ | ||
34 | #include "../src/pep/impls/krylov/pepkrylov.h" | ||
35 | #include <slepcblaslapack.h> | ||
36 | |||
37 | static PetscBool cited = PETSC_FALSE; | ||
38 | static const char citation[] = | ||
39 | "@Article{slepc-pep,\n" | ||
40 | " author = \"C. Campos and J. E. Roman\",\n" | ||
41 | " title = \"Parallel {Krylov} solvers for the polynomial eigenvalue problem in {SLEPc}\",\n" | ||
42 | " journal = \"{SIAM} J. Sci. Comput.\",\n" | ||
43 | " volume = \"38\",\n" | ||
44 | " number = \"5\",\n" | ||
45 | " pages = \"S385--S411\",\n" | ||
46 | " year = \"2016,\"\n" | ||
47 | " doi = \"https://doi.org/10.1137/15M1022458\"\n" | ||
48 | "}\n"; | ||
49 | |||
50 | 974 | static PetscErrorCode PEPSetUp_TOAR(PEP pep) | |
51 | { | ||
52 | 974 | PEP_TOAR *ctx = (PEP_TOAR*)pep->data; | |
53 | 974 | PetscBool sinv,flg; | |
54 | 974 | PetscInt i; | |
55 | |||
56 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
974 | PetscFunctionBegin; |
57 |
6/10✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
974 | PEPCheckShiftSinvert(pep); |
58 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd)); |
59 |
3/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
974 | PetscCheck(ctx->lock || pep->mpd>=pep->ncv,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant"); |
60 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
974 | if (pep->max_it==PETSC_DETERMINE) pep->max_it = PetscMax(100,2*(pep->nmat-1)*pep->n/pep->ncv); |
61 |
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.
|
974 | if (!pep->which) PetscCall(PEPSetWhichEigenpairs_Default(pep)); |
62 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
974 | PetscCheck(pep->which!=PEP_ALL,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues"); |
63 |
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.
|
974 | if (pep->problem_type!=PEP_GENERAL) PetscCall(PetscInfo(pep,"Problem type ignored, performing a non-symmetric linearization\n")); |
64 | |||
65 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
974 | if (!ctx->keep) ctx->keep = 0.5; |
66 | |||
67 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(PEPAllocateSolution(pep,pep->nmat-1)); |
68 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(PEPSetWorkVecs(pep,3)); |
69 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(DSSetType(pep->ds,DSNHEP)); |
70 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(DSSetExtraRow(pep->ds,PETSC_TRUE)); |
71 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(DSAllocate(pep->ds,pep->ncv+1)); |
72 | |||
73 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(PEPBasisCoefficients(pep,pep->pbc)); |
74 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(STGetTransform(pep->st,&flg)); |
75 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
974 | if (!flg) { |
76 |
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.
|
864 | PetscCall(PetscFree(pep->solvematcoeffs)); |
77 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
864 | PetscCall(PetscMalloc1(pep->nmat,&pep->solvematcoeffs)); |
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.
|
864 | PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv)); |
79 |
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.
|
864 | if (sinv) PetscCall(PEPEvaluateBasis(pep,pep->target,0,pep->solvematcoeffs,NULL)); |
80 | else { | ||
81 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
955 | for (i=0;i<pep->nmat-1;i++) pep->solvematcoeffs[i] = 0.0; |
82 | 305 | pep->solvematcoeffs[pep->nmat-1] = 1.0; | |
83 | } | ||
84 | } | ||
85 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(BVDestroy(&ctx->V)); |
86 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(BVCreateTensor(pep->V,pep->nmat-1,&ctx->V)); |
87 |
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.
|
188 | PetscFunctionReturn(PETSC_SUCCESS); |
88 | } | ||
89 | |||
90 | /* | ||
91 | Extend the TOAR basis by applying the matrix operator | ||
92 | over a vector which is decomposed in the TOAR way | ||
93 | Input: | ||
94 | - pbc: array containing the polynomial basis coefficients | ||
95 | - S,V: define the latest Arnoldi vector (nv vectors in V) | ||
96 | Output: | ||
97 | - t: new vector extending the TOAR basis | ||
98 | - r: temporary coefficients to compute the TOAR coefficients | ||
99 | for the new Arnoldi vector | ||
100 | Workspace: t_ (two vectors) | ||
101 | */ | ||
102 | 94142 | static PetscErrorCode PEPTOARExtendBasis(PEP pep,PetscBool sinvert,PetscScalar sigma,PetscScalar *S,PetscInt ls,PetscInt nv,BV V,Vec t,PetscScalar *r,PetscInt lr,Vec *t_) | |
103 | { | ||
104 | 94142 | PetscInt nmat=pep->nmat,deg=nmat-1,k,j,off=0,lss; | |
105 | 94142 | Vec v=t_[0],ve=t_[1],q=t_[2]; | |
106 | 94142 | PetscScalar alpha=1.0,*ss,a; | |
107 | 94142 | PetscReal *ca=pep->pbc,*cb=pep->pbc+nmat,*cg=pep->pbc+2*nmat; | |
108 | 94142 | PetscBool flg; | |
109 | |||
110 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
94142 | PetscFunctionBegin; |
111 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(BVSetActiveColumns(pep->V,0,nv)); |
112 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(STGetTransform(pep->st,&flg)); |
113 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
94142 | if (sinvert) { |
114 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1397545 | for (j=0;j<nv;j++) { |
115 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
1332343 | if (deg>1) r[lr+j] = S[j]/ca[0]; |
116 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1332343 | if (deg>2) r[2*lr+j] = (S[ls+j]+(sigma-cb[1])*r[lr+j])/ca[1]; |
117 | } | ||
118 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
319487 | for (k=2;k<deg-1;k++) { |
119 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6772146 | for (j=0;j<nv;j++) r[(k+1)*lr+j] = (S[k*ls+j]+(sigma-cb[k])*r[k*lr+j]-cg[k]*r[(k-1)*lr+j])/ca[k]; |
120 | } | ||
121 | 1397545 | k = deg-1; | |
122 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1397545 | for (j=0;j<nv;j++) r[j] = (S[k*ls+j]+(sigma-cb[k])*r[k*lr+j]-cg[k]*r[(k-1)*lr+j])/ca[k]; |
123 | 65202 | ss = r; lss = lr; off = 1; alpha = -1.0; a = pep->sfactor; | |
124 | } else { | ||
125 | 28940 | ss = S; lss = ls; off = 0; alpha = -ca[deg-1]; a = 1.0; | |
126 | } | ||
127 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(BVMultVec(V,1.0,0.0,v,ss+off*lss)); |
128 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
94142 | if (PetscUnlikely(pep->Dr)) { /* balancing */ |
129 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
14983 | PetscCall(VecPointwiseMult(v,v,pep->Dr)); |
130 | } | ||
131 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(STMatMult(pep->st,off,v,q)); |
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.
|
94142 | PetscCall(VecScale(q,a)); |
133 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
382000 | for (j=1+off;j<deg+off-1;j++) { |
134 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
287858 | PetscCall(BVMultVec(V,1.0,0.0,v,ss+j*lss)); |
135 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
287858 | if (PetscUnlikely(pep->Dr)) PetscCall(VecPointwiseMult(v,v,pep->Dr)); |
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.
|
287858 | PetscCall(STMatMult(pep->st,j,v,t)); |
137 | 287858 | a *= pep->sfactor; | |
138 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
287858 | PetscCall(VecAXPY(q,a,t)); |
139 | } | ||
140 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
94142 | if (sinvert) { |
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.
|
65202 | PetscCall(BVMultVec(V,1.0,0.0,v,ss)); |
142 |
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.
|
65202 | if (PetscUnlikely(pep->Dr)) PetscCall(VecPointwiseMult(v,v,pep->Dr)); |
143 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
65202 | PetscCall(STMatMult(pep->st,deg,v,t)); |
144 | 65202 | a *= pep->sfactor; | |
145 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
65202 | PetscCall(VecAXPY(q,a,t)); |
146 | } else { | ||
147 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28940 | PetscCall(BVMultVec(V,1.0,0.0,ve,ss+(deg-1)*lss)); |
148 |
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.
|
28940 | if (PetscUnlikely(pep->Dr)) PetscCall(VecPointwiseMult(ve,ve,pep->Dr)); |
149 | 28940 | a *= pep->sfactor; | |
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.
|
28940 | PetscCall(STMatMult(pep->st,deg-1,ve,t)); |
151 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28940 | PetscCall(VecAXPY(q,a,t)); |
152 | 28940 | a *= pep->sfactor; | |
153 | } | ||
154 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
94142 | if (flg || !sinvert) alpha /= a; |
155 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(STMatSolve(pep->st,q,t)); |
156 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(VecScale(t,alpha)); |
157 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
94142 | if (!sinvert) { |
158 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28940 | PetscCall(VecAXPY(t,cg[deg-1],v)); |
159 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28940 | PetscCall(VecAXPY(t,cb[deg-1],ve)); |
160 | } | ||
161 |
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.
|
94142 | if (PetscUnlikely(pep->Dr)) PetscCall(VecPointwiseDivide(t,t,pep->Dr)); |
162 |
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.
|
18011 | PetscFunctionReturn(PETSC_SUCCESS); |
163 | } | ||
164 | |||
165 | /* | ||
166 | Compute TOAR coefficients of the blocks of the new Arnoldi vector computed | ||
167 | */ | ||
168 | 94142 | static PetscErrorCode PEPTOARCoefficients(PEP pep,PetscBool sinvert,PetscScalar sigma,PetscInt nv,PetscScalar *S,PetscInt ls,PetscScalar *r,PetscInt lr,PetscScalar *x) | |
169 | { | ||
170 | 94142 | PetscInt k,j,nmat=pep->nmat,d=nmat-1; | |
171 | 94142 | PetscReal *ca=pep->pbc,*cb=pep->pbc+nmat,*cg=pep->pbc+2*nmat; | |
172 | 94142 | PetscScalar t=1.0,tp=0.0,tt; | |
173 | |||
174 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
94142 | PetscFunctionBegin; |
175 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
94142 | if (sinvert) { |
176 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
415642 | for (k=1;k<d;k++) { |
177 | 350440 | tt = t; | |
178 | 350440 | t = ((sigma-cb[k-1])*t-cg[k-1]*tp)/ca[k-1]; /* k-th basis polynomial */ | |
179 | 350440 | tp = tt; | |
180 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9239420 | for (j=0;j<=nv;j++) r[k*lr+j] += t*x[j]; |
181 | } | ||
182 | } else { | ||
183 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
503064 | for (j=0;j<=nv;j++) r[j] = (cb[0]-sigma)*S[j]+ca[0]*S[ls+j]; |
184 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
31560 | for (k=1;k<d-1;k++) { |
185 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
48300 | for (j=0;j<=nv;j++) r[k*lr+j] = (cb[k]-sigma)*S[k*ls+j]+ca[k]*S[(k+1)*ls+j]+cg[k]*S[(k-1)*ls+j]; |
186 | } | ||
187 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
28940 | if (sigma!=0.0) for (j=0;j<=nv;j++) r[(d-1)*lr+j] -= sigma*S[(d-1)*ls+j]; |
188 | } | ||
189 |
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.
|
94142 | PetscFunctionReturn(PETSC_SUCCESS); |
190 | } | ||
191 | |||
192 | /* | ||
193 | Compute a run of Arnoldi iterations dim(work)=ld | ||
194 | */ | ||
195 | 8719 | static PetscErrorCode PEPTOARrun(PEP pep,PetscScalar sigma,Mat A,PetscInt k,PetscInt *M,PetscReal *beta,PetscBool *breakdown,Vec *t_) | |
196 | { | ||
197 | 8719 | PEP_TOAR *ctx = (PEP_TOAR*)pep->data; | |
198 | 8719 | PetscInt j,m=*M,deg=pep->nmat-1,ld; | |
199 | 8719 | PetscInt ldh,lds,nqt,l; | |
200 | 8719 | Vec t; | |
201 | 8719 | PetscReal norm=0.0; | |
202 | 8719 | PetscBool flg,sinvert=PETSC_FALSE,lindep; | |
203 | 8719 | PetscScalar *H,*x,*S; | |
204 | 8719 | Mat MS; | |
205 | |||
206 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
8719 | PetscFunctionBegin; |
207 | 8719 | *beta = 0.0; | |
208 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(MatDenseGetArray(A,&H)); |
209 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(MatDenseGetLDA(A,&ldh)); |
210 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(BVTensorGetFactors(ctx->V,NULL,&MS)); |
211 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(MatDenseGetArray(MS,&S)); |
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.
|
8719 | PetscCall(BVGetSizes(pep->V,NULL,NULL,&ld)); |
213 | 8719 | lds = ld*deg; | |
214 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(BVGetActiveColumns(pep->V,&l,&nqt)); |
215 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(STGetTransform(pep->st,&flg)); |
216 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
8719 | if (!flg) { |
217 | /* spectral transformation handled by the solver */ | ||
218 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
7809 | PetscCall(PetscObjectTypeCompareAny((PetscObject)pep->st,&flg,STSINVERT,STSHIFT,"")); |
219 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
7809 | PetscCheck(flg,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"ST type not supported for TOAR without transforming matrices"); |
220 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
7809 | PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinvert)); |
221 | } | ||
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.
|
8719 | PetscCall(BVSetActiveColumns(ctx->V,0,m)); |
223 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
102841 | for (j=k;j<m;j++) { |
224 | /* apply operator */ | ||
225 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(BVGetColumn(pep->V,nqt,&t)); |
226 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(PEPTOARExtendBasis(pep,sinvert,sigma,S+j*lds,ld,nqt,pep->V,t,S+(j+1)*lds,ld,t_)); |
227 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(BVRestoreColumn(pep->V,nqt,&t)); |
228 | |||
229 | /* orthogonalize */ | ||
230 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
94142 | if (sinvert) x = S+(j+1)*lds; |
231 | 28940 | else x = S+(deg-1)*ld+(j+1)*lds; | |
232 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(BVOrthogonalizeColumn(pep->V,nqt,x,&norm,&lindep)); |
233 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
94142 | if (!lindep) { |
234 | 86095 | x[nqt] = norm; | |
235 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
86095 | PetscCall(BVScaleColumn(pep->V,nqt,1.0/norm)); |
236 | 86095 | nqt++; | |
237 | } | ||
238 | |||
239 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(PEPTOARCoefficients(pep,sinvert,sigma,nqt-1,S+j*lds,ld,S+(j+1)*lds,ld,x)); |
240 | |||
241 | /* level-2 orthogonalization */ | ||
242 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94142 | PetscCall(BVOrthogonalizeColumn(ctx->V,j+1,H+j*ldh,&norm,breakdown)); |
243 | 94142 | H[j+1+ldh*j] = norm; | |
244 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
94142 | if (PetscUnlikely(*breakdown)) { |
245 | 20 | *M = j+1; | |
246 | 20 | break; | |
247 | } | ||
248 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94122 | PetscCall(BVScaleColumn(ctx->V,j+1,1.0/norm)); |
249 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
94122 | PetscCall(BVSetActiveColumns(pep->V,l,nqt)); |
250 | } | ||
251 | 8719 | *beta = norm; | |
252 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(BVSetActiveColumns(ctx->V,0,*M)); |
253 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(MatDenseRestoreArray(MS,&S)); |
254 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(BVTensorRestoreFactors(ctx->V,NULL,&MS)); |
255 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(MatDenseRestoreArray(A,&H)); |
256 |
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.
|
1662 | PetscFunctionReturn(PETSC_SUCCESS); |
257 | } | ||
258 | |||
259 | /* | ||
260 | Computes T_j = phi_idx(T). In T_j and T_p are phi_{idx-1}(T) | ||
261 | and phi_{idx-2}(T) respectively or null if idx=0,1. | ||
262 | Tp and Tj are input/output arguments | ||
263 | */ | ||
264 | 50 | static PetscErrorCode PEPEvaluateBasisM(PEP pep,PetscInt k,PetscScalar *T,PetscInt ldt,PetscInt idx,PetscScalar **Tp,PetscScalar **Tj) | |
265 | { | ||
266 | 50 | PetscInt i; | |
267 | 50 | PetscReal *ca,*cb,*cg; | |
268 | 50 | PetscScalar *pt,g,a; | |
269 | 50 | PetscBLASInt k_,ldt_; | |
270 | |||
271 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
50 | PetscFunctionBegin; |
272 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
50 | if (idx==0) { |
273 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscArrayzero(*Tj,k*k)); |
274 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscArrayzero(*Tp,k*k)); |
275 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
145 | for (i=0;i<k;i++) (*Tj)[i+i*k] = 1.0; |
276 | } else { | ||
277 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscBLASIntCast(ldt,&ldt_)); |
278 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscBLASIntCast(k,&k_)); |
279 | 40 | ca = pep->pbc; cb = pep->pbc+pep->nmat; cg = pep->pbc+2*pep->nmat; | |
280 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
580 | for (i=0;i<k;i++) T[i*ldt+i] -= cb[idx-1]; |
281 | 40 | a = 1/ca[idx-1]; | |
282 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
40 | g = (idx==1)?0.0:-cg[idx-1]/ca[idx-1]; |
283 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
40 | PetscCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&a,T,&ldt_,*Tj,&k_,&g,*Tp,&k_)); |
284 | 40 | pt = *Tj; *Tj = *Tp; *Tp = pt; | |
285 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
580 | for (i=0;i<k;i++) T[i*ldt+i] += cb[idx-1]; |
286 | } | ||
287 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
288 | } | ||
289 | |||
290 | 140 | static PetscErrorCode PEPExtractInvariantPair(PEP pep,PetscScalar sigma,PetscInt sr,PetscInt k,PetscScalar *S,PetscInt ld,PetscInt deg,Mat HH) | |
291 | { | ||
292 | 140 | PetscInt i,j,jj,ldh,lds,ldt,d=pep->nmat-1,idxcpy=0; | |
293 | 140 | PetscScalar *H,*At,*Bt,*Hj,*Hp,*T,sone=1.0,g,a,*pM,*work; | |
294 | 140 | PetscBLASInt k_,sr_,lds_,ldh_,info,*p,lwork,ldt_; | |
295 | 140 | PetscBool transf=PETSC_FALSE,flg; | |
296 | 140 | PetscReal norm,maxnrm,*rwork; | |
297 | 140 | BV *R,Y; | |
298 | 140 | Mat M,*A; | |
299 | |||
300 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
140 | PetscFunctionBegin; |
301 |
2/14✓ Branch 0 taken 8 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
140 | if (k==0) PetscFunctionReturn(PETSC_SUCCESS); |
302 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(MatDenseGetArray(HH,&H)); |
303 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(MatDenseGetLDA(HH,&ldh)); |
304 | 140 | lds = deg*ld; | |
305 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(PetscCalloc6(k,&p,sr*k,&At,k*k,&Bt,k*k,&Hj,k*k,&Hp,sr*k,&work)); |
306 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(PetscBLASIntCast(sr,&sr_)); |
307 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(PetscBLASIntCast(k,&k_)); |
308 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(PetscBLASIntCast(lds,&lds_)); |
309 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(PetscBLASIntCast(ldh,&ldh_)); |
310 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(STGetTransform(pep->st,&flg)); |
311 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
140 | if (!flg) { |
312 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
110 | PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&flg)); |
313 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
110 | if (flg || sigma!=0.0) transf=PETSC_TRUE; |
314 | } | ||
315 | 70 | if (transf) { | |
316 |
4/6✓ Branch 0 taken 2 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(PetscMalloc1(k*k,&T)); |
317 | 420 | ldt = k; | |
318 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
420 | for (i=0;i<k;i++) PetscCall(PetscArraycpy(T+k*i,H+i*ldh,k)); |
319 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
70 | if (flg) { |
320 |
4/6✓ Branch 0 taken 2 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(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); |
321 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
70 | PetscCallBLAS("LAPACKgetrf",LAPACKgetrf_(&k_,&k_,T,&k_,p,&info)); |
322 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
70 | SlepcCheckLapackInfo("getrf",info); |
323 |
4/6✓ Branch 0 taken 2 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(PetscBLASIntCast(sr*k,&lwork)); |
324 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
70 | PetscCallBLAS("LAPACKgetri",LAPACKgetri_(&k_,T,&k_,p,work,&lwork,&info)); |
325 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
70 | SlepcCheckLapackInfo("getri",info); |
326 |
4/6✓ Branch 0 taken 2 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(PetscFPTrapPop()); |
327 | } | ||
328 |
3/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
420 | if (sigma!=0.0) for (i=0;i<k;i++) T[i+k*i] += sigma; |
329 | } else { | ||
330 | 70 | T = H; ldt = ldh; | |
331 | } | ||
332 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(PetscBLASIntCast(ldt,&ldt_)); |
333 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
140 | switch (pep->extract) { |
334 | case PEP_EXTRACT_NONE: | ||
335 | break; | ||
336 | 110 | case PEP_EXTRACT_NORM: | |
337 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
110 | if (pep->basis == PEP_BASIS_MONOMIAL) { |
338 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
110 | PetscCall(PetscBLASIntCast(ldt,&ldt_)); |
339 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
110 | PetscCall(PetscMalloc1(k,&rwork)); |
340 | 110 | norm = LAPACKlange_("F",&k_,&k_,T,&ldt_,rwork); | |
341 |
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.
|
110 | PetscCall(PetscFree(rwork)); |
342 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
110 | if (norm>1.0) idxcpy = d-1; |
343 | } else { | ||
344 | ✗ | PetscCall(PetscBLASIntCast(ldt,&ldt_)); | |
345 | ✗ | PetscCall(PetscMalloc1(k,&rwork)); | |
346 | maxnrm = 0.0; | ||
347 | ✗ | for (i=0;i<pep->nmat-1;i++) { | |
348 | ✗ | PetscCall(PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj)); | |
349 | ✗ | norm = LAPACKlange_("F",&k_,&k_,Hj,&k_,rwork); | |
350 | ✗ | if (norm > maxnrm) { | |
351 | ✗ | idxcpy = i; | |
352 | ✗ | maxnrm = norm; | |
353 | } | ||
354 | } | ||
355 | ✗ | PetscCall(PetscFree(rwork)); | |
356 | } | ||
357 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
110 | if (idxcpy>0) { |
358 | /* copy block idxcpy of S to the first one */ | ||
359 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
790 | for (j=0;j<k;j++) PetscCall(PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr)); |
360 | } | ||
361 | break; | ||
362 | 10 | case PEP_EXTRACT_RESIDUAL: | |
363 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(STGetTransform(pep->st,&flg)); |
364 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | if (flg) { |
365 | ✗ | PetscCall(PetscMalloc1(pep->nmat,&A)); | |
366 | ✗ | for (i=0;i<pep->nmat;i++) PetscCall(STGetMatrixTransformed(pep->st,i,A+i)); | |
367 | 10 | } else A = pep->A; | |
368 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscMalloc1(pep->nmat-1,&R)); |
369 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
30 | for (i=0;i<pep->nmat-1;i++) PetscCall(BVDuplicateResize(pep->V,k,R+i)); |
370 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(BVDuplicateResize(pep->V,sr,&Y)); |
371 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,sr,k,NULL,&M)); |
372 | 10 | g = 0.0; a = 1.0; | |
373 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(BVSetActiveColumns(pep->V,0,sr)); |
374 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
40 | for (j=0;j<pep->nmat;j++) { |
375 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(BVMatMult(pep->V,A[j],Y)); |
376 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj)); |
377 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
90 | for (i=0;i<pep->nmat-1;i++) { |
378 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
60 | PetscCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&k_,&k_,&a,S+i*ld,&lds_,Hj,&k_,&g,At,&sr_)); |
379 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(MatDenseGetArray(M,&pM)); |
380 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
870 | for (jj=0;jj<k;jj++) PetscCall(PetscArraycpy(pM+jj*sr,At+jj*sr,sr)); |
381 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(MatDenseRestoreArray(M,&pM)); |
382 |
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.
|
90 | PetscCall(BVMult(R[i],1.0,(i==0)?0.0:1.0,Y,M)); |
383 | } | ||
384 | } | ||
385 | |||
386 | /* frobenius norm */ | ||
387 | maxnrm = 0.0; | ||
388 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
30 | for (i=0;i<pep->nmat-1;i++) { |
389 |
4/6✓ Branch 0 taken 2 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(BVNorm(R[i],NORM_FROBENIUS,&norm)); |
390 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
20 | if (maxnrm > norm) { |
391 | ✗ | maxnrm = norm; | |
392 | ✗ | idxcpy = i; | |
393 | } | ||
394 | } | ||
395 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | if (idxcpy>0) { |
396 | /* copy block idxcpy of S to the first one */ | ||
397 | ✗ | for (j=0;j<k;j++) PetscCall(PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr)); | |
398 | } | ||
399 |
1/10✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
10 | if (flg) PetscCall(PetscFree(A)); |
400 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
30 | for (i=0;i<pep->nmat-1;i++) PetscCall(BVDestroy(&R[i])); |
401 |
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.
|
10 | PetscCall(PetscFree(R)); |
402 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(BVDestroy(&Y)); |
403 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatDestroy(&M)); |
404 | break; | ||
405 | case PEP_EXTRACT_STRUCTURED: | ||
406 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
145 | for (j=0;j<k;j++) Bt[j+j*k] = 1.0; |
407 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
145 | for (j=0;j<sr;j++) { |
408 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1960 | for (i=0;i<k;i++) At[j*k+i] = PetscConj(S[i*lds+j]); |
409 | } | ||
410 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PEPEvaluateBasisM(pep,k,T,ldt,0,&Hp,&Hj)); |
411 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20 | for (i=1;i<deg;i++) { |
412 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj)); |
413 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
10 | PetscCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&sr_,&k_,&sone,Hj,&k_,S+i*ld,&lds_,&sone,At,&k_)); |
414 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
10 | PetscCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&k_,&k_,&sone,Hj,&k_,Hj,&k_,&sone,Bt,&k_)); |
415 | } | ||
416 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); |
417 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
10 | PetscCallBLAS("LAPACKgesv",LAPACKgesv_(&k_,&sr_,Bt,&k_,p,At,&k_,&info)); |
418 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscFPTrapPop()); |
419 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
10 | SlepcCheckLapackInfo("gesv",info); |
420 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
145 | for (j=0;j<sr;j++) { |
421 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1960 | for (i=0;i<k;i++) S[i*lds+j] = PetscConj(At[j*k+i]); |
422 | } | ||
423 | break; | ||
424 | } | ||
425 |
7/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
140 | if (transf) PetscCall(PetscFree(T)); |
426 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(PetscFree6(p,At,Bt,Hj,Hp,work)); |
427 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(MatDenseRestoreArray(HH,&H)); |
428 |
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.
|
28 | PetscFunctionReturn(PETSC_SUCCESS); |
429 | } | ||
430 | |||
431 | 974 | static PetscErrorCode PEPSolve_TOAR(PEP pep) | |
432 | { | ||
433 | 974 | PEP_TOAR *ctx = (PEP_TOAR*)pep->data; | |
434 | 974 | PetscInt i,j,k,l,nv=0,ld,lds,nq=0,nconv=0; | |
435 | 974 | PetscInt nmat=pep->nmat,deg=nmat-1; | |
436 | 974 | PetscScalar *S,sigma; | |
437 | 974 | PetscReal beta; | |
438 | 974 | PetscBool breakdown=PETSC_FALSE,flg,falselock=PETSC_FALSE,sinv=PETSC_FALSE; | |
439 | 974 | Mat H,MS,MQ; | |
440 | |||
441 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
974 | PetscFunctionBegin; |
442 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(PetscCitationsRegister(citation,&cited)); |
443 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
974 | if (ctx->lock) { |
444 | /* undocumented option to use a cheaper locking instead of the true locking */ | ||
445 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
964 | PetscCall(PetscOptionsGetBool(NULL,NULL,"-pep_toar_falselocking",&falselock,NULL)); |
446 | } | ||
447 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(STGetShift(pep->st,&sigma)); |
448 | |||
449 | /* update polynomial basis coefficients */ | ||
450 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(STGetTransform(pep->st,&flg)); |
451 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
974 | if (pep->sfactor!=1.0) { |
452 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
144 | for (i=0;i<nmat;i++) { |
453 | 108 | pep->pbc[nmat+i] /= pep->sfactor; | |
454 | 108 | pep->pbc[2*nmat+i] /= pep->sfactor*pep->sfactor; | |
455 | } | ||
456 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
36 | if (!flg) { |
457 | 36 | pep->target /= pep->sfactor; | |
458 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(RGPushScale(pep->rg,1.0/pep->sfactor)); |
459 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(STScaleShift(pep->st,1.0/pep->sfactor)); |
460 | 36 | sigma /= pep->sfactor; | |
461 | } else { | ||
462 | ✗ | PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv)); | |
463 | ✗ | pep->target = sinv?pep->target*pep->sfactor:pep->target/pep->sfactor; | |
464 | ✗ | PetscCall(RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor)); | |
465 | ✗ | PetscCall(STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor)); | |
466 | } | ||
467 | } | ||
468 | |||
469 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
974 | if (flg) sigma = 0.0; |
470 | |||
471 | /* clean projected matrix (including the extra-arrow) */ | ||
472 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(DSSetDimensions(pep->ds,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_DETERMINE)); |
473 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(DSGetMat(pep->ds,DS_MAT_A,&H)); |
474 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(MatZeroEntries(H)); |
475 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(DSRestoreMat(pep->ds,DS_MAT_A,&H)); |
476 | |||
477 | /* Get the starting Arnoldi vector */ | ||
478 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(BVTensorBuildFirstColumn(ctx->V,pep->nini)); |
479 | |||
480 | /* restart loop */ | ||
481 | 974 | l = 0; | |
482 | 974 | while (pep->reason == PEP_CONVERGED_ITERATING) { | |
483 | 8719 | pep->its++; | |
484 | |||
485 | /* compute an nv-step Lanczos factorization */ | ||
486 | 8719 | nv = PetscMax(PetscMin(nconv+pep->mpd,pep->ncv),nv); | |
487 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(DSGetMat(pep->ds,DS_MAT_A,&H)); |
488 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(PEPTOARrun(pep,sigma,H,pep->nconv+l,&nv,&beta,&breakdown,pep->work)); |
489 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(DSRestoreMat(pep->ds,DS_MAT_A,&H)); |
490 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(DSSetDimensions(pep->ds,nv,pep->nconv,pep->nconv+l)); |
491 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(DSSetState(pep->ds,l?DS_STATE_RAW:DS_STATE_INTERMEDIATE)); |
492 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(BVSetActiveColumns(ctx->V,pep->nconv,nv)); |
493 | |||
494 | /* solve projected problem */ | ||
495 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(DSSolve(pep->ds,pep->eigr,pep->eigi)); |
496 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL)); |
497 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(DSUpdateExtraRow(pep->ds)); |
498 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(DSSynchronize(pep->ds,pep->eigr,pep->eigi)); |
499 | |||
500 | /* check convergence */ | ||
501 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k)); |
502 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall((*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx)); |
503 | |||
504 | /* update l */ | ||
505 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
8719 | if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0; |
506 | else { | ||
507 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
7745 | l = (nv==k)?0:PetscMax(1,(PetscInt)((nv-k)*ctx->keep)); |
508 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
7745 | PetscCall(DSGetTruncateSize(pep->ds,k,nv,&l)); |
509 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
7745 | if (!breakdown) { |
510 | /* prepare the Rayleigh quotient for restart */ | ||
511 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
7745 | PetscCall(DSTruncate(pep->ds,k+l,PETSC_FALSE)); |
512 | } | ||
513 | } | ||
514 | 8719 | nconv = k; | |
515 |
5/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 10 times.
✗ Branch 5 not taken.
|
8719 | if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */ |
516 |
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.
|
8719 | if (l) PetscCall(PetscInfo(pep,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l)); |
517 | |||
518 | /* update S */ | ||
519 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(DSGetMat(pep->ds,DS_MAT_Q,&MQ)); |
520 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(BVMultInPlace(ctx->V,MQ,pep->nconv,k+l)); |
521 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(DSRestoreMat(pep->ds,DS_MAT_Q,&MQ)); |
522 | |||
523 | /* copy last column of S */ | ||
524 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(BVCopyColumn(ctx->V,nv,k+l)); |
525 | |||
526 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
8719 | if (PetscUnlikely(breakdown && pep->reason == PEP_CONVERGED_ITERATING)) { |
527 | /* stop if breakdown */ | ||
528 | ✗ | PetscCall(PetscInfo(pep,"Breakdown TOAR method (it=%" PetscInt_FMT " norm=%g)\n",pep->its,(double)beta)); | |
529 | ✗ | pep->reason = PEP_DIVERGED_BREAKDOWN; | |
530 | } | ||
531 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
8719 | if (pep->reason != PEP_CONVERGED_ITERATING) l--; |
532 | /* truncate S */ | ||
533 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8719 | PetscCall(BVGetActiveColumns(pep->V,NULL,&nq)); |
534 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
8719 | if (k+l+deg<=nq) { |
535 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8176 | PetscCall(BVSetActiveColumns(ctx->V,pep->nconv,k+l+1)); |
536 |
7/10✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
8176 | if (!falselock && ctx->lock) PetscCall(BVTensorCompress(ctx->V,k-pep->nconv)); |
537 |
4/6✓ Branch 0 taken 2 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 | else PetscCall(BVTensorCompress(ctx->V,0)); |
538 | } | ||
539 | 8719 | pep->nconv = k; | |
540 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
9693 | PetscCall(PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv)); |
541 | } | ||
542 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
974 | if (pep->nconv>0) { |
543 | /* {V*S_nconv^i}_{i=0}^{d-1} has rank nconv instead of nconv+d-1. Force zeros in each S_nconv^i block */ | ||
544 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(BVSetActiveColumns(ctx->V,0,pep->nconv)); |
545 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(BVGetActiveColumns(pep->V,NULL,&nq)); |
546 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(BVSetActiveColumns(pep->V,0,nq)); |
547 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
974 | if (nq>pep->nconv) { |
548 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
25 | PetscCall(BVTensorCompress(ctx->V,pep->nconv)); |
549 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
25 | PetscCall(BVSetActiveColumns(pep->V,0,pep->nconv)); |
550 | 25 | nq = pep->nconv; | |
551 | } | ||
552 | |||
553 | /* perform Newton refinement if required */ | ||
554 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
974 | if (pep->refine==PEP_REFINE_MULTIPLE && pep->rits>0) { |
555 | /* extract invariant pair */ | ||
556 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(BVTensorGetFactors(ctx->V,NULL,&MS)); |
557 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(MatDenseGetArray(MS,&S)); |
558 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(DSGetMat(pep->ds,DS_MAT_A,&H)); |
559 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(BVGetSizes(pep->V,NULL,NULL,&ld)); |
560 | 140 | lds = deg*ld; | |
561 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(PEPExtractInvariantPair(pep,sigma,nq,pep->nconv,S,ld,deg,H)); |
562 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(DSRestoreMat(pep->ds,DS_MAT_A,&H)); |
563 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(DSSetDimensions(pep->ds,pep->nconv,0,0)); |
564 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(DSSetState(pep->ds,DS_STATE_RAW)); |
565 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(PEPNewtonRefinement_TOAR(pep,sigma,&pep->rits,NULL,pep->nconv,S,lds)); |
566 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(DSSolve(pep->ds,pep->eigr,pep->eigi)); |
567 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL)); |
568 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(DSSynchronize(pep->ds,pep->eigr,pep->eigi)); |
569 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(DSGetMat(pep->ds,DS_MAT_Q,&MQ)); |
570 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(BVMultInPlace(ctx->V,MQ,0,pep->nconv)); |
571 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(DSRestoreMat(pep->ds,DS_MAT_Q,&MQ)); |
572 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(MatDenseRestoreArray(MS,&S)); |
573 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
140 | PetscCall(BVTensorRestoreFactors(ctx->V,NULL,&MS)); |
574 | } | ||
575 | } | ||
576 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(STGetTransform(pep->st,&flg)); |
577 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
974 | if (pep->refine!=PEP_REFINE_MULTIPLE || pep->rits==0) { |
578 |
7/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
834 | if (!flg) PetscTryTypeMethod(pep,backtransform); |
579 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
834 | if (pep->sfactor!=1.0) { |
580 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
220 | for (j=0;j<pep->nconv;j++) { |
581 | 184 | pep->eigr[j] *= pep->sfactor; | |
582 | 184 | pep->eigi[j] *= pep->sfactor; | |
583 | } | ||
584 | /* restore original values */ | ||
585 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
144 | for (i=0;i<pep->nmat;i++) { |
586 | 108 | pep->pbc[pep->nmat+i] *= pep->sfactor; | |
587 | 108 | pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor; | |
588 | } | ||
589 | } | ||
590 | } | ||
591 | /* restore original values */ | ||
592 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
974 | if (!flg) { |
593 | 864 | pep->target *= pep->sfactor; | |
594 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
864 | PetscCall(STScaleShift(pep->st,pep->sfactor)); |
595 | } else { | ||
596 |
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.
|
110 | PetscCall(STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor)); |
597 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
110 | pep->target = (sinv)?pep->target/pep->sfactor:pep->target*pep->sfactor; |
598 | } | ||
599 |
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.
|
974 | if (pep->sfactor!=1.0) PetscCall(RGPopScale(pep->rg)); |
600 | |||
601 | /* change the state to raw so that DSVectors() computes eigenvectors from scratch */ | ||
602 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(DSSetDimensions(pep->ds,pep->nconv,0,0)); |
603 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
974 | PetscCall(DSSetState(pep->ds,DS_STATE_RAW)); |
604 |
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.
|
188 | PetscFunctionReturn(PETSC_SUCCESS); |
605 | } | ||
606 | |||
607 | 20 | static PetscErrorCode PEPTOARSetRestart_TOAR(PEP pep,PetscReal keep) | |
608 | { | ||
609 | 20 | PEP_TOAR *ctx = (PEP_TOAR*)pep->data; | |
610 | |||
611 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBegin; |
612 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
20 | if (keep==(PetscReal)PETSC_DEFAULT || keep==(PetscReal)PETSC_DECIDE) ctx->keep = 0.5; |
613 | else { | ||
614 |
2/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
20 | PetscCheck(keep>=0.1 && keep<=0.9,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]"); |
615 | 20 | ctx->keep = keep; | |
616 | } | ||
617 |
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.
|
4 | PetscFunctionReturn(PETSC_SUCCESS); |
618 | } | ||
619 | |||
620 | /*@ | ||
621 | PEPTOARSetRestart - Sets the restart parameter for the TOAR | ||
622 | method, in particular the proportion of basis vectors that must be kept | ||
623 | after restart. | ||
624 | |||
625 | Logically Collective | ||
626 | |||
627 | Input Parameters: | ||
628 | + pep - the eigenproblem solver context | ||
629 | - keep - the number of vectors to be kept at restart | ||
630 | |||
631 | Options Database Key: | ||
632 | . -pep_toar_restart - Sets the restart parameter | ||
633 | |||
634 | Notes: | ||
635 | Allowed values are in the range [0.1,0.9]. The default is 0.5. | ||
636 | |||
637 | Level: advanced | ||
638 | |||
639 | .seealso: PEPTOARGetRestart() | ||
640 | @*/ | ||
641 | 20 | PetscErrorCode PEPTOARSetRestart(PEP pep,PetscReal keep) | |
642 | { | ||
643 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBegin; |
644 |
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.
|
20 | PetscValidHeaderSpecific(pep,PEP_CLASSID,1); |
645 |
29/66✗ Branch 0 not taken.
✓ 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 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✓ Branch 26 taken 2 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
✓ Branch 33 taken 2 times.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✓ Branch 36 taken 2 times.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✓ Branch 39 taken 2 times.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 42 not taken.
✓ Branch 43 taken 2 times.
✗ Branch 44 not taken.
✗ Branch 45 not taken.
✗ Branch 46 not taken.
✓ Branch 47 taken 2 times.
✓ Branch 48 taken 2 times.
✗ Branch 49 not taken.
✗ Branch 50 not taken.
✓ Branch 51 taken 2 times.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✓ Branch 54 taken 2 times.
✗ Branch 55 not taken.
✗ Branch 56 not taken.
✓ Branch 57 taken 2 times.
✗ Branch 58 not taken.
✗ Branch 59 not taken.
✓ Branch 60 taken 2 times.
✗ Branch 61 not taken.
✗ Branch 62 not taken.
✓ Branch 63 taken 2 times.
✗ Branch 64 not taken.
✗ Branch 65 not taken.
|
20 | PetscValidLogicalCollectiveReal(pep,keep,2); |
646 |
8/14✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
|
20 | PetscTryMethod(pep,"PEPTOARSetRestart_C",(PEP,PetscReal),(pep,keep)); |
647 |
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.
|
20 | PetscFunctionReturn(PETSC_SUCCESS); |
648 | } | ||
649 | |||
650 | 10 | static PetscErrorCode PEPTOARGetRestart_TOAR(PEP pep,PetscReal *keep) | |
651 | { | ||
652 | 10 | PEP_TOAR *ctx = (PEP_TOAR*)pep->data; | |
653 | |||
654 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
655 | 10 | *keep = ctx->keep; | |
656 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
657 | } | ||
658 | |||
659 | /*@ | ||
660 | PEPTOARGetRestart - Gets the restart parameter used in the TOAR method. | ||
661 | |||
662 | Not Collective | ||
663 | |||
664 | Input Parameter: | ||
665 | . pep - the eigenproblem solver context | ||
666 | |||
667 | Output Parameter: | ||
668 | . keep - the restart parameter | ||
669 | |||
670 | Level: advanced | ||
671 | |||
672 | .seealso: PEPTOARSetRestart() | ||
673 | @*/ | ||
674 | 10 | PetscErrorCode PEPTOARGetRestart(PEP pep,PetscReal *keep) | |
675 | { | ||
676 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
677 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(pep,PEP_CLASSID,1); |
678 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | PetscAssertPointer(keep,2); |
679 |
9/16✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
|
10 | PetscUseMethod(pep,"PEPTOARGetRestart_C",(PEP,PetscReal*),(pep,keep)); |
680 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
681 | } | ||
682 | |||
683 | 10 | static PetscErrorCode PEPTOARSetLocking_TOAR(PEP pep,PetscBool lock) | |
684 | { | ||
685 | 10 | PEP_TOAR *ctx = (PEP_TOAR*)pep->data; | |
686 | |||
687 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
688 | 10 | ctx->lock = lock; | |
689 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
690 | } | ||
691 | |||
692 | /*@ | ||
693 | PEPTOARSetLocking - Choose between locking and non-locking variants of | ||
694 | the TOAR method. | ||
695 | |||
696 | Logically Collective | ||
697 | |||
698 | Input Parameters: | ||
699 | + pep - the eigenproblem solver context | ||
700 | - lock - true if the locking variant must be selected | ||
701 | |||
702 | Options Database Key: | ||
703 | . -pep_toar_locking - Sets the locking flag | ||
704 | |||
705 | Notes: | ||
706 | The default is to lock converged eigenpairs when the method restarts. | ||
707 | This behaviour can be changed so that all directions are kept in the | ||
708 | working subspace even if already converged to working accuracy (the | ||
709 | non-locking variant). | ||
710 | |||
711 | Level: advanced | ||
712 | |||
713 | .seealso: PEPTOARGetLocking() | ||
714 | @*/ | ||
715 | 10 | PetscErrorCode PEPTOARSetLocking(PEP pep,PetscBool lock) | |
716 | { | ||
717 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
718 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(pep,PEP_CLASSID,1); |
719 |
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.
|
10 | PetscValidLogicalCollectiveBool(pep,lock,2); |
720 |
8/14✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
|
10 | PetscTryMethod(pep,"PEPTOARSetLocking_C",(PEP,PetscBool),(pep,lock)); |
721 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
722 | } | ||
723 | |||
724 | 10 | static PetscErrorCode PEPTOARGetLocking_TOAR(PEP pep,PetscBool *lock) | |
725 | { | ||
726 | 10 | PEP_TOAR *ctx = (PEP_TOAR*)pep->data; | |
727 | |||
728 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
729 | 10 | *lock = ctx->lock; | |
730 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
731 | } | ||
732 | |||
733 | /*@ | ||
734 | PEPTOARGetLocking - Gets the locking flag used in the TOAR method. | ||
735 | |||
736 | Not Collective | ||
737 | |||
738 | Input Parameter: | ||
739 | . pep - the eigenproblem solver context | ||
740 | |||
741 | Output Parameter: | ||
742 | . lock - the locking flag | ||
743 | |||
744 | Level: advanced | ||
745 | |||
746 | .seealso: PEPTOARSetLocking() | ||
747 | @*/ | ||
748 | 10 | PetscErrorCode PEPTOARGetLocking(PEP pep,PetscBool *lock) | |
749 | { | ||
750 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
751 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(pep,PEP_CLASSID,1); |
752 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | PetscAssertPointer(lock,2); |
753 |
9/16✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
|
10 | PetscUseMethod(pep,"PEPTOARGetLocking_C",(PEP,PetscBool*),(pep,lock)); |
754 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
755 | } | ||
756 | |||
757 | 930 | static PetscErrorCode PEPSetFromOptions_TOAR(PEP pep,PetscOptionItems PetscOptionsObject) | |
758 | { | ||
759 | 930 | PetscBool flg,lock; | |
760 | 930 | PetscReal keep; | |
761 | |||
762 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
930 | PetscFunctionBegin; |
763 |
1/12✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
930 | PetscOptionsHeadBegin(PetscOptionsObject,"PEP TOAR Options"); |
764 | |||
765 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
930 | PetscCall(PetscOptionsReal("-pep_toar_restart","Proportion of vectors kept after restart","PEPTOARSetRestart",0.5,&keep,&flg)); |
766 |
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.
|
930 | if (flg) PetscCall(PEPTOARSetRestart(pep,keep)); |
767 | |||
768 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
930 | PetscCall(PetscOptionsBool("-pep_toar_locking","Choose between locking and non-locking variants","PEPTOARSetLocking",PETSC_FALSE,&lock,&flg)); |
769 |
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.
|
930 | if (flg) PetscCall(PEPTOARSetLocking(pep,lock)); |
770 | |||
771 |
1/14✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
930 | PetscOptionsHeadEnd(); |
772 |
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.
|
175 | PetscFunctionReturn(PETSC_SUCCESS); |
773 | } | ||
774 | |||
775 | 15 | static PetscErrorCode PEPView_TOAR(PEP pep,PetscViewer viewer) | |
776 | { | ||
777 | 15 | PEP_TOAR *ctx = (PEP_TOAR*)pep->data; | |
778 | 15 | PetscBool isascii; | |
779 | |||
780 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
15 | PetscFunctionBegin; |
781 |
4/6✓ Branch 0 taken 2 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 | PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); |
782 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
15 | if (isascii) { |
783 |
4/6✓ Branch 0 taken 2 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 | PetscCall(PetscViewerASCIIPrintf(viewer," %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep))); |
784 |
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.
|
15 | PetscCall(PetscViewerASCIIPrintf(viewer," using the %slocking variant\n",ctx->lock?"":"non-")); |
785 | } | ||
786 |
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.
|
3 | PetscFunctionReturn(PETSC_SUCCESS); |
787 | } | ||
788 | |||
789 | 958 | static PetscErrorCode PEPDestroy_TOAR(PEP pep) | |
790 | { | ||
791 | 958 | PEP_TOAR *ctx = (PEP_TOAR*)pep->data; | |
792 | |||
793 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
958 | PetscFunctionBegin; |
794 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
958 | PetscCall(BVDestroy(&ctx->V)); |
795 |
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.
|
958 | PetscCall(PetscFree(pep->data)); |
796 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
958 | PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",NULL)); |
797 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
958 | PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",NULL)); |
798 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
958 | PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",NULL)); |
799 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
958 | PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",NULL)); |
800 |
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.
|
182 | PetscFunctionReturn(PETSC_SUCCESS); |
801 | } | ||
802 | |||
803 | 958 | SLEPC_EXTERN PetscErrorCode PEPCreate_TOAR(PEP pep) | |
804 | { | ||
805 | 958 | PEP_TOAR *ctx; | |
806 | |||
807 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
958 | PetscFunctionBegin; |
808 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
958 | PetscCall(PetscNew(&ctx)); |
809 | 958 | pep->data = (void*)ctx; | |
810 | |||
811 | 958 | pep->lineariz = PETSC_TRUE; | |
812 | 958 | ctx->lock = PETSC_TRUE; | |
813 | |||
814 | 958 | pep->ops->solve = PEPSolve_TOAR; | |
815 | 958 | pep->ops->setup = PEPSetUp_TOAR; | |
816 | 958 | pep->ops->setfromoptions = PEPSetFromOptions_TOAR; | |
817 | 958 | pep->ops->destroy = PEPDestroy_TOAR; | |
818 | 958 | pep->ops->view = PEPView_TOAR; | |
819 | 958 | pep->ops->backtransform = PEPBackTransform_Default; | |
820 | 958 | pep->ops->computevectors = PEPComputeVectors_Default; | |
821 | 958 | pep->ops->extractvectors = PEPExtractVectors_TOAR; | |
822 | |||
823 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
958 | PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",PEPTOARSetRestart_TOAR)); |
824 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
958 | PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",PEPTOARGetRestart_TOAR)); |
825 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
958 | PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",PEPTOARSetLocking_TOAR)); |
826 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
958 | PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",PEPTOARGetLocking_TOAR)); |
827 |
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.
|
182 | PetscFunctionReturn(PETSC_SUCCESS); |
828 | } | ||
829 |