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 eigensolvers in ChASE. | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/ | ||
15 | #include <petsc/private/petscscalapack.h> | ||
16 | |||
17 | #if !defined(PETSC_USE_COMPLEX) | ||
18 | #if defined(PETSC_USE_REAL_SINGLE) | ||
19 | /* s */ | ||
20 | #define CHASEchase_init_blockcyclic pschase_init_blockcyclic_ | ||
21 | #define CHASEchase pschase_ | ||
22 | #define CHASEchase_finalize pschase_finalize_ | ||
23 | #elif defined(PETSC_USE_REAL_DOUBLE) | ||
24 | /* d */ | ||
25 | #define CHASEchase_init_blockcyclic pdchase_init_blockcyclic_ | ||
26 | #define CHASEchase pdchase_ | ||
27 | #define CHASEchase_finalize pdchase_finalize_ | ||
28 | #endif | ||
29 | #else | ||
30 | #if defined(PETSC_USE_REAL_SINGLE) | ||
31 | /* c */ | ||
32 | #define CHASEchase_init_blockcyclic pcchase_init_blockcyclic_ | ||
33 | #define CHASEchase pcchase_ | ||
34 | #define CHASEchase_finalize pcchase_finalize_ | ||
35 | #elif defined(PETSC_USE_REAL_DOUBLE) | ||
36 | /* z */ | ||
37 | #define CHASEchase_init_blockcyclic pzchase_init_blockcyclic_ | ||
38 | #define CHASEchase pzchase_ | ||
39 | #define CHASEchase_finalize pzchase_finalize_ | ||
40 | #endif | ||
41 | #endif | ||
42 | |||
43 | SLEPC_EXTERN void CHASEchase_init_blockcyclic(PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscScalar*,PetscReal*,PetscInt*,PetscInt*,char*,PetscInt*,PetscInt*,MPI_Comm*,PetscInt*); | ||
44 | SLEPC_EXTERN void CHASEchase(PetscInt*,PetscReal*,char*,char*,char*); | ||
45 | SLEPC_EXTERN void CHASEchase_finalize(PetscInt*); | ||
46 | |||
47 | typedef struct { | ||
48 | Mat As; /* converted matrix */ | ||
49 | PetscInt nex; /* extra searching space size */ | ||
50 | PetscInt deg; /* initial degree of Chebyshev polynomial filter */ | ||
51 | PetscBool opt; /* internal optimization of polynomial degree */ | ||
52 | } EPS_ChASE; | ||
53 | |||
54 | 6 | static PetscErrorCode EPSSetUp_ChASE(EPS eps) | |
55 | { | ||
56 | 6 | EPS_ChASE *ctx = (EPS_ChASE*)eps->data; | |
57 | 6 | Mat A; | |
58 | 6 | PetscBool isshift; | |
59 | 6 | PetscScalar shift; | |
60 | |||
61 | 6 | PetscFunctionBegin; | |
62 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | EPSCheckStandard(eps); |
63 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | EPSCheckHermitian(eps); |
64 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | EPSCheckNotStructured(eps); |
65 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift)); |
66 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | PetscCheck(isshift,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support spectral transformations"); |
67 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | if (eps->nev==0) eps->nev = 1; |
68 |
2/4✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
|
6 | if (eps->ncv==PETSC_DETERMINE) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15)); |
69 | ✗ | else PetscCheck(eps->ncv>=eps->nev+1 || (eps->ncv==eps->nev && eps->ncv==eps->n),PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev+1"); | |
70 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n")); |
71 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6 | if (eps->max_it==PETSC_DETERMINE) eps->max_it = 1; |
72 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | if (!eps->which) eps->which = EPS_SMALLEST_REAL; |
73 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | PetscCheck(eps->which==EPS_SMALLEST_REAL,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver only supports computation of leftmost eigenvalues"); |
74 |
4/14✗ 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 taken 2 times.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
6 | EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION); |
75 |
3/12✗ 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.
|
6 | EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE | EPS_FEATURE_STOPPING); |
76 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(EPSAllocateSolution(eps,0)); |
77 | |||
78 | 6 | ctx->nex = eps->ncv-eps->nev; | |
79 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | if (ctx->deg==PETSC_DEFAULT) ctx->deg = 10; |
80 | |||
81 | /* convert matrices */ | ||
82 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(MatDestroy(&ctx->As)); |
83 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(STGetMatrix(eps->st,0,&A)); |
84 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(MatConvert(A,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->As)); |
85 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(STGetShift(eps->st,&shift)); |
86 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | if (shift != 0.0) PetscCall(MatShift(ctx->As,-shift)); |
87 | PetscFunctionReturn(PETSC_SUCCESS); | ||
88 | } | ||
89 | |||
90 | 6 | static PetscErrorCode EPSSolve_ChASE(EPS eps) | |
91 | { | ||
92 | 6 | EPS_ChASE *ctx = (EPS_ChASE*)eps->data; | |
93 | 6 | Mat As = ctx->As,Vs,V; | |
94 | 6 | Mat_ScaLAPACK *a = (Mat_ScaLAPACK*)As->data,*v; | |
95 | 6 | MPI_Comm comm; | |
96 | 6 | PetscReal *w = eps->errest; /* used to store real eigenvalues */ | |
97 | 6 | PetscInt i,init=0,flg=1; | |
98 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | char grid_major='R',mode='X',opt=ctx->opt?'S':'X',qr='C'; |
99 | |||
100 | 6 | PetscFunctionBegin; | |
101 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | PetscCheck(a->grid->npcol==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Current implementation of the ChASE interface only supports process grids with one column; use -mat_scalapack_grid_height <p> where <p> is the number of processes"); |
102 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(BVGetMat(eps->V,&V)); |
103 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(MatConvert(V,MATSCALAPACK,MAT_INITIAL_MATRIX,&Vs)); |
104 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(BVRestoreMat(eps->V,&V)); |
105 | 6 | v = (Mat_ScaLAPACK*)Vs->data; | |
106 | |||
107 | /* initialization */ | ||
108 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(PetscObjectGetComm((PetscObject)As,&comm)); |
109 | 6 | CHASEchase_init_blockcyclic(&a->N,&eps->nev,&ctx->nex,&a->mb,&a->nb,a->loc,&a->lld,v->loc,w,&a->grid->nprow,&a->grid->npcol,&grid_major,&a->rsrc,&a->csrc,&comm,&init); | |
110 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | PetscCheck(init==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Problem initializing ChASE"); |
111 | |||
112 | /* solve */ | ||
113 | 6 | CHASEchase(&ctx->deg,&eps->tol,&mode,&opt,&qr); | |
114 | 6 | CHASEchase_finalize(&flg); | |
115 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | PetscCheck(flg==0,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Problem solving with ChASE"); |
116 | |||
117 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
132 | for (i=0;i<eps->ncv;i++) { |
118 | 126 | eps->eigr[i] = eps->errest[i]; | |
119 | 126 | eps->errest[i] = eps->tol; | |
120 | } | ||
121 | |||
122 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(BVGetMat(eps->V,&V)); |
123 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(MatConvert(Vs,MATDENSE,MAT_REUSE_MATRIX,&V)); |
124 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(BVRestoreMat(eps->V,&V)); |
125 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(MatDestroy(&Vs)); |
126 | |||
127 | 6 | eps->nconv = eps->nev; | |
128 | 6 | eps->its = 1; | |
129 | 6 | eps->reason = EPS_CONVERGED_TOL; | |
130 | 6 | PetscFunctionReturn(PETSC_SUCCESS); | |
131 | } | ||
132 | |||
133 | 6 | static PetscErrorCode EPSCHASESetDegree_ChASE(EPS eps,PetscInt deg,PetscBool opt) | |
134 | { | ||
135 | 6 | EPS_ChASE *ctx = (EPS_ChASE*)eps->data; | |
136 | |||
137 | 6 | PetscFunctionBegin; | |
138 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | if (deg==PETSC_DEFAULT || deg==PETSC_DECIDE) { |
139 | ✗ | ctx->deg = PETSC_DEFAULT; | |
140 | ✗ | eps->state = EPS_STATE_INITIAL; | |
141 | } else { | ||
142 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
6 | PetscCheck(deg>0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"The degree must be >0"); |
143 | 6 | ctx->deg = deg; | |
144 | } | ||
145 | 6 | ctx->opt = opt; | |
146 | 6 | PetscFunctionReturn(PETSC_SUCCESS); | |
147 | } | ||
148 | |||
149 | /*@ | ||
150 | EPSCHASESetDegree - Sets the degree of the Chebyshev polynomial filter in the ChASE solver. | ||
151 | |||
152 | Logically Collective | ||
153 | |||
154 | Input Parameters: | ||
155 | + eps - the eigenproblem solver context | ||
156 | . deg - initial degree of Chebyshev polynomial filter | ||
157 | - opt - internal optimization of polynomial degree | ||
158 | |||
159 | Options Database Keys: | ||
160 | + -eps_chase_degree - Sets the initial degree | ||
161 | - -eps_chase_degree_opt - Enables/disables the optimization | ||
162 | |||
163 | Level: advanced | ||
164 | |||
165 | .seealso: EPSCHASEGetDegree() | ||
166 | @*/ | ||
167 | 6 | PetscErrorCode EPSCHASESetDegree(EPS eps,PetscInt deg,PetscBool opt) | |
168 | { | ||
169 | 6 | PetscFunctionBegin; | |
170 | 6 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); | |
171 | 6 | PetscValidLogicalCollectiveInt(eps,deg,2); | |
172 | 6 | PetscValidLogicalCollectiveBool(eps,opt,3); | |
173 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
6 | PetscTryMethod(eps,"EPSCHASESetDegree_C",(EPS,PetscInt,PetscBool),(eps,deg,opt)); |
174 | 6 | PetscFunctionReturn(PETSC_SUCCESS); | |
175 | } | ||
176 | |||
177 | 6 | static PetscErrorCode EPSCHASEGetDegree_ChASE(EPS eps,PetscInt *deg,PetscBool *opt) | |
178 | { | ||
179 | 6 | EPS_ChASE *ctx = (EPS_ChASE*)eps->data; | |
180 | |||
181 | 6 | PetscFunctionBegin; | |
182 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6 | if (deg) *deg = ctx->deg; |
183 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6 | if (opt) *opt = ctx->opt; |
184 | 6 | PetscFunctionReturn(PETSC_SUCCESS); | |
185 | } | ||
186 | |||
187 | /*@ | ||
188 | EPSCHASEGetDegree - Gets the degree of the Chebyshev polynomial filter used in the ChASE solver. | ||
189 | |||
190 | Not Collective | ||
191 | |||
192 | Input Parameter: | ||
193 | . eps - the eigenproblem solver context | ||
194 | |||
195 | Output Parameters: | ||
196 | + deg - initial degree of Chebyshev polynomial filter | ||
197 | - opt - internal optimization of polynomial degree | ||
198 | |||
199 | Level: advanced | ||
200 | |||
201 | .seealso: EPSCHASESetDegree() | ||
202 | @*/ | ||
203 | 6 | PetscErrorCode EPSCHASEGetDegree(EPS eps,PetscInt *deg,PetscBool *opt) | |
204 | { | ||
205 | 6 | PetscFunctionBegin; | |
206 | 6 | PetscValidHeaderSpecific(eps,EPS_CLASSID,1); | |
207 |
3/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
6 | PetscUseMethod(eps,"EPSCHASEGetDegree_C",(EPS,PetscInt*,PetscBool*),(eps,deg,opt)); |
208 | 6 | PetscFunctionReturn(PETSC_SUCCESS); | |
209 | } | ||
210 | |||
211 | 6 | static PetscErrorCode EPSSetFromOptions_ChASE(EPS eps,PetscOptionItems PetscOptionsObject) | |
212 | { | ||
213 | 6 | EPS_ChASE *ctx = (EPS_ChASE*)eps->data; | |
214 | 6 | PetscBool flg1,flg2; | |
215 | 6 | PetscInt deg; | |
216 | 6 | PetscBool opt; | |
217 | |||
218 | 6 | PetscFunctionBegin; | |
219 |
1/8✗ 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.
|
6 | PetscOptionsHeadBegin(PetscOptionsObject,"EPS ChASE Options"); |
220 | |||
221 | 6 | deg = ctx->deg; | |
222 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(PetscOptionsInt("-eps_chase_degree","Initial degree of Chebyshev polynomial filter","EPSCHASESetDegree",deg,°,&flg1)); |
223 | 6 | opt = ctx->opt; | |
224 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(PetscOptionsBool("-eps_chase_degree_opt","Internal optimization of polynomial degree","EPSCHASESetDegree",opt,&opt,&flg2)); |
225 |
2/6✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
6 | if (flg1 || flg2) PetscCall(EPSCHASESetDegree(eps,deg,opt)); |
226 | |||
227 | 6 | PetscOptionsHeadEnd(); | |
228 | PetscFunctionReturn(PETSC_SUCCESS); | ||
229 | } | ||
230 | |||
231 | 6 | static PetscErrorCode EPSDestroy_ChASE(EPS eps) | |
232 | { | ||
233 | 6 | PetscFunctionBegin; | |
234 |
2/4✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
|
6 | PetscCall(PetscFree(eps->data)); |
235 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCHASESetDegree_C",NULL)); |
236 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCHASEGetDegree_C",NULL)); |
237 | PetscFunctionReturn(PETSC_SUCCESS); | ||
238 | } | ||
239 | |||
240 | 6 | static PetscErrorCode EPSReset_ChASE(EPS eps) | |
241 | { | ||
242 | 6 | EPS_ChASE *ctx = (EPS_ChASE*)eps->data; | |
243 | |||
244 | 6 | PetscFunctionBegin; | |
245 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(MatDestroy(&ctx->As)); |
246 | PetscFunctionReturn(PETSC_SUCCESS); | ||
247 | } | ||
248 | |||
249 | ✗ | static PetscErrorCode EPSView_ChASE(EPS eps,PetscViewer viewer) | |
250 | { | ||
251 | ✗ | EPS_ChASE *ctx = (EPS_ChASE*)eps->data; | |
252 | ✗ | PetscBool isascii; | |
253 | |||
254 | ✗ | PetscFunctionBegin; | |
255 | ✗ | PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); | |
256 | ✗ | if (isascii) { | |
257 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer," initial degree of Chebyshev polynomial filter: %" PetscInt_FMT "\n",ctx->deg)); | |
258 | ✗ | if (ctx->opt) PetscCall(PetscViewerASCIIPrintf(viewer," internal optimization of polynomial degree enabled\n")); | |
259 | } | ||
260 | PetscFunctionReturn(PETSC_SUCCESS); | ||
261 | } | ||
262 | |||
263 | 6 | SLEPC_EXTERN PetscErrorCode EPSCreate_ChASE(EPS eps) | |
264 | { | ||
265 | 6 | EPS_ChASE *ctx; | |
266 | |||
267 | 6 | PetscFunctionBegin; | |
268 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(PetscNew(&ctx)); |
269 | 6 | eps->data = (void*)ctx; | |
270 | 6 | ctx->deg = PETSC_DEFAULT; | |
271 | 6 | ctx->opt = PETSC_TRUE; | |
272 | |||
273 | 6 | eps->categ = EPS_CATEGORY_OTHER; | |
274 | |||
275 | 6 | eps->ops->solve = EPSSolve_ChASE; | |
276 | 6 | eps->ops->setup = EPSSetUp_ChASE; | |
277 | 6 | eps->ops->setupsort = EPSSetUpSort_Basic; | |
278 | 6 | eps->ops->setfromoptions = EPSSetFromOptions_ChASE; | |
279 | 6 | eps->ops->destroy = EPSDestroy_ChASE; | |
280 | 6 | eps->ops->reset = EPSReset_ChASE; | |
281 | 6 | eps->ops->view = EPSView_ChASE; | |
282 | 6 | eps->ops->backtransform = EPSBackTransform_Default; | |
283 | 6 | eps->ops->setdefaultst = EPSSetDefaultST_NoFactor; | |
284 | |||
285 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCHASESetDegree_C",EPSCHASESetDegree_ChASE)); |
286 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSCHASEGetDegree_C",EPSCHASEGetDegree_ChASE)); |
287 | PetscFunctionReturn(PETSC_SUCCESS); | ||
288 | } | ||
289 |