Line data Source code
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: "jd"
12 :
13 : Method: Jacobi-Davidson
14 :
15 : Algorithm:
16 :
17 : Jacobi-Davidson with various subspace extraction and
18 : restart techniques.
19 :
20 : References:
21 :
22 : [1] G.L.G. Sleijpen and H.A. van der Vorst, "A Jacobi-Davidson
23 : iteration method for linear eigenvalue problems", SIAM J.
24 : Matrix Anal. Appl. 17(2):401-425, 1996.
25 :
26 : [2] E. Romero and J.E. Roman, "A parallel implementation of
27 : Davidson methods for large-scale eigenvalue problems in
28 : SLEPc", ACM Trans. Math. Software 40(2), Article 13, 2014.
29 : */
30 :
31 : #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
32 : #include <../src/eps/impls/davidson/davidson.h>
33 :
34 20 : static PetscErrorCode EPSSetFromOptions_JD(EPS eps,PetscOptionItems *PetscOptionsObject)
35 : {
36 20 : PetscBool flg,flg2,op,orth;
37 20 : PetscInt opi,opi0;
38 20 : PetscReal opf;
39 :
40 20 : PetscFunctionBegin;
41 20 : PetscOptionsHeadBegin(PetscOptionsObject,"EPS Jacobi-Davidson (JD) Options");
42 :
43 20 : PetscCall(EPSJDGetKrylovStart(eps,&op));
44 20 : PetscCall(PetscOptionsBool("-eps_jd_krylov_start","Start the search subspace with a Krylov basis","EPSJDSetKrylovStart",op,&op,&flg));
45 20 : if (flg) PetscCall(EPSJDSetKrylovStart(eps,op));
46 :
47 20 : PetscCall(EPSJDGetBOrth(eps,&orth));
48 20 : PetscCall(PetscOptionsBool("-eps_jd_borth","Use B-orthogonalization in the search subspace","EPSJDSetBOrth",op,&op,&flg));
49 20 : if (flg) PetscCall(EPSJDSetBOrth(eps,op));
50 :
51 20 : PetscCall(EPSJDGetBlockSize(eps,&opi));
52 20 : PetscCall(PetscOptionsInt("-eps_jd_blocksize","Number of vectors to add to the search subspace","EPSJDSetBlockSize",opi,&opi,&flg));
53 20 : if (flg) PetscCall(EPSJDSetBlockSize(eps,opi));
54 :
55 20 : PetscCall(EPSJDGetRestart(eps,&opi,&opi0));
56 20 : PetscCall(PetscOptionsInt("-eps_jd_minv","Size of the search subspace after restarting","EPSJDSetRestart",opi,&opi,&flg));
57 20 : PetscCall(PetscOptionsInt("-eps_jd_plusk","Number of eigenvectors saved from the previous iteration when restarting","EPSJDSetRestart",opi0,&opi0,&flg2));
58 20 : if (flg || flg2) PetscCall(EPSJDSetRestart(eps,opi,opi0));
59 :
60 20 : PetscCall(EPSJDGetInitialSize(eps,&opi));
61 20 : PetscCall(PetscOptionsInt("-eps_jd_initial_size","Initial size of the search subspace","EPSJDSetInitialSize",opi,&opi,&flg));
62 20 : if (flg) PetscCall(EPSJDSetInitialSize(eps,opi));
63 :
64 20 : PetscCall(EPSJDGetFix(eps,&opf));
65 20 : PetscCall(PetscOptionsReal("-eps_jd_fix","Tolerance for changing the target in the correction equation","EPSJDSetFix",opf,&opf,&flg));
66 20 : if (flg) PetscCall(EPSJDSetFix(eps,opf));
67 :
68 20 : PetscCall(EPSJDGetConstCorrectionTol(eps,&op));
69 20 : PetscCall(PetscOptionsBool("-eps_jd_const_correction_tol","Disable the dynamic stopping criterion when solving the correction equation","EPSJDSetConstCorrectionTol",op,&op,&flg));
70 20 : if (flg) PetscCall(EPSJDSetConstCorrectionTol(eps,op));
71 :
72 20 : PetscOptionsHeadEnd();
73 20 : PetscFunctionReturn(PETSC_SUCCESS);
74 : }
75 :
76 51 : static PetscErrorCode EPSSetDefaultST_JD(EPS eps)
77 : {
78 51 : KSP ksp;
79 :
80 51 : PetscFunctionBegin;
81 51 : if (!((PetscObject)eps->st)->type_name) {
82 20 : PetscCall(STSetType(eps->st,STPRECOND));
83 20 : PetscCall(STPrecondSetKSPHasMat(eps->st,PETSC_TRUE));
84 : }
85 51 : PetscCall(STGetKSP(eps->st,&ksp));
86 51 : if (!((PetscObject)ksp)->type_name) {
87 21 : PetscCall(KSPSetType(ksp,KSPBCGSL));
88 21 : PetscCall(KSPSetTolerances(ksp,1e-4,PETSC_CURRENT,PETSC_CURRENT,90));
89 : }
90 51 : PetscFunctionReturn(PETSC_SUCCESS);
91 : }
92 :
93 31 : static PetscErrorCode EPSSetUp_JD(EPS eps)
94 : {
95 31 : PetscBool t;
96 31 : KSP ksp;
97 :
98 31 : PetscFunctionBegin;
99 : /* Setup common for all davidson solvers */
100 31 : PetscCall(EPSSetUp_XD(eps));
101 :
102 : /* Check some constraints */
103 31 : PetscCall(STGetKSP(eps->st,&ksp));
104 31 : PetscCall(PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&t));
105 31 : PetscCheck(!t,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"EPSJD does not work with KSPPREONLY");
106 31 : PetscFunctionReturn(PETSC_SUCCESS);
107 : }
108 :
109 0 : static PetscErrorCode EPSView_JD(EPS eps,PetscViewer viewer)
110 : {
111 0 : PetscBool isascii,opb;
112 0 : PetscReal opf;
113 0 : PetscInt opi,opi0;
114 :
115 0 : PetscFunctionBegin;
116 0 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
117 0 : if (isascii) {
118 0 : PetscCall(EPSXDGetBOrth_XD(eps,&opb));
119 0 : if (opb) PetscCall(PetscViewerASCIIPrintf(viewer," search subspace is B-orthogonalized\n"));
120 0 : else PetscCall(PetscViewerASCIIPrintf(viewer," search subspace is orthogonalized\n"));
121 0 : PetscCall(EPSXDGetBlockSize_XD(eps,&opi));
122 0 : PetscCall(PetscViewerASCIIPrintf(viewer," block size=%" PetscInt_FMT "\n",opi));
123 0 : PetscCall(EPSXDGetKrylovStart_XD(eps,&opb));
124 0 : if (!opb) PetscCall(PetscViewerASCIIPrintf(viewer," type of the initial subspace: non-Krylov\n"));
125 0 : else PetscCall(PetscViewerASCIIPrintf(viewer," type of the initial subspace: Krylov\n"));
126 0 : PetscCall(EPSXDGetRestart_XD(eps,&opi,&opi0));
127 0 : PetscCall(PetscViewerASCIIPrintf(viewer," size of the subspace after restarting: %" PetscInt_FMT "\n",opi));
128 0 : PetscCall(PetscViewerASCIIPrintf(viewer," number of vectors after restarting from the previous iteration: %" PetscInt_FMT "\n",opi0));
129 :
130 0 : PetscCall(EPSJDGetFix_JD(eps,&opf));
131 0 : PetscCall(PetscViewerASCIIPrintf(viewer," threshold for changing the target in the correction equation (fix): %g\n",(double)opf));
132 :
133 0 : PetscCall(EPSJDGetConstCorrectionTol_JD(eps,&opb));
134 0 : if (!opb) PetscCall(PetscViewerASCIIPrintf(viewer," using dynamic tolerance for the correction equation\n"));
135 : }
136 0 : PetscFunctionReturn(PETSC_SUCCESS);
137 : }
138 :
139 21 : static PetscErrorCode EPSDestroy_JD(EPS eps)
140 : {
141 21 : PetscFunctionBegin;
142 21 : PetscCall(PetscFree(eps->data));
143 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetKrylovStart_C",NULL));
144 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetKrylovStart_C",NULL));
145 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetBlockSize_C",NULL));
146 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetBlockSize_C",NULL));
147 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetRestart_C",NULL));
148 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetRestart_C",NULL));
149 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetInitialSize_C",NULL));
150 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetInitialSize_C",NULL));
151 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetFix_C",NULL));
152 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetFix_C",NULL));
153 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetConstCorrectionTol_C",NULL));
154 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetConstCorrectionTol_C",NULL));
155 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetBOrth_C",NULL));
156 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetBOrth_C",NULL));
157 21 : PetscFunctionReturn(PETSC_SUCCESS);
158 : }
159 :
160 : /*@
161 : EPSJDSetKrylovStart - Activates or deactivates starting the searching
162 : subspace with a Krylov basis.
163 :
164 : Logically Collective
165 :
166 : Input Parameters:
167 : + eps - the eigenproblem solver context
168 : - krylovstart - boolean flag
169 :
170 : Options Database Key:
171 : . -eps_jd_krylov_start - Activates starting the searching subspace with a
172 : Krylov basis
173 :
174 : Level: advanced
175 :
176 : .seealso: EPSJDGetKrylovStart()
177 : @*/
178 1 : PetscErrorCode EPSJDSetKrylovStart(EPS eps,PetscBool krylovstart)
179 : {
180 1 : PetscFunctionBegin;
181 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
182 3 : PetscValidLogicalCollectiveBool(eps,krylovstart,2);
183 1 : PetscTryMethod(eps,"EPSJDSetKrylovStart_C",(EPS,PetscBool),(eps,krylovstart));
184 1 : PetscFunctionReturn(PETSC_SUCCESS);
185 : }
186 :
187 : /*@
188 : EPSJDGetKrylovStart - Returns a flag indicating if the searching subspace is started with a
189 : Krylov basis.
190 :
191 : Not Collective
192 :
193 : Input Parameter:
194 : . eps - the eigenproblem solver context
195 :
196 : Output Parameters:
197 : . krylovstart - boolean flag indicating if the searching subspace is started
198 : with a Krylov basis
199 :
200 : Level: advanced
201 :
202 : .seealso: EPSJDSetKrylovStart()
203 : @*/
204 20 : PetscErrorCode EPSJDGetKrylovStart(EPS eps,PetscBool *krylovstart)
205 : {
206 20 : PetscFunctionBegin;
207 20 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
208 20 : PetscAssertPointer(krylovstart,2);
209 20 : PetscUseMethod(eps,"EPSJDGetKrylovStart_C",(EPS,PetscBool*),(eps,krylovstart));
210 20 : PetscFunctionReturn(PETSC_SUCCESS);
211 : }
212 :
213 : /*@
214 : EPSJDSetBlockSize - Sets the number of vectors to be added to the searching space
215 : in every iteration.
216 :
217 : Logically Collective
218 :
219 : Input Parameters:
220 : + eps - the eigenproblem solver context
221 : - blocksize - number of vectors added to the search space in every iteration
222 :
223 : Options Database Key:
224 : . -eps_jd_blocksize - number of vectors added to the searching space every iteration
225 :
226 : Level: advanced
227 :
228 : .seealso: EPSJDSetKrylovStart()
229 : @*/
230 1 : PetscErrorCode EPSJDSetBlockSize(EPS eps,PetscInt blocksize)
231 : {
232 1 : PetscFunctionBegin;
233 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
234 3 : PetscValidLogicalCollectiveInt(eps,blocksize,2);
235 1 : PetscTryMethod(eps,"EPSJDSetBlockSize_C",(EPS,PetscInt),(eps,blocksize));
236 1 : PetscFunctionReturn(PETSC_SUCCESS);
237 : }
238 :
239 : /*@
240 : EPSJDGetBlockSize - Returns the number of vectors to be added to the searching space
241 : in every iteration.
242 :
243 : Not Collective
244 :
245 : Input Parameter:
246 : . eps - the eigenproblem solver context
247 :
248 : Output Parameter:
249 : . blocksize - number of vectors added to the search space in every iteration
250 :
251 : Level: advanced
252 :
253 : .seealso: EPSJDSetBlockSize()
254 : @*/
255 20 : PetscErrorCode EPSJDGetBlockSize(EPS eps,PetscInt *blocksize)
256 : {
257 20 : PetscFunctionBegin;
258 20 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
259 20 : PetscAssertPointer(blocksize,2);
260 20 : PetscUseMethod(eps,"EPSJDGetBlockSize_C",(EPS,PetscInt*),(eps,blocksize));
261 20 : PetscFunctionReturn(PETSC_SUCCESS);
262 : }
263 :
264 : /*@
265 : EPSJDSetRestart - Sets the number of vectors of the searching space after
266 : restarting and the number of vectors saved from the previous iteration.
267 :
268 : Logically Collective
269 :
270 : Input Parameters:
271 : + eps - the eigenproblem solver context
272 : . minv - number of vectors of the searching subspace after restarting
273 : - plusk - number of vectors saved from the previous iteration
274 :
275 : Options Database Keys:
276 : + -eps_jd_minv - number of vectors of the searching subspace after restarting
277 : - -eps_jd_plusk - number of vectors saved from the previous iteration
278 :
279 : Note:
280 : PETSC_CURRENT can be used to preserve the current value of any of the
281 : arguments, and PETSC_DETERMINE to set them to a default value.
282 :
283 : Level: advanced
284 :
285 : .seealso: EPSJDGetRestart()
286 : @*/
287 2 : PetscErrorCode EPSJDSetRestart(EPS eps,PetscInt minv,PetscInt plusk)
288 : {
289 2 : PetscFunctionBegin;
290 2 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
291 6 : PetscValidLogicalCollectiveInt(eps,minv,2);
292 6 : PetscValidLogicalCollectiveInt(eps,plusk,3);
293 2 : PetscTryMethod(eps,"EPSJDSetRestart_C",(EPS,PetscInt,PetscInt),(eps,minv,plusk));
294 2 : PetscFunctionReturn(PETSC_SUCCESS);
295 : }
296 :
297 : /*@
298 : EPSJDGetRestart - Gets the number of vectors of the searching space after
299 : restarting and the number of vectors saved from the previous iteration.
300 :
301 : Not Collective
302 :
303 : Input Parameter:
304 : . eps - the eigenproblem solver context
305 :
306 : Output Parameters:
307 : + minv - number of vectors of the searching subspace after restarting
308 : - plusk - number of vectors saved from the previous iteration
309 :
310 : Level: advanced
311 :
312 : .seealso: EPSJDSetRestart()
313 : @*/
314 20 : PetscErrorCode EPSJDGetRestart(EPS eps,PetscInt *minv,PetscInt *plusk)
315 : {
316 20 : PetscFunctionBegin;
317 20 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
318 20 : PetscUseMethod(eps,"EPSJDGetRestart_C",(EPS,PetscInt*,PetscInt*),(eps,minv,plusk));
319 20 : PetscFunctionReturn(PETSC_SUCCESS);
320 : }
321 :
322 : /*@
323 : EPSJDSetInitialSize - Sets the initial size of the searching space.
324 :
325 : Logically Collective
326 :
327 : Input Parameters:
328 : + eps - the eigenproblem solver context
329 : - initialsize - number of vectors of the initial searching subspace
330 :
331 : Options Database Key:
332 : . -eps_jd_initial_size - number of vectors of the initial searching subspace
333 :
334 : Notes:
335 : If EPSJDGetKrylovStart() is PETSC_FALSE and the user provides vectors with
336 : EPSSetInitialSpace(), up to initialsize vectors will be used; and if the
337 : provided vectors are not enough, the solver completes the subspace with
338 : random vectors. In the case of EPSJDGetKrylovStart() being PETSC_TRUE, the solver
339 : gets the first vector provided by the user or, if not available, a random vector,
340 : and expands the Krylov basis up to initialsize vectors.
341 :
342 : Level: advanced
343 :
344 : .seealso: EPSJDGetInitialSize(), EPSJDGetKrylovStart()
345 : @*/
346 1 : PetscErrorCode EPSJDSetInitialSize(EPS eps,PetscInt initialsize)
347 : {
348 1 : PetscFunctionBegin;
349 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
350 3 : PetscValidLogicalCollectiveInt(eps,initialsize,2);
351 1 : PetscTryMethod(eps,"EPSJDSetInitialSize_C",(EPS,PetscInt),(eps,initialsize));
352 1 : PetscFunctionReturn(PETSC_SUCCESS);
353 : }
354 :
355 : /*@
356 : EPSJDGetInitialSize - Returns the initial size of the searching space.
357 :
358 : Not Collective
359 :
360 : Input Parameter:
361 : . eps - the eigenproblem solver context
362 :
363 : Output Parameter:
364 : . initialsize - number of vectors of the initial searching subspace
365 :
366 : Notes:
367 : If EPSJDGetKrylovStart() is PETSC_FALSE and the user provides vectors with
368 : EPSSetInitialSpace(), up to initialsize vectors will be used; and if the
369 : provided vectors are not enough, the solver completes the subspace with
370 : random vectors. In the case of EPSJDGetKrylovStart() being PETSC_TRUE, the solver
371 : gets the first vector provided by the user or, if not available, a random vector,
372 : and expands the Krylov basis up to initialsize vectors.
373 :
374 : Level: advanced
375 :
376 : .seealso: EPSJDSetInitialSize(), EPSJDGetKrylovStart()
377 : @*/
378 20 : PetscErrorCode EPSJDGetInitialSize(EPS eps,PetscInt *initialsize)
379 : {
380 20 : PetscFunctionBegin;
381 20 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
382 20 : PetscAssertPointer(initialsize,2);
383 20 : PetscUseMethod(eps,"EPSJDGetInitialSize_C",(EPS,PetscInt*),(eps,initialsize));
384 20 : PetscFunctionReturn(PETSC_SUCCESS);
385 : }
386 :
387 1 : static PetscErrorCode EPSJDSetFix_JD(EPS eps,PetscReal fix)
388 : {
389 1 : EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
390 :
391 1 : PetscFunctionBegin;
392 1 : if (fix == (PetscReal)PETSC_DEFAULT || fix == (PetscReal)PETSC_DECIDE) fix = 0.01;
393 1 : PetscCheck(fix>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value, must be >0");
394 1 : data->fix = fix;
395 1 : PetscFunctionReturn(PETSC_SUCCESS);
396 : }
397 :
398 : /*@
399 : EPSJDSetFix - Sets the threshold for changing the target in the correction
400 : equation.
401 :
402 : Logically Collective
403 :
404 : Input Parameters:
405 : + eps - the eigenproblem solver context
406 : - fix - threshold for changing the target
407 :
408 : Options Database Key:
409 : . -eps_jd_fix - the fix value
410 :
411 : Note:
412 : The target in the correction equation is fixed at the first iterations.
413 : When the norm of the residual vector is lower than the fix value,
414 : the target is set to the corresponding eigenvalue.
415 :
416 : Level: advanced
417 :
418 : .seealso: EPSJDGetFix()
419 : @*/
420 1 : PetscErrorCode EPSJDSetFix(EPS eps,PetscReal fix)
421 : {
422 1 : PetscFunctionBegin;
423 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
424 3 : PetscValidLogicalCollectiveReal(eps,fix,2);
425 1 : PetscTryMethod(eps,"EPSJDSetFix_C",(EPS,PetscReal),(eps,fix));
426 1 : PetscFunctionReturn(PETSC_SUCCESS);
427 : }
428 :
429 20 : PetscErrorCode EPSJDGetFix_JD(EPS eps,PetscReal *fix)
430 : {
431 20 : EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
432 :
433 20 : PetscFunctionBegin;
434 20 : *fix = data->fix;
435 20 : PetscFunctionReturn(PETSC_SUCCESS);
436 : }
437 :
438 : /*@
439 : EPSJDGetFix - Returns the threshold for changing the target in the correction
440 : equation.
441 :
442 : Not Collective
443 :
444 : Input Parameter:
445 : . eps - the eigenproblem solver context
446 :
447 : Output Parameter:
448 : . fix - threshold for changing the target
449 :
450 : Note:
451 : The target in the correction equation is fixed at the first iterations.
452 : When the norm of the residual vector is lower than the fix value,
453 : the target is set to the corresponding eigenvalue.
454 :
455 : Level: advanced
456 :
457 : .seealso: EPSJDSetFix()
458 : @*/
459 20 : PetscErrorCode EPSJDGetFix(EPS eps,PetscReal *fix)
460 : {
461 20 : PetscFunctionBegin;
462 20 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
463 20 : PetscAssertPointer(fix,2);
464 20 : PetscUseMethod(eps,"EPSJDGetFix_C",(EPS,PetscReal*),(eps,fix));
465 20 : PetscFunctionReturn(PETSC_SUCCESS);
466 : }
467 :
468 1 : static PetscErrorCode EPSJDSetConstCorrectionTol_JD(EPS eps,PetscBool constant)
469 : {
470 1 : EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
471 :
472 1 : PetscFunctionBegin;
473 1 : data->dynamic = PetscNot(constant);
474 1 : PetscFunctionReturn(PETSC_SUCCESS);
475 : }
476 :
477 : /*@
478 : EPSJDSetConstCorrectionTol - If true, deactivates the dynamic stopping criterion
479 : (also called Newton) that sets the KSP relative tolerance
480 : to 0.5**i, where i is the number of EPS iterations from the last converged value.
481 :
482 : Logically Collective
483 :
484 : Input Parameters:
485 : + eps - the eigenproblem solver context
486 : - constant - if false, the KSP relative tolerance is set to 0.5**i.
487 :
488 : Options Database Key:
489 : . -eps_jd_const_correction_tol - Deactivates the dynamic stopping criterion.
490 :
491 : Level: advanced
492 :
493 : .seealso: EPSJDGetConstCorrectionTol()
494 : @*/
495 1 : PetscErrorCode EPSJDSetConstCorrectionTol(EPS eps,PetscBool constant)
496 : {
497 1 : PetscFunctionBegin;
498 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
499 3 : PetscValidLogicalCollectiveBool(eps,constant,2);
500 1 : PetscTryMethod(eps,"EPSJDSetConstCorrectionTol_C",(EPS,PetscBool),(eps,constant));
501 1 : PetscFunctionReturn(PETSC_SUCCESS);
502 : }
503 :
504 20 : PetscErrorCode EPSJDGetConstCorrectionTol_JD(EPS eps,PetscBool *constant)
505 : {
506 20 : EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
507 :
508 20 : PetscFunctionBegin;
509 20 : *constant = PetscNot(data->dynamic);
510 20 : PetscFunctionReturn(PETSC_SUCCESS);
511 : }
512 :
513 : /*@
514 : EPSJDGetConstCorrectionTol - Returns a flag indicating if the dynamic stopping is being used for
515 : solving the correction equation.
516 :
517 : Not Collective
518 :
519 : Input Parameter:
520 : . eps - the eigenproblem solver context
521 :
522 : Output Parameter:
523 : . constant - boolean flag indicating if the dynamic stopping criterion is not being used.
524 :
525 : Notes:
526 : If the flag is false the KSP relative tolerance is set to 0.5**i, where i is the number
527 : of EPS iterations from the last converged value.
528 :
529 : Level: advanced
530 :
531 : .seealso: EPSJDSetConstCorrectionTol()
532 : @*/
533 20 : PetscErrorCode EPSJDGetConstCorrectionTol(EPS eps,PetscBool *constant)
534 : {
535 20 : PetscFunctionBegin;
536 20 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
537 20 : PetscAssertPointer(constant,2);
538 20 : PetscUseMethod(eps,"EPSJDGetConstCorrectionTol_C",(EPS,PetscBool*),(eps,constant));
539 20 : PetscFunctionReturn(PETSC_SUCCESS);
540 : }
541 :
542 : /*@
543 : EPSJDSetBOrth - Selects the orthogonalization that will be used in the search
544 : subspace in case of generalized Hermitian problems.
545 :
546 : Logically Collective
547 :
548 : Input Parameters:
549 : + eps - the eigenproblem solver context
550 : - borth - whether to B-orthogonalize the search subspace
551 :
552 : Options Database Key:
553 : . -eps_jd_borth - Set the orthogonalization used in the search subspace
554 :
555 : Level: advanced
556 :
557 : .seealso: EPSJDGetBOrth()
558 : @*/
559 1 : PetscErrorCode EPSJDSetBOrth(EPS eps,PetscBool borth)
560 : {
561 1 : PetscFunctionBegin;
562 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
563 3 : PetscValidLogicalCollectiveBool(eps,borth,2);
564 1 : PetscTryMethod(eps,"EPSJDSetBOrth_C",(EPS,PetscBool),(eps,borth));
565 1 : PetscFunctionReturn(PETSC_SUCCESS);
566 : }
567 :
568 : /*@
569 : EPSJDGetBOrth - Returns the orthogonalization used in the search
570 : subspace in case of generalized Hermitian problems.
571 :
572 : Not Collective
573 :
574 : Input Parameter:
575 : . eps - the eigenproblem solver context
576 :
577 : Output Parameters:
578 : . borth - whether to B-orthogonalize the search subspace
579 :
580 : Level: advanced
581 :
582 : .seealso: EPSJDSetBOrth()
583 : @*/
584 20 : PetscErrorCode EPSJDGetBOrth(EPS eps,PetscBool *borth)
585 : {
586 20 : PetscFunctionBegin;
587 20 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
588 20 : PetscAssertPointer(borth,2);
589 20 : PetscUseMethod(eps,"EPSJDGetBOrth_C",(EPS,PetscBool*),(eps,borth));
590 20 : PetscFunctionReturn(PETSC_SUCCESS);
591 : }
592 :
593 21 : SLEPC_EXTERN PetscErrorCode EPSCreate_JD(EPS eps)
594 : {
595 21 : EPS_DAVIDSON *data;
596 :
597 21 : PetscFunctionBegin;
598 21 : PetscCall(PetscNew(&data));
599 21 : eps->data = (void*)data;
600 :
601 21 : data->blocksize = 1;
602 21 : data->initialsize = 0;
603 21 : data->minv = 0;
604 21 : data->plusk = PETSC_DETERMINE;
605 21 : data->ipB = PETSC_TRUE;
606 21 : data->fix = 0.01;
607 21 : data->krylovstart = PETSC_FALSE;
608 21 : data->dynamic = PETSC_FALSE;
609 :
610 21 : eps->useds = PETSC_TRUE;
611 21 : eps->categ = EPS_CATEGORY_PRECOND;
612 :
613 21 : eps->ops->solve = EPSSolve_XD;
614 21 : eps->ops->setup = EPSSetUp_JD;
615 21 : eps->ops->setupsort = EPSSetUpSort_Default;
616 21 : eps->ops->setfromoptions = EPSSetFromOptions_JD;
617 21 : eps->ops->destroy = EPSDestroy_JD;
618 21 : eps->ops->reset = EPSReset_XD;
619 21 : eps->ops->view = EPSView_JD;
620 21 : eps->ops->backtransform = EPSBackTransform_Default;
621 21 : eps->ops->computevectors = EPSComputeVectors_XD;
622 21 : eps->ops->setdefaultst = EPSSetDefaultST_JD;
623 :
624 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetKrylovStart_C",EPSXDSetKrylovStart_XD));
625 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetKrylovStart_C",EPSXDGetKrylovStart_XD));
626 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetBlockSize_C",EPSXDSetBlockSize_XD));
627 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetBlockSize_C",EPSXDGetBlockSize_XD));
628 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetRestart_C",EPSXDSetRestart_XD));
629 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetRestart_C",EPSXDGetRestart_XD));
630 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetInitialSize_C",EPSXDSetInitialSize_XD));
631 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetInitialSize_C",EPSXDGetInitialSize_XD));
632 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetFix_C",EPSJDSetFix_JD));
633 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetFix_C",EPSJDGetFix_JD));
634 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetConstCorrectionTol_C",EPSJDSetConstCorrectionTol_JD));
635 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetConstCorrectionTol_C",EPSJDGetConstCorrectionTol_JD));
636 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetBOrth_C",EPSXDSetBOrth_XD));
637 21 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetBOrth_C",EPSXDGetBOrth_XD));
638 21 : PetscFunctionReturn(PETSC_SUCCESS);
639 : }
|