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 | This file implements a wrapper to the BLOPEX package | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/ | ||
15 | #include "blopex.h" | ||
16 | #include <lobpcg.h> | ||
17 | #include <interpreter.h> | ||
18 | #include <multivector.h> | ||
19 | #include <temp_multivector.h> | ||
20 | |||
21 | PetscInt slepc_blopex_useconstr = -1; | ||
22 | |||
23 | typedef struct { | ||
24 | lobpcg_Tolerance tol; | ||
25 | lobpcg_BLASLAPACKFunctions blap_fn; | ||
26 | mv_InterfaceInterpreter ii; | ||
27 | ST st; | ||
28 | Vec w; | ||
29 | PetscInt bs; /* block size */ | ||
30 | } EPS_BLOPEX; | ||
31 | |||
32 | 4710 | static void Precond_FnSingleVector(void *data,void *x,void *y) | |
33 | { | ||
34 | 4710 | EPS_BLOPEX *blopex = (EPS_BLOPEX*)data; | |
35 | 4710 | MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st); | |
36 | 4710 | KSP ksp; | |
37 | |||
38 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4710 | PetscFunctionBegin; |
39 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4710 | PetscCallAbort(comm,STGetKSP(blopex->st,&ksp)); |
40 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4710 | PetscCallAbort(comm,KSPSolve(ksp,(Vec)x,(Vec)y)); |
41 |
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.
|
4710 | PetscFunctionReturnVoid(); |
42 | } | ||
43 | |||
44 | 1580 | static void Precond_FnMultiVector(void *data,void *x,void *y) | |
45 | { | ||
46 | 1580 | EPS_BLOPEX *blopex = (EPS_BLOPEX*)data; | |
47 | |||
48 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1580 | PetscFunctionBegin; |
49 | 1580 | blopex->ii.Eval(Precond_FnSingleVector,data,x,y); | |
50 |
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.
|
1580 | PetscFunctionReturnVoid(); |
51 | } | ||
52 | |||
53 | 4920 | static void OperatorASingleVector(void *data,void *x,void *y) | |
54 | { | ||
55 | 4920 | EPS_BLOPEX *blopex = (EPS_BLOPEX*)data; | |
56 | 4920 | MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st); | |
57 | 4920 | Mat A,B; | |
58 | 4920 | PetscScalar sigma; | |
59 | 4920 | PetscInt nmat; | |
60 | |||
61 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4920 | PetscFunctionBegin; |
62 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4920 | PetscCallAbort(comm,STGetNumMatrices(blopex->st,&nmat)); |
63 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4920 | PetscCallAbort(comm,STGetMatrix(blopex->st,0,&A)); |
64 |
6/8✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
4920 | if (nmat>1) PetscCallAbort(comm,STGetMatrix(blopex->st,1,&B)); |
65 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4920 | PetscCallAbort(comm,MatMult(A,(Vec)x,(Vec)y)); |
66 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
4920 | PetscCallAbort(comm,STGetShift(blopex->st,&sigma)); |
67 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
4920 | if (sigma != 0.0) { |
68 |
6/8✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1068 | if (nmat>1) PetscCallAbort(comm,MatMult(B,(Vec)x,blopex->w)); |
69 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
825 | else PetscCallAbort(comm,VecCopy((Vec)x,blopex->w)); |
70 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1068 | PetscCallAbort(comm,VecAXPY((Vec)y,-sigma,blopex->w)); |
71 | } | ||
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.
|
4920 | PetscFunctionReturnVoid(); |
73 | } | ||
74 | |||
75 | 1634 | static void OperatorAMultiVector(void *data,void *x,void *y) | |
76 | { | ||
77 | 1634 | EPS_BLOPEX *blopex = (EPS_BLOPEX*)data; | |
78 | |||
79 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1634 | PetscFunctionBegin; |
80 | 1634 | blopex->ii.Eval(OperatorASingleVector,data,x,y); | |
81 |
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.
|
1634 | PetscFunctionReturnVoid(); |
82 | } | ||
83 | |||
84 | 2796 | static void OperatorBSingleVector(void *data,void *x,void *y) | |
85 | { | ||
86 | 2796 | EPS_BLOPEX *blopex = (EPS_BLOPEX*)data; | |
87 | 2796 | MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st); | |
88 | 2796 | Mat B; | |
89 | |||
90 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2796 | PetscFunctionBegin; |
91 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2796 | PetscCallAbort(comm,STGetMatrix(blopex->st,1,&B)); |
92 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2796 | PetscCallAbort(comm,MatMult(B,(Vec)x,(Vec)y)); |
93 |
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.
|
2796 | PetscFunctionReturnVoid(); |
94 | } | ||
95 | |||
96 | 969 | static void OperatorBMultiVector(void *data,void *x,void *y) | |
97 | { | ||
98 | 969 | EPS_BLOPEX *blopex = (EPS_BLOPEX*)data; | |
99 | |||
100 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
969 | PetscFunctionBegin; |
101 | 969 | blopex->ii.Eval(OperatorBSingleVector,data,x,y); | |
102 |
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.
|
969 | PetscFunctionReturnVoid(); |
103 | } | ||
104 | |||
105 | 42 | static PetscErrorCode EPSSetDimensions_BLOPEX(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd) | |
106 | { | ||
107 | 42 | EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data; | |
108 | 42 | PetscInt k; | |
109 | |||
110 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
42 | PetscFunctionBegin; |
111 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
42 | if (*nev==0) *nev = 1; |
112 | 42 | k = ((*nev-1)/ctx->bs+1)*ctx->bs; | |
113 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
42 | if (*ncv!=PETSC_DETERMINE) { /* ncv set */ |
114 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | PetscCheck(*ncv>=k,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv is not sufficiently large"); |
115 | 36 | } else *ncv = k; | |
116 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
42 | if (*mpd==PETSC_DETERMINE) *mpd = *ncv; |
117 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6 | else PetscCall(PetscInfo(eps,"Warning: given value of mpd ignored\n")); |
118 |
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.
|
14 | PetscFunctionReturn(PETSC_SUCCESS); |
119 | } | ||
120 | |||
121 | 42 | static PetscErrorCode EPSSetUp_BLOPEX(EPS eps) | |
122 | { | ||
123 | 42 | EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data; | |
124 | 42 | PetscBool flg; | |
125 | 42 | KSP ksp; | |
126 | |||
127 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
42 | PetscFunctionBegin; |
128 |
4/10✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 6 times.
✓ Branch 5 taken 6 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 6 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
42 | EPSCheckHermitianDefinite(eps); |
129 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
42 | EPSCheckNotStructured(eps); |
130 |
4/6✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 6 times.
✗ Branch 5 not taken.
|
66 | if (!blopex->bs) blopex->bs = PetscMin(16,eps->nev?eps->nev:1); |
131 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(EPSSetDimensions_BLOPEX(eps,&eps->nev,&eps->ncv,&eps->mpd)); |
132 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
42 | if (eps->max_it==PETSC_DETERMINE) eps->max_it = PetscMax(100,2*eps->n/eps->ncv); |
133 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
42 | if (!eps->which) eps->which = EPS_SMALLEST_REAL; |
134 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
42 | PetscCheck(eps->which==EPS_SMALLEST_REAL,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only smallest real eigenvalues"); |
135 |
8/18✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 4 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 4 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 6 times.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
|
42 | EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING); |
136 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 4 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.
|
42 | EPSCheckIgnored(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_EXTRACTION); |
137 | |||
138 | 42 | blopex->st = eps->st; | |
139 | |||
140 |
3/6✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 6 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
42 | PetscCheck(eps->converged==EPSConvergedRelative || eps->converged==EPSConvergedAbsolute,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Convergence test not supported in this solver"); |
141 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
42 | if (eps->converged == EPSConvergedRelative) { |
142 | 24 | blopex->tol.absolute = 0.0; | |
143 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
30 | blopex->tol.relative = SlepcDefaultTol(eps->tol); |
144 | } else { /* EPSConvergedAbsolute */ | ||
145 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | blopex->tol.absolute = SlepcDefaultTol(eps->tol); |
146 | 18 | blopex->tol.relative = 0.0; | |
147 | } | ||
148 | |||
149 | 42 | SLEPCSetupInterpreter(&blopex->ii); | |
150 | |||
151 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(STGetKSP(eps->st,&ksp)); |
152 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg)); |
153 |
6/8✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
42 | if (!flg) PetscCall(PetscInfo(eps,"Warning: ignoring KSP, should use KSPPREONLY\n")); |
154 | |||
155 | /* allocate memory */ | ||
156 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 6 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.
|
42 | if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V)); |
157 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(PetscObjectTypeCompareAny((PetscObject)eps->V,&flg,BVVECS,BVCONTIGUOUS,"")); |
158 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
42 | if (!flg) { /* blopex only works with BVVECS or BVCONTIGUOUS */ |
159 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(BVSetType(eps->V,BVCONTIGUOUS)); |
160 | } | ||
161 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(EPSAllocateSolution(eps,0)); |
162 |
6/8✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
42 | if (!blopex->w) PetscCall(BVCreateVec(eps->V,&blopex->w)); |
163 | |||
164 | #if defined(PETSC_USE_COMPLEX) | ||
165 | 21 | blopex->blap_fn.zpotrf = PETSC_zpotrf_interface; | |
166 | 21 | blopex->blap_fn.zhegv = PETSC_zsygv_interface; | |
167 | #else | ||
168 | 21 | blopex->blap_fn.dpotrf = PETSC_dpotrf_interface; | |
169 | 21 | blopex->blap_fn.dsygv = PETSC_dsygv_interface; | |
170 | #endif | ||
171 |
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.
|
42 | PetscFunctionReturn(PETSC_SUCCESS); |
172 | } | ||
173 | |||
174 | 42 | static PetscErrorCode EPSSolve_BLOPEX(EPS eps) | |
175 | { | ||
176 | 42 | EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data; | |
177 | 42 | PetscScalar sigma,*eigr=NULL; | |
178 | 42 | PetscReal *errest=NULL; | |
179 | 42 | int i,j,info,its,nconv; | |
180 | 42 | double *residhist=NULL; | |
181 | 42 | mv_MultiVectorPtr eigenvectors,constraints; | |
182 | #if defined(PETSC_USE_COMPLEX) | ||
183 | 21 | komplex *lambda=NULL,*lambdahist=NULL; | |
184 | #else | ||
185 | 21 | double *lambda=NULL,*lambdahist=NULL; | |
186 | #endif | ||
187 | |||
188 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
42 | PetscFunctionBegin; |
189 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(STGetShift(eps->st,&sigma)); |
190 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
42 | PetscCall(PetscMalloc1(blopex->bs,&lambda)); |
191 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 6 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.
|
42 | if (eps->numbermonitors>0) PetscCall(PetscMalloc4(blopex->bs*(eps->max_it+1),&lambdahist,eps->ncv,&eigr,blopex->bs*(eps->max_it+1),&residhist,eps->ncv,&errest)); |
192 | |||
193 | /* Complete the initial basis with random vectors */ | ||
194 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
48 | for (i=0;i<eps->nini;i++) { /* in case the initial vectors were also set with VecSetRandom */ |
195 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6 | PetscCall(BVSetRandomColumn(eps->V,eps->nini)); |
196 | } | ||
197 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 6 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
246 | for (i=eps->nini;i<eps->ncv;i++) PetscCall(BVSetRandomColumn(eps->V,i)); |
198 | |||
199 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
96 | while (eps->reason == EPS_CONVERGED_ITERATING) { |
200 | |||
201 | /* Create multivector of constraints from leading columns of V */ | ||
202 |
7/10✓ Branch 0 taken 6 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 6 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
54 | PetscCall(PetscObjectComposedDataSetInt((PetscObject)eps->V,slepc_blopex_useconstr,1)); |
203 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
54 | PetscCall(BVSetActiveColumns(eps->V,0,eps->nconv)); |
204 | 54 | constraints = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->nds+eps->nconv,eps->V); | |
205 | |||
206 | /* Create multivector where eigenvectors of this run will be stored */ | ||
207 |
5/10✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
54 | PetscCall(PetscObjectComposedDataSetInt((PetscObject)eps->V,slepc_blopex_useconstr,0)); |
208 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
54 | PetscCall(BVSetActiveColumns(eps->V,eps->nconv,eps->nconv+blopex->bs)); |
209 | 54 | eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,blopex->bs,eps->V); | |
210 | |||
211 | #if defined(PETSC_USE_COMPLEX) | ||
212 |
2/2✓ Branch 0 taken 3 times.
✓ Branch 1 taken 3 times.
|
45 | info = lobpcg_solve_complex(eigenvectors,blopex,OperatorAMultiVector, |
213 | 27 | eps->isgeneralized?blopex:NULL,eps->isgeneralized?OperatorBMultiVector:NULL, | |
214 | blopex,Precond_FnMultiVector,constraints, | ||
215 | blopex->blap_fn,blopex->tol,eps->max_it,0,&its, | ||
216 |
2/2✓ Branch 0 taken 3 times.
✓ Branch 1 taken 3 times.
|
27 | lambda,lambdahist,blopex->bs,eps->errest+eps->nconv,residhist,blopex->bs); |
217 | #else | ||
218 |
2/2✓ Branch 0 taken 3 times.
✓ Branch 1 taken 3 times.
|
45 | info = lobpcg_solve_double(eigenvectors,blopex,OperatorAMultiVector, |
219 | 27 | eps->isgeneralized?blopex:NULL,eps->isgeneralized?OperatorBMultiVector:NULL, | |
220 | blopex,Precond_FnMultiVector,constraints, | ||
221 | blopex->blap_fn,blopex->tol,eps->max_it,0,&its, | ||
222 |
2/2✓ Branch 0 taken 3 times.
✓ Branch 1 taken 3 times.
|
27 | lambda,lambdahist,blopex->bs,eps->errest+eps->nconv,residhist,blopex->bs); |
223 | #endif | ||
224 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
54 | PetscCheck(info==0,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"BLOPEX failed with exit code=%d",info); |
225 | 54 | mv_MultiVectorDestroy(constraints); | |
226 | 54 | mv_MultiVectorDestroy(eigenvectors); | |
227 | |||
228 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
318 | for (j=0;j<blopex->bs;j++) { |
229 | #if defined(PETSC_USE_COMPLEX) | ||
230 | 105 | eps->eigr[eps->nconv+j] = PetscCMPLX(lambda[j].real,lambda[j].imag); | |
231 | #else | ||
232 | 105 | eps->eigr[eps->nconv+j] = lambda[j]; | |
233 | #endif | ||
234 | } | ||
235 | |||
236 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
54 | if (eps->numbermonitors>0) { |
237 | ✗ | for (i=0;i<its;i++) { | |
238 | nconv = 0; | ||
239 | ✗ | for (j=0;j<blopex->bs;j++) { | |
240 | #if defined(PETSC_USE_COMPLEX) | ||
241 | ✗ | eigr[eps->nconv+j] = PetscCMPLX(lambdahist[j+i*blopex->bs].real,lambdahist[j+i*blopex->bs].imag); | |
242 | #else | ||
243 | ✗ | eigr[eps->nconv+j] = lambdahist[j+i*blopex->bs]; | |
244 | #endif | ||
245 | ✗ | errest[eps->nconv+j] = residhist[j+i*blopex->bs]; | |
246 | ✗ | if (residhist[j+i*blopex->bs]<=eps->tol) nconv++; | |
247 | } | ||
248 | ✗ | PetscCall(EPSMonitor(eps,eps->its+i,eps->nconv+nconv,eigr,eps->eigi,errest,eps->nconv+blopex->bs)); | |
249 | } | ||
250 | } | ||
251 | |||
252 | 54 | eps->its += its; | |
253 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
54 | if (info==-1) { |
254 | ✗ | eps->reason = EPS_DIVERGED_ITS; | |
255 | ✗ | break; | |
256 | } else { | ||
257 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
264 | for (i=0;i<blopex->bs;i++) { |
258 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
210 | if (sigma != 0.0) eps->eigr[eps->nconv+i] += sigma; |
259 | } | ||
260 | 54 | eps->nconv += blopex->bs; | |
261 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
54 | if (eps->nconv>=eps->nev) eps->reason = EPS_CONVERGED_TOL; |
262 | } | ||
263 | } | ||
264 | |||
265 |
5/8✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
42 | PetscCall(PetscFree(lambda)); |
266 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 6 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.
|
42 | if (eps->numbermonitors>0) PetscCall(PetscFree4(lambdahist,eigr,residhist,errest)); |
267 |
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.
|
14 | PetscFunctionReturn(PETSC_SUCCESS); |
268 | } | ||
269 | |||
270 | 12 | static PetscErrorCode EPSBLOPEXSetBlockSize_BLOPEX(EPS eps,PetscInt bs) | |
271 | { | ||
272 | 12 | EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data; | |
273 | |||
274 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
12 | PetscFunctionBegin; |
275 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | if (bs==PETSC_DEFAULT || bs==PETSC_DECIDE) { |
276 | ✗ | ctx->bs = 0; | |
277 | ✗ | eps->state = EPS_STATE_INITIAL; | |
278 | } else { | ||
279 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
12 | PetscCheck(bs>0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block size must be >0"); |
280 | 12 | ctx->bs = bs; | |
281 | } | ||
282 |
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); |
283 | } | ||
284 | |||
285 | /*@ | ||
286 | EPSBLOPEXSetBlockSize - Sets the block size of the BLOPEX solver. | ||
287 | |||
288 | Logically Collective | ||
289 | |||
290 | Input Parameters: | ||
291 | + eps - the eigenproblem solver context | ||
292 | - bs - the block size | ||
293 | |||
294 | Options Database Key: | ||
295 | . -eps_blopex_blocksize - Sets the block size | ||
296 | |||
297 | Level: advanced | ||
298 | |||
299 | .seealso: EPSBLOPEXGetBlockSize() | ||
300 | @*/ | ||
301 | 12 | PetscErrorCode EPSBLOPEXSetBlockSize(EPS eps,PetscInt bs) | |
302 | { | ||
303 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
12 | PetscFunctionBegin; |
304 |
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.
|
12 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
305 |
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.
|
12 | PetscValidLogicalCollectiveInt(eps,bs,2); |
306 |
8/14✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 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.
|
12 | PetscTryMethod(eps,"EPSBLOPEXSetBlockSize_C",(EPS,PetscInt),(eps,bs)); |
307 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
308 | } | ||
309 | |||
310 | 6 | static PetscErrorCode EPSBLOPEXGetBlockSize_BLOPEX(EPS eps,PetscInt *bs) | |
311 | { | ||
312 | 6 | EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data; | |
313 | |||
314 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6 | PetscFunctionBegin; |
315 | 6 | *bs = ctx->bs; | |
316 |
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.
|
6 | PetscFunctionReturn(PETSC_SUCCESS); |
317 | } | ||
318 | |||
319 | /*@ | ||
320 | EPSBLOPEXGetBlockSize - Gets the block size used in the BLOPEX solver. | ||
321 | |||
322 | Not Collective | ||
323 | |||
324 | Input Parameter: | ||
325 | . eps - the eigenproblem solver context | ||
326 | |||
327 | Output Parameter: | ||
328 | . bs - the block size | ||
329 | |||
330 | Level: advanced | ||
331 | |||
332 | .seealso: EPSBLOPEXSetBlockSize() | ||
333 | @*/ | ||
334 | 6 | PetscErrorCode EPSBLOPEXGetBlockSize(EPS eps,PetscInt *bs) | |
335 | { | ||
336 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6 | PetscFunctionBegin; |
337 |
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.
|
6 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); |
338 |
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.
|
6 | PetscAssertPointer(bs,2); |
339 |
9/16✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 6 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.
|
6 | PetscUseMethod(eps,"EPSBLOPEXGetBlockSize_C",(EPS,PetscInt*),(eps,bs)); |
340 |
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.
|
6 | PetscFunctionReturn(PETSC_SUCCESS); |
341 | } | ||
342 | |||
343 | 36 | static PetscErrorCode EPSReset_BLOPEX(EPS eps) | |
344 | { | ||
345 | 36 | EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data; | |
346 | |||
347 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | PetscFunctionBegin; |
348 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(VecDestroy(&blopex->w)); |
349 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
350 | } | ||
351 | |||
352 | 36 | static PetscErrorCode EPSDestroy_BLOPEX(EPS eps) | |
353 | { | ||
354 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | PetscFunctionBegin; |
355 |
5/8✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
36 | PetscCall(PetscFree(eps->data)); |
356 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXSetBlockSize_C",NULL)); |
357 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXGetBlockSize_C",NULL)); |
358 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
359 | } | ||
360 | |||
361 | 6 | static PetscErrorCode EPSView_BLOPEX(EPS eps,PetscViewer viewer) | |
362 | { | ||
363 | 6 | EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data; | |
364 | 6 | PetscBool isascii; | |
365 | |||
366 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6 | PetscFunctionBegin; |
367 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6 | PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); |
368 |
5/8✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
6 | if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer," block size %" PetscInt_FMT "\n",ctx->bs)); |
369 |
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.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
370 | } | ||
371 | |||
372 | 36 | static PetscErrorCode EPSSetFromOptions_BLOPEX(EPS eps,PetscOptionItems PetscOptionsObject) | |
373 | { | ||
374 | 36 | PetscBool flg; | |
375 | 36 | PetscInt bs; | |
376 | |||
377 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | PetscFunctionBegin; |
378 |
1/12✗ Branch 0 not taken.
✓ Branch 1 taken 6 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.
|
36 | PetscOptionsHeadBegin(PetscOptionsObject,"EPS BLOPEX Options"); |
379 | |||
380 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscOptionsInt("-eps_blopex_blocksize","Block size","EPSBLOPEXSetBlockSize",20,&bs,&flg)); |
381 |
6/8✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
36 | if (flg) PetscCall(EPSBLOPEXSetBlockSize(eps,bs)); |
382 | |||
383 |
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.
|
36 | PetscOptionsHeadEnd(); |
384 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
385 | } | ||
386 | |||
387 | 36 | SLEPC_EXTERN PetscErrorCode EPSCreate_BLOPEX(EPS eps) | |
388 | { | ||
389 | 36 | EPS_BLOPEX *ctx; | |
390 | |||
391 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | PetscFunctionBegin; |
392 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscNew(&ctx)); |
393 | 36 | eps->data = (void*)ctx; | |
394 | |||
395 | 36 | eps->categ = EPS_CATEGORY_PRECOND; | |
396 | |||
397 | 36 | eps->ops->solve = EPSSolve_BLOPEX; | |
398 | 36 | eps->ops->setup = EPSSetUp_BLOPEX; | |
399 | 36 | eps->ops->setupsort = EPSSetUpSort_Basic; | |
400 | 36 | eps->ops->setfromoptions = EPSSetFromOptions_BLOPEX; | |
401 | 36 | eps->ops->destroy = EPSDestroy_BLOPEX; | |
402 | 36 | eps->ops->reset = EPSReset_BLOPEX; | |
403 | 36 | eps->ops->view = EPSView_BLOPEX; | |
404 | 36 | eps->ops->backtransform = EPSBackTransform_Default; | |
405 | 36 | eps->ops->setdefaultst = EPSSetDefaultST_GMRES; | |
406 | |||
407 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXSetBlockSize_C",EPSBLOPEXSetBlockSize_BLOPEX)); |
408 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
36 | PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXGetBlockSize_C",EPSBLOPEXGetBlockSize_BLOPEX)); |
409 |
5/8✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
36 | if (slepc_blopex_useconstr < 0) PetscCall(PetscObjectComposedDataRegister(&slepc_blopex_useconstr)); |
410 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
411 | } | ||
412 |