Line | Branch | Exec | Source |
---|---|---|---|
1 | /* | ||
2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
5 | |||
6 | This file is part of SLEPc. | ||
7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
9 | */ | ||
10 | /* | ||
11 | SLEPc eigensolver: "krylovschur" | ||
12 | |||
13 | Method: Two-sided Arnoldi with Krylov-Schur restart (for left eigenvectors) | ||
14 | |||
15 | References: | ||
16 | |||
17 | [1] I.N. Zwaan and M.E. Hochstenbach, "Krylov-Schur-type restarts | ||
18 | for the two-sided Arnoldi method", SIAM J. Matrix Anal. Appl. | ||
19 | 38(2):297-321, 2017. | ||
20 | |||
21 | */ | ||
22 | |||
23 | #include <slepc/private/epsimpl.h> | ||
24 | #include "krylovschur.h" | ||
25 | #include <slepcblaslapack.h> | ||
26 | |||
27 | 1006 | static PetscErrorCode EPSTwoSidedRQUpdate1(EPS eps,Mat M,PetscInt nv,PetscReal beta,PetscReal betat) | |
28 | { | ||
29 | 1006 | PetscScalar *T,*S,*A,*w; | |
30 | 1006 | const PetscScalar *pM; | |
31 | 1006 | Vec u; | |
32 | 1006 | PetscInt ld,ncv=eps->ncv,i,l,nnv; | |
33 | 1006 | PetscBLASInt info,n_,ncv_,*p,one=1; | |
34 | |||
35 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1006 | PetscFunctionBegin; |
36 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSGetLeadingDimension(eps->ds,&ld)); |
37 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(PetscMalloc3(nv,&p,ncv*ncv,&A,ncv,&w)); |
38 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVGetActiveColumns(eps->V,&l,&nnv)); |
39 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVSetActiveColumns(eps->V,0,nv)); |
40 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVSetActiveColumns(eps->W,0,nv)); |
41 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVGetColumn(eps->V,nv,&u)); |
42 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVDotVec(eps->W,u,w)); |
43 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVRestoreColumn(eps->V,nv,&u)); |
44 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(MatDenseGetArrayRead(M,&pM)); |
45 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(PetscArraycpy(A,pM,ncv*ncv)); |
46 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(MatDenseRestoreArrayRead(M,&pM)); |
47 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(PetscBLASIntCast(nv,&n_)); |
48 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(PetscBLASIntCast(ncv,&ncv_)); |
49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); |
50 |
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.
|
1006 | PetscCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n_,&n_,A,&ncv_,p,&info)); |
51 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
1006 | SlepcCheckLapackInfo("getrf",info); |
52 |
3/6✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(PetscLogFlops(2.0*n_*n_*n_/3.0)); |
53 |
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.
|
1006 | PetscCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n_,&one,A,&ncv_,p,w,&ncv_,&info)); |
54 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
1006 | SlepcCheckLapackInfo("getrs",info); |
55 |
3/6✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(PetscLogFlops(2.0*n_*n_-n_)); |
56 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVMultColumn(eps->V,-1.0,1.0,nv,w)); |
57 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSGetArray(eps->ds,DS_MAT_A,&S)); |
58 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20702 | for (i=0;i<nv;i++) S[(nv-1)*ld+i] += beta*w[i]; |
59 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSRestoreArray(eps->ds,DS_MAT_A,&S)); |
60 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVGetColumn(eps->W,nv,&u)); |
61 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVDotVec(eps->V,u,w)); |
62 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVRestoreColumn(eps->W,nv,&u)); |
63 |
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.
|
1006 | PetscCallBLAS("LAPACKgetrs",LAPACKgetrs_("C",&n_,&one,A,&ncv_,p,w,&ncv_,&info)); |
64 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(PetscFPTrapPop()); |
65 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVMultColumn(eps->W,-1.0,1.0,nv,w)); |
66 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSGetArray(eps->ds,DS_MAT_B,&T)); |
67 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20702 | for (i=0;i<nv;i++) T[(nv-1)*ld+i] += betat*w[i]; |
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.
|
1006 | PetscCall(DSRestoreArray(eps->ds,DS_MAT_B,&T)); |
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.
|
1006 | PetscCall(PetscFree3(p,A,w)); |
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.
|
1006 | PetscCall(BVSetActiveColumns(eps->V,l,nnv)); |
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.
|
1006 | PetscCall(BVSetActiveColumns(eps->W,l,nnv)); |
72 |
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.
|
200 | PetscFunctionReturn(PETSC_SUCCESS); |
73 | } | ||
74 | |||
75 | 847 | static PetscErrorCode EPSTwoSidedRQUpdate2(EPS eps,Mat M,PetscInt k) | |
76 | { | ||
77 | 847 | PetscScalar *Q,*pM,*w,zero=0.0,sone=1.0,*c,*A; | |
78 | 847 | PetscBLASInt n_,ncv_,ld_; | |
79 | 847 | PetscReal norm; | |
80 | 847 | PetscInt l,nv,ncv=eps->ncv,ld,i,j; | |
81 | |||
82 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
847 | PetscFunctionBegin; |
83 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(DSGetLeadingDimension(eps->ds,&ld)); |
84 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(BVGetActiveColumns(eps->V,&l,&nv)); |
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.
|
847 | PetscCall(BVSetActiveColumns(eps->V,0,nv)); |
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.
|
847 | PetscCall(BVSetActiveColumns(eps->W,0,nv)); |
87 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(PetscMalloc2(ncv*ncv,&w,ncv,&c)); |
88 | /* u = u - V*V'*u */ | ||
89 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(BVOrthogonalizeColumn(eps->V,k,c,&norm,NULL)); |
90 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(BVScaleColumn(eps->V,k,1.0/norm)); |
91 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(DSGetArray(eps->ds,DS_MAT_A,&A)); |
92 | /* H = H + V'*u*b' */ | ||
93 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9048 | for (j=l;j<k;j++) { |
94 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
93334 | for (i=0;i<k;i++) A[i+j*ld] += c[i]*A[k+j*ld]; |
95 | 8201 | A[k+j*ld] *= norm; | |
96 | } | ||
97 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(DSRestoreArray(eps->ds,DS_MAT_A,&A)); |
98 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(BVOrthogonalizeColumn(eps->W,k,c,&norm,NULL)); |
99 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(BVScaleColumn(eps->W,k,1.0/norm)); |
100 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(DSGetArray(eps->ds,DS_MAT_B,&A)); |
101 | /* H = H + V'*u*b' */ | ||
102 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9048 | for (j=l;j<k;j++) { |
103 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
93334 | for (i=0;i<k;i++) A[i+j*ld] += c[i]*A[k+j*ld]; |
104 | 8201 | A[k+j*ld] *= norm; | |
105 | } | ||
106 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(DSRestoreArray(eps->ds,DS_MAT_B,&A)); |
107 | |||
108 | /* M = Q'*M*Q */ | ||
109 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(MatDenseGetArray(M,&pM)); |
110 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(PetscBLASIntCast(ncv,&ncv_)); |
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.
|
847 | PetscCall(PetscBLASIntCast(nv,&n_)); |
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.
|
847 | PetscCall(PetscBLASIntCast(ld,&ld_)); |
113 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(DSGetArray(eps->ds,DS_MAT_Q,&Q)); |
114 |
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.
|
847 | PetscCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pM,&ncv_,Q,&ld_,&zero,w,&ncv_)); |
115 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(DSRestoreArray(eps->ds,DS_MAT_Q,&Q)); |
116 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(DSGetArray(eps->ds,DS_MAT_Z,&Q)); |
117 |
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.
|
847 | PetscCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,Q,&ld_,w,&ncv_,&zero,pM,&ncv_)); |
118 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(DSRestoreArray(eps->ds,DS_MAT_Z,&Q)); |
119 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(MatDenseRestoreArray(M,&pM)); |
120 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(PetscFree2(w,c)); |
121 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(BVSetActiveColumns(eps->V,l,nv)); |
122 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(BVSetActiveColumns(eps->W,l,nv)); |
123 |
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.
|
169 | PetscFunctionReturn(PETSC_SUCCESS); |
124 | } | ||
125 | |||
126 | 159 | PetscErrorCode EPSSolve_KrylovSchur_TwoSided(EPS eps) | |
127 | { | ||
128 | 159 | EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data; | |
129 | 159 | Mat M,U,Op,OpHT,S,T; | |
130 | 159 | PetscReal norm,norm2,beta,betat; | |
131 | 159 | PetscInt ld,l,nv,nvt,k,nconv,dsn,dsk; | |
132 | 159 | PetscBool breakdownt,breakdown,breakdownl; | |
133 | |||
134 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
159 | PetscFunctionBegin; |
135 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
159 | PetscCall(DSGetLeadingDimension(eps->ds,&ld)); |
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.
|
159 | PetscCall(EPSGetStartVector(eps,0,NULL)); |
137 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
159 | PetscCall(EPSGetLeftStartVector(eps,0,NULL)); |
138 | 159 | l = 0; | |
139 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
159 | PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,eps->ncv,eps->ncv,NULL,&M)); |
140 | |||
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.
|
159 | PetscCall(STGetOperator(eps->st,&Op)); |
142 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
159 | PetscCall(MatCreateHermitianTranspose(Op,&OpHT)); |
143 | |||
144 | /* Restart loop */ | ||
145 | 1165 | while (eps->reason == EPS_CONVERGED_ITERATING) { | |
146 | 1006 | eps->its++; | |
147 | |||
148 | /* Compute an nv-step Arnoldi factorization for Op */ | ||
149 | 1006 | nv = PetscMin(eps->nconv+eps->mpd,eps->ncv); | |
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.
|
1006 | PetscCall(DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l)); |
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.
|
1006 | PetscCall(DSGetMat(eps->ds,DS_MAT_A,&S)); |
152 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVMatArnoldi(eps->V,Op,S,eps->nconv+l,&nv,&beta,&breakdown)); |
153 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&S)); |
154 | |||
155 | /* Compute an nv-step Arnoldi factorization for Op' */ | ||
156 | 1006 | nvt = nv; | |
157 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l)); |
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.
|
1006 | PetscCall(DSGetMat(eps->ds,DS_MAT_B,&T)); |
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.
|
1006 | PetscCall(BVMatArnoldi(eps->W,OpHT,T,eps->nconv+l,&nvt,&betat,&breakdownt)); |
160 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSRestoreMat(eps->ds,DS_MAT_B,&T)); |
161 | |||
162 | /* Make sure both factorizations have the same length */ | ||
163 | 1006 | nv = PetscMin(nv,nvt); | |
164 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l)); |
165 |
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.
|
1006 | if (l==0) PetscCall(DSSetState(eps->ds,DS_STATE_INTERMEDIATE)); |
166 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | else PetscCall(DSSetState(eps->ds,DS_STATE_RAW)); |
167 |
2/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
1006 | breakdown = (breakdown || breakdownt)? PETSC_TRUE: PETSC_FALSE; |
168 | |||
169 | /* Update M, modify Rayleigh quotients S and T */ | ||
170 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVSetActiveColumns(eps->V,eps->nconv+l,nv)); |
171 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVSetActiveColumns(eps->W,eps->nconv+l,nv)); |
172 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVMatProject(eps->V,NULL,eps->W,M)); |
173 | |||
174 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(EPSTwoSidedRQUpdate1(eps,M,nv,beta,betat)); |
175 | |||
176 | /* Solve projected problem */ | ||
177 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi)); |
178 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL)); |
179 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi)); |
180 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSUpdateExtraRow(eps->ds)); |
181 | |||
182 | /* Check convergence */ | ||
183 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVNormColumn(eps->V,nv,NORM_2,&norm)); |
184 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVNormColumn(eps->W,nv,NORM_2,&norm2)); |
185 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta*norm,betat*norm2,1.0,&k)); |
186 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx)); |
187 | 1006 | nconv = k; | |
188 | |||
189 | /* Update l */ | ||
190 |
4/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
|
1006 | if (eps->reason != EPS_CONVERGED_ITERATING || breakdown || k==nv) l = 0; |
191 | else { | ||
192 | 847 | l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep)); | |
193 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(DSGetTruncateSize(eps->ds,k,nv,&l)); |
194 | } | ||
195 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
1006 | if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */ |
196 |
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.
|
1006 | if (l) PetscCall(PetscInfo(eps,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l)); |
197 | |||
198 | /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */ | ||
199 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVSetActiveColumns(eps->V,eps->nconv,nv)); |
200 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVSetActiveColumns(eps->W,eps->nconv,nv)); |
201 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSGetMat(eps->ds,DS_MAT_Q,&U)); |
202 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVMultInPlace(eps->V,U,eps->nconv,k+l)); |
203 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSRestoreMat(eps->ds,DS_MAT_Q,&U)); |
204 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSGetMat(eps->ds,DS_MAT_Z,&U)); |
205 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(BVMultInPlace(eps->W,U,eps->nconv,k+l)); |
206 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1006 | PetscCall(DSRestoreMat(eps->ds,DS_MAT_Z,&U)); |
207 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
1006 | if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) { |
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.
|
847 | PetscCall(BVCopyColumn(eps->V,nv,k+l)); |
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.
|
847 | PetscCall(BVCopyColumn(eps->W,nv,k+l)); |
210 | } | ||
211 | |||
212 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1006 | if (eps->reason == EPS_CONVERGED_ITERATING) { |
213 |
2/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
847 | if (breakdown || k==nv) { |
214 | /* Start a new Arnoldi factorization */ | ||
215 | ✗ | PetscCall(PetscInfo(eps,"Breakdown in Krylov-Schur method (it=%" PetscInt_FMT " norm=%g)\n",eps->its,(double)beta)); | |
216 | ✗ | if (k<eps->nev) { | |
217 | ✗ | PetscCall(EPSGetStartVector(eps,k,&breakdown)); | |
218 | ✗ | PetscCall(EPSGetLeftStartVector(eps,k,&breakdownl)); | |
219 | ✗ | if (breakdown || breakdownl) { | |
220 | ✗ | eps->reason = EPS_DIVERGED_BREAKDOWN; | |
221 | ✗ | PetscCall(PetscInfo(eps,"Unable to generate more start vectors\n")); | |
222 | } | ||
223 | } | ||
224 | } else { | ||
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.
|
847 | PetscCall(DSGetDimensions(eps->ds,&dsn,NULL,&dsk,NULL)); |
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.
|
847 | PetscCall(DSSetDimensions(eps->ds,dsn,k,dsk)); |
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.
|
847 | PetscCall(DSTruncate(eps->ds,k+l,PETSC_FALSE)); |
228 | } | ||
229 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
847 | PetscCall(EPSTwoSidedRQUpdate2(eps,M,k+l)); |
230 | } | ||
231 | 1006 | eps->nconv = k; | |
232 |
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.
|
1165 | PetscCall(EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv)); |
233 | } | ||
234 | |||
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.
|
159 | PetscCall(STRestoreOperator(eps->st,&Op)); |
236 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
159 | PetscCall(MatDestroy(&OpHT)); |
237 | |||
238 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
159 | PetscCall(DSTruncate(eps->ds,eps->nconv,PETSC_TRUE)); |
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.
|
159 | PetscCall(MatDestroy(&M)); |
240 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
31 | PetscFunctionReturn(PETSC_SUCCESS); |
241 | } | ||
242 |