Actual source code: epsbasic.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: Basic EPS routines
12: */
14: #include <slepc/private/epsimpl.h>
16: /* Logging support */
17: PetscClassId EPS_CLASSID = 0;
18: PetscLogEvent EPS_SetUp = 0,EPS_Solve = 0,EPS_CISS_SVD = 0;
20: /* List of registered EPS routines */
21: PetscFunctionList EPSList = NULL;
22: PetscBool EPSRegisterAllCalled = PETSC_FALSE;
24: /* List of registered EPS monitors */
25: PetscFunctionList EPSMonitorList = NULL;
26: PetscFunctionList EPSMonitorCreateList = NULL;
27: PetscFunctionList EPSMonitorDestroyList = NULL;
28: PetscBool EPSMonitorRegisterAllCalled = PETSC_FALSE;
30: /*@
31: EPSCreate - Creates the `EPS` context.
33: Collective
35: Input Parameter:
36: . comm - MPI communicator
38: Output Parameter:
39: . outeps - location to put the `EPS` context
41: Note:
42: The default `EPS` type is `EPSKRYLOVSCHUR`.
44: Level: beginner
46: .seealso: [](ch:eps), `EPSSetUp()`, `EPSSolve()`, `EPSDestroy()`, `EPS`
47: @*/
48: PetscErrorCode EPSCreate(MPI_Comm comm,EPS *outeps)
49: {
50: EPS eps;
52: PetscFunctionBegin;
53: PetscAssertPointer(outeps,2);
54: PetscCall(EPSInitializePackage());
55: PetscCall(SlepcHeaderCreate(eps,EPS_CLASSID,"EPS","Eigenvalue Problem Solver","EPS",comm,EPSDestroy,EPSView));
57: eps->max_it = PETSC_DETERMINE;
58: eps->nev = 0;
59: eps->ncv = PETSC_DETERMINE;
60: eps->mpd = PETSC_DETERMINE;
61: eps->nini = 0;
62: eps->nds = 0;
63: eps->target = 0.0;
64: eps->tol = PETSC_DETERMINE;
65: eps->thres = PETSC_MIN_REAL;
66: eps->threlative = PETSC_FALSE;
67: eps->conv = EPS_CONV_REL;
68: eps->stop = EPS_STOP_BASIC;
69: eps->which = (EPSWhich)0;
70: eps->inta = 0.0;
71: eps->intb = 0.0;
72: eps->problem_type = (EPSProblemType)0;
73: eps->extraction = EPS_RITZ;
74: eps->balance = EPS_BALANCE_NONE;
75: eps->balance_its = 5;
76: eps->balance_cutoff = 1e-8;
77: eps->trueres = PETSC_FALSE;
78: eps->trackall = PETSC_FALSE;
79: eps->purify = PETSC_TRUE;
80: eps->twosided = PETSC_FALSE;
82: eps->converged = EPSConvergedRelative;
83: eps->convergeduser = NULL;
84: eps->convergeddestroy= NULL;
85: eps->stopping = EPSStoppingBasic;
86: eps->stoppinguser = NULL;
87: eps->stoppingdestroy = NULL;
88: eps->arbitrary = NULL;
89: eps->arbitrarydestroy= NULL;
90: eps->convergedctx = NULL;
91: eps->stoppingctx = NULL;
92: eps->arbitraryctx = NULL;
93: eps->numbermonitors = 0;
95: eps->st = NULL;
96: eps->ds = NULL;
97: eps->V = NULL;
98: eps->W = NULL;
99: eps->rg = NULL;
100: eps->D = NULL;
101: eps->IS = NULL;
102: eps->ISL = NULL;
103: eps->defl = NULL;
104: eps->eigr = NULL;
105: eps->eigi = NULL;
106: eps->errest = NULL;
107: eps->rr = NULL;
108: eps->ri = NULL;
109: eps->perm = NULL;
110: eps->nwork = 0;
111: eps->work = NULL;
112: eps->data = NULL;
114: eps->state = EPS_STATE_INITIAL;
115: eps->categ = EPS_CATEGORY_KRYLOV;
116: eps->nconv = 0;
117: eps->its = 0;
118: eps->nloc = 0;
119: eps->nrma = 0.0;
120: eps->nrmb = 0.0;
121: eps->useds = PETSC_FALSE;
122: eps->isgeneralized = PETSC_FALSE;
123: eps->ispositive = PETSC_FALSE;
124: eps->ishermitian = PETSC_FALSE;
125: eps->isstructured = PETSC_FALSE;
126: eps->reason = EPS_CONVERGED_ITERATING;
128: PetscCall(PetscNew(&eps->sc));
129: *outeps = eps;
130: PetscFunctionReturn(PETSC_SUCCESS);
131: }
133: /*@
134: EPSSetType - Selects the particular solver to be used in the `EPS` object.
136: Logically Collective
138: Input Parameters:
139: + eps - the linear eigensolver context
140: - type - a known method
142: Options Database Key:
143: . -eps_type type - sets the method; use `-help` for a list of available methods
145: Notes:
146: See `EPSType` for available methods. The default is `EPSKRYLOVSCHUR`.
148: Normally, it is best to use the `EPSSetFromOptions()` command and
149: then set the `EPS` type from the options database rather than by using
150: this routine. Using the options database provides the user with
151: maximum flexibility in evaluating the different available methods.
152: The `EPSSetType()` routine is provided for those situations where it
153: is necessary to set the iterative solver independently of the command
154: line or options database.
156: Level: intermediate
158: .seealso: [](ch:eps), `STSetType()`, `EPSType`
159: @*/
160: PetscErrorCode EPSSetType(EPS eps,EPSType type)
161: {
162: PetscErrorCode (*r)(EPS);
163: PetscBool match;
165: PetscFunctionBegin;
167: PetscAssertPointer(type,2);
169: PetscCall(PetscObjectTypeCompare((PetscObject)eps,type,&match));
170: if (match) PetscFunctionReturn(PETSC_SUCCESS);
172: PetscCall(PetscFunctionListFind(EPSList,type,&r));
173: PetscCheck(r,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown EPS type given: %s",type);
175: PetscTryTypeMethod(eps,destroy);
176: PetscCall(PetscMemzero(eps->ops,sizeof(struct _EPSOps)));
178: eps->state = EPS_STATE_INITIAL;
179: PetscCall(PetscObjectChangeTypeName((PetscObject)eps,type));
180: PetscCall((*r)(eps));
181: PetscFunctionReturn(PETSC_SUCCESS);
182: }
184: /*@
185: EPSGetType - Gets the `EPS` type as a string from the `EPS` object.
187: Not Collective
189: Input Parameter:
190: . eps - the linear eigensolver context
192: Output Parameter:
193: . type - name of `EPS` method
195: Level: intermediate
197: .seealso: [](ch:eps), `EPSSetType()`
198: @*/
199: PetscErrorCode EPSGetType(EPS eps,EPSType *type)
200: {
201: PetscFunctionBegin;
203: PetscAssertPointer(type,2);
204: *type = ((PetscObject)eps)->type_name;
205: PetscFunctionReturn(PETSC_SUCCESS);
206: }
208: /*@C
209: EPSRegister - Adds a method to the eigenproblem solver package.
211: Not Collective
213: Input Parameters:
214: + name - name of a new user-defined solver
215: - function - routine to create the solver context
217: Note:
218: `EPSRegister()` may be called multiple times to add several user-defined solvers.
220: Example Usage:
221: .vb
222: EPSRegister("my_solver",MySolverCreate);
223: .ve
225: Then, your solver can be chosen with the procedural interface via
226: .vb
227: EPSSetType(eps,"my_solver")
228: .ve
229: or at runtime via the option `-eps_type my_solver`.
231: Level: advanced
233: .seealso: [](ch:eps), `EPSRegisterAll()`
234: @*/
235: PetscErrorCode EPSRegister(const char *name,PetscErrorCode (*function)(EPS))
236: {
237: PetscFunctionBegin;
238: PetscCall(EPSInitializePackage());
239: PetscCall(PetscFunctionListAdd(&EPSList,name,function));
240: PetscFunctionReturn(PETSC_SUCCESS);
241: }
243: /*@C
244: EPSMonitorRegister - Registers an `EPS` monitor routine that may be accessed with
245: `EPSMonitorSetFromOptions()`.
247: Not Collective
249: Input Parameters:
250: + name - name of a new monitor routine
251: . vtype - a `PetscViewerType` for the output
252: . format - a `PetscViewerFormat` for the output
253: . monitor - monitor routine, see `EPSMonitorRegisterFn`
254: . create - creation routine, or `NULL`
255: - destroy - destruction routine, or `NULL`
257: Notes:
258: `EPSMonitorRegister()` may be called multiple times to add several user-defined monitors.
260: The calling sequence for the given function matches the calling sequence of `EPSMonitorFn`
261: functions passed to `EPSMonitorSet()` with the additional requirement that its final argument
262: be a `PetscViewerAndFormat`.
264: Example Usage:
265: .vb
266: EPSMonitorRegister("my_monitor",PETSCVIEWERASCII,PETSC_VIEWER_ASCII_INFO_DETAIL,MyMonitor,NULL,NULL);
267: .ve
269: Then, your monitor can be chosen with the procedural interface via
270: .vb
271: EPSMonitorSetFromOptions(eps,"-eps_monitor_my_monitor","my_monitor",NULL)
272: .ve
273: or at runtime via the option `-eps_monitor_my_monitor`.
275: Level: advanced
277: .seealso: [](ch:eps), `EPSMonitorSet()`, `EPSMonitorRegisterAll()`, `EPSMonitorSetFromOptions()`
278: @*/
279: PetscErrorCode EPSMonitorRegister(const char name[],PetscViewerType vtype,PetscViewerFormat format,EPSMonitorRegisterFn *monitor,EPSMonitorRegisterCreateFn *create,EPSMonitorRegisterDestroyFn *destroy)
280: {
281: char key[PETSC_MAX_PATH_LEN];
283: PetscFunctionBegin;
284: PetscCall(EPSInitializePackage());
285: PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
286: PetscCall(PetscFunctionListAdd(&EPSMonitorList,key,monitor));
287: if (create) PetscCall(PetscFunctionListAdd(&EPSMonitorCreateList,key,create));
288: if (destroy) PetscCall(PetscFunctionListAdd(&EPSMonitorDestroyList,key,destroy));
289: PetscFunctionReturn(PETSC_SUCCESS);
290: }
292: /*@
293: EPSReset - Resets the `EPS` context to the initial state (prior to setup)
294: and destroys any allocated `Vec`s and `Mat`s.
296: Collective
298: Input Parameter:
299: . eps - the linear eigensolver context
301: Note:
302: This can be used when a problem of different matrix size wants to be solved.
303: All options that have previously been set are preserved, so in a next use
304: the solver configuration is the same, but new sizes for matrices and vectors
305: are allowed.
307: Level: advanced
309: .seealso: [](ch:eps), `EPSDestroy()`
310: @*/
311: PetscErrorCode EPSReset(EPS eps)
312: {
313: PetscFunctionBegin;
315: if (!eps) PetscFunctionReturn(PETSC_SUCCESS);
316: PetscTryTypeMethod(eps,reset);
317: if (eps->st) PetscCall(STReset(eps->st));
318: PetscCall(VecDestroy(&eps->D));
319: PetscCall(BVDestroy(&eps->V));
320: PetscCall(BVDestroy(&eps->W));
321: PetscCall(VecDestroyVecs(eps->nwork,&eps->work));
322: eps->nwork = 0;
323: eps->state = EPS_STATE_INITIAL;
324: PetscFunctionReturn(PETSC_SUCCESS);
325: }
327: /*@
328: EPSDestroy - Destroys the `EPS` context.
330: Collective
332: Input Parameter:
333: . eps - the linear eigensolver context
335: Level: beginner
337: .seealso: [](ch:eps), `EPSCreate()`, `EPSSetUp()`, `EPSSolve()`
338: @*/
339: PetscErrorCode EPSDestroy(EPS *eps)
340: {
341: PetscFunctionBegin;
342: if (!*eps) PetscFunctionReturn(PETSC_SUCCESS);
344: if (--((PetscObject)*eps)->refct > 0) { *eps = NULL; PetscFunctionReturn(PETSC_SUCCESS); }
345: PetscCall(EPSReset(*eps));
346: PetscTryTypeMethod(*eps,destroy);
347: if ((*eps)->eigr) PetscCall(PetscFree4((*eps)->eigr,(*eps)->eigi,(*eps)->errest,(*eps)->perm));
348: if ((*eps)->rr) PetscCall(PetscFree2((*eps)->rr,(*eps)->ri));
349: PetscCall(STDestroy(&(*eps)->st));
350: PetscCall(RGDestroy(&(*eps)->rg));
351: PetscCall(DSDestroy(&(*eps)->ds));
352: PetscCall(PetscFree((*eps)->sc));
353: /* just in case the initial vectors have not been used */
354: PetscCall(SlepcBasisDestroy_Private(&(*eps)->nds,&(*eps)->defl));
355: PetscCall(SlepcBasisDestroy_Private(&(*eps)->nini,&(*eps)->IS));
356: PetscCall(SlepcBasisDestroy_Private(&(*eps)->ninil,&(*eps)->ISL));
357: if ((*eps)->convergeddestroy) PetscCall((*(*eps)->convergeddestroy)(&(*eps)->convergedctx));
358: if ((*eps)->stoppingdestroy) PetscCall((*(*eps)->stoppingdestroy)(&(*eps)->stoppingctx));
359: if ((*eps)->arbitrarydestroy) PetscCall((*(*eps)->arbitrarydestroy)(&(*eps)->arbitraryctx));
360: PetscCall(EPSMonitorCancel(*eps));
361: PetscCall(PetscHeaderDestroy(eps));
362: PetscFunctionReturn(PETSC_SUCCESS);
363: }
365: /*@
366: EPSSetTarget - Sets the value of the target.
368: Logically Collective
370: Input Parameters:
371: + eps - the linear eigensolver context
372: - target - the value of the target
374: Options Database Key:
375: . -eps_target target - the value of the target
377: Notes:
378: The target is a scalar value used to determine the portion of the spectrum
379: of interest. It is used in combination with `EPSSetWhichEigenpairs()`.
381: When PETSc is built with real scalars, it is not possible to specify a
382: complex target.
384: In the case of complex scalars, a complex value can be provided in the
385: command line with `[+/-][realnumber][+/-]realnumberi` with no spaces, e.g.
386: `-eps_target 1.0+2.0i`.
388: Level: intermediate
390: .seealso: [](ch:eps), `EPSGetTarget()`, `EPSSetWhichEigenpairs()`
391: @*/
392: PetscErrorCode EPSSetTarget(EPS eps,PetscScalar target)
393: {
394: PetscFunctionBegin;
397: eps->target = target;
398: if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
399: PetscCall(STSetDefaultShift(eps->st,target));
400: PetscFunctionReturn(PETSC_SUCCESS);
401: }
403: /*@
404: EPSGetTarget - Gets the value of the target.
406: Not Collective
408: Input Parameter:
409: . eps - the linear eigensolver context
411: Output Parameter:
412: . target - the value of the target
414: Note:
415: If the target was not set by the user, then zero is returned.
417: Level: intermediate
419: .seealso: [](ch:eps), `EPSSetTarget()`
420: @*/
421: PetscErrorCode EPSGetTarget(EPS eps,PetscScalar* target)
422: {
423: PetscFunctionBegin;
425: PetscAssertPointer(target,2);
426: *target = eps->target;
427: PetscFunctionReturn(PETSC_SUCCESS);
428: }
430: /*@
431: EPSSetInterval - Defines the computational interval for spectrum slicing.
433: Logically Collective
435: Input Parameters:
436: + eps - the linear eigensolver context
437: . inta - left end of the interval
438: - intb - right end of the interval
440: Options Database Key:
441: . -eps_interval a,b - set $[a,b]$ as the interval of interest
443: Notes:
444: Spectrum slicing is a technique employed for computing all eigenvalues of
445: symmetric eigenproblems in a given interval, see section [](#sec:slice).
446: This function provides the interval to be considered. It must be used in
447: combination with `EPS_ALL`, see `EPSSetWhichEigenpairs()`.
449: A computational interval is also needed when using polynomial filters,
450: see `STFILTER` and section [](#sec:filter).
452: In the command-line option, two values must be provided. For an open interval,
453: one can give an infinite, e.g., `-eps_interval 1.0,inf` or `-eps_interval -inf,1.0`.
454: An open interval in the programmatic interface can be specified with
455: `PETSC_MAX_REAL` and -`PETSC_MAX_REAL`.
457: Level: intermediate
459: .seealso: [](ch:eps), [](#sec:slice), [](#sec:filter), `EPSGetInterval()`, `EPSSetWhichEigenpairs()`, `STFILTER`
460: @*/
461: PetscErrorCode EPSSetInterval(EPS eps,PetscReal inta,PetscReal intb)
462: {
463: PetscFunctionBegin;
467: PetscCheck(inta<intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be inta<intb");
468: if (eps->inta != inta || eps->intb != intb) {
469: eps->inta = inta;
470: eps->intb = intb;
471: eps->state = EPS_STATE_INITIAL;
472: }
473: PetscFunctionReturn(PETSC_SUCCESS);
474: }
476: /*@
477: EPSGetInterval - Gets the computational interval for spectrum slicing.
479: Not Collective
481: Input Parameter:
482: . eps - the linear eigensolver context
484: Output Parameters:
485: + inta - left end of the interval
486: - intb - right end of the interval
488: Level: intermediate
490: Note:
491: If the interval was not set by the user, then zeros are returned.
493: .seealso: [](ch:eps), `EPSSetInterval()`
494: @*/
495: PetscErrorCode EPSGetInterval(EPS eps,PetscReal* inta,PetscReal* intb)
496: {
497: PetscFunctionBegin;
499: if (inta) *inta = eps->inta;
500: if (intb) *intb = eps->intb;
501: PetscFunctionReturn(PETSC_SUCCESS);
502: }
504: /*@
505: EPSSetST - Associates a spectral transformation object to the eigensolver.
507: Collective
509: Input Parameters:
510: + eps - the linear eigensolver context
511: - st - the spectral transformation object
513: Note:
514: Use `EPSGetST()` to retrieve the spectral transformation context at a later time
515: (for example, to free it at the end of the computations).
517: Level: advanced
519: .seealso: [](ch:eps), `EPSGetST()`
520: @*/
521: PetscErrorCode EPSSetST(EPS eps,ST st)
522: {
523: PetscFunctionBegin;
526: PetscCheckSameComm(eps,1,st,2);
527: PetscCall(PetscObjectReference((PetscObject)st));
528: PetscCall(STDestroy(&eps->st));
529: eps->st = st;
530: PetscFunctionReturn(PETSC_SUCCESS);
531: }
533: /*@
534: EPSGetST - Obtain the spectral transformation (`ST`) object associated
535: to the eigensolver object.
537: Not Collective
539: Input Parameter:
540: . eps - the linear eigensolver context
542: Output Parameter:
543: . st - spectral transformation context
545: Level: intermediate
547: .seealso: [](ch:eps), `EPSSetST()`
548: @*/
549: PetscErrorCode EPSGetST(EPS eps,ST *st)
550: {
551: PetscFunctionBegin;
553: PetscAssertPointer(st,2);
554: if (!eps->st) {
555: PetscCall(STCreate(PetscObjectComm((PetscObject)eps),&eps->st));
556: PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->st,(PetscObject)eps,0));
557: PetscCall(PetscObjectSetOptions((PetscObject)eps->st,((PetscObject)eps)->options));
558: }
559: *st = eps->st;
560: PetscFunctionReturn(PETSC_SUCCESS);
561: }
563: /*@
564: EPSSetBV - Associates a basis vectors object to the eigensolver.
566: Collective
568: Input Parameters:
569: + eps - the linear eigensolver context
570: - V - the basis vectors object
572: Level: advanced
574: .seealso: [](ch:eps), `EPSGetBV()`
575: @*/
576: PetscErrorCode EPSSetBV(EPS eps,BV V)
577: {
578: PetscFunctionBegin;
581: PetscCheckSameComm(eps,1,V,2);
582: PetscCall(PetscObjectReference((PetscObject)V));
583: PetscCall(BVDestroy(&eps->V));
584: eps->V = V;
585: PetscFunctionReturn(PETSC_SUCCESS);
586: }
588: /*@
589: EPSGetBV - Obtain the basis vectors object associated to the eigensolver object.
591: Not Collective
593: Input Parameter:
594: . eps - the linear eigensolver context
596: Output Parameter:
597: . V - basis vectors context
599: Level: advanced
601: .seealso: [](ch:eps), `EPSSetBV()`
602: @*/
603: PetscErrorCode EPSGetBV(EPS eps,BV *V)
604: {
605: PetscFunctionBegin;
607: PetscAssertPointer(V,2);
608: if (!eps->V) {
609: PetscCall(BVCreate(PetscObjectComm((PetscObject)eps),&eps->V));
610: PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->V,(PetscObject)eps,0));
611: PetscCall(PetscObjectSetOptions((PetscObject)eps->V,((PetscObject)eps)->options));
612: }
613: *V = eps->V;
614: PetscFunctionReturn(PETSC_SUCCESS);
615: }
617: /*@
618: EPSSetRG - Associates a region object to the eigensolver.
620: Collective
622: Input Parameters:
623: + eps - the linear eigensolver context
624: - rg - the region object
626: Note:
627: Use `EPSGetRG()` to retrieve the region context at a later time (for example,
628: to free it at the end of the computations).
630: Level: advanced
632: .seealso: [](ch:eps), `EPSGetRG()`
633: @*/
634: PetscErrorCode EPSSetRG(EPS eps,RG rg)
635: {
636: PetscFunctionBegin;
638: if (rg) {
640: PetscCheckSameComm(eps,1,rg,2);
641: }
642: PetscCall(PetscObjectReference((PetscObject)rg));
643: PetscCall(RGDestroy(&eps->rg));
644: eps->rg = rg;
645: PetscFunctionReturn(PETSC_SUCCESS);
646: }
648: /*@
649: EPSGetRG - Obtain the region object associated to the eigensolver.
651: Not Collective
653: Input Parameter:
654: . eps - the linear eigensolver context
656: Output Parameter:
657: . rg - region context
659: Level: advanced
661: .seealso: [](ch:eps), `EPSSetRG()`
662: @*/
663: PetscErrorCode EPSGetRG(EPS eps,RG *rg)
664: {
665: PetscFunctionBegin;
667: PetscAssertPointer(rg,2);
668: if (!eps->rg) {
669: PetscCall(RGCreate(PetscObjectComm((PetscObject)eps),&eps->rg));
670: PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->rg,(PetscObject)eps,0));
671: PetscCall(PetscObjectSetOptions((PetscObject)eps->rg,((PetscObject)eps)->options));
672: }
673: *rg = eps->rg;
674: PetscFunctionReturn(PETSC_SUCCESS);
675: }
677: /*@
678: EPSSetDS - Associates a direct solver object to the eigensolver.
680: Collective
682: Input Parameters:
683: + eps - the linear eigensolver context
684: - ds - the direct solver object
686: Note:
687: Use `EPSGetDS()` to retrieve the direct solver context at a later time (for example,
688: to free it at the end of the computations).
690: Level: advanced
692: .seealso: [](ch:eps), `EPSGetDS()`
693: @*/
694: PetscErrorCode EPSSetDS(EPS eps,DS ds)
695: {
696: PetscFunctionBegin;
699: PetscCheckSameComm(eps,1,ds,2);
700: PetscCall(PetscObjectReference((PetscObject)ds));
701: PetscCall(DSDestroy(&eps->ds));
702: eps->ds = ds;
703: PetscFunctionReturn(PETSC_SUCCESS);
704: }
706: /*@
707: EPSGetDS - Obtain the direct solver object associated to the eigensolver object.
709: Not Collective
711: Input Parameter:
712: . eps - the linear eigensolver context
714: Output Parameter:
715: . ds - direct solver context
717: Level: advanced
719: .seealso: [](ch:eps), `EPSSetDS()`
720: @*/
721: PetscErrorCode EPSGetDS(EPS eps,DS *ds)
722: {
723: PetscFunctionBegin;
725: PetscAssertPointer(ds,2);
726: if (!eps->ds) {
727: PetscCall(DSCreate(PetscObjectComm((PetscObject)eps),&eps->ds));
728: PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->ds,(PetscObject)eps,0));
729: PetscCall(PetscObjectSetOptions((PetscObject)eps->ds,((PetscObject)eps)->options));
730: }
731: *ds = eps->ds;
732: PetscFunctionReturn(PETSC_SUCCESS);
733: }
735: /*@
736: EPSIsGeneralized - Ask if the `EPS` object corresponds to a generalized
737: eigenvalue problem.
739: Not Collective
741: Input Parameter:
742: . eps - the linear eigensolver context
744: Output Parameter:
745: . is - `PETSC_TRUE` if the problem is generalized
747: Level: intermediate
749: .seealso: [](ch:eps), `EPSIsHermitian()`, `EPSIsPositive()`, `EPSIsStructured()`
750: @*/
751: PetscErrorCode EPSIsGeneralized(EPS eps,PetscBool* is)
752: {
753: PetscFunctionBegin;
755: PetscAssertPointer(is,2);
756: *is = eps->isgeneralized;
757: PetscFunctionReturn(PETSC_SUCCESS);
758: }
760: /*@
761: EPSIsHermitian - Ask if the `EPS` object corresponds to a Hermitian
762: eigenvalue problem.
764: Not Collective
766: Input Parameter:
767: . eps - the linear eigensolver context
769: Output Parameter:
770: . is - `PETSC_TRUE` if the problem is Hermitian
772: Level: intermediate
774: .seealso: [](ch:eps), `EPSIsGeneralized()`, `EPSIsPositive()`, `EPSIsStructured()`
775: @*/
776: PetscErrorCode EPSIsHermitian(EPS eps,PetscBool* is)
777: {
778: PetscFunctionBegin;
780: PetscAssertPointer(is,2);
781: *is = eps->ishermitian;
782: PetscFunctionReturn(PETSC_SUCCESS);
783: }
785: /*@
786: EPSIsPositive - Ask if the `EPS` object corresponds to an eigenvalue
787: problem type that requires a positive (semi-) definite matrix $B$.
789: Not Collective
791: Input Parameter:
792: . eps - the linear eigensolver context
794: Output Parameter:
795: . is - `PETSC_TRUE` if the problem is positive (semi-) definite
797: Level: intermediate
799: .seealso: [](ch:eps), `EPSIsGeneralized()`, `EPSIsHermitian()`, `EPSIsStructured()`
800: @*/
801: PetscErrorCode EPSIsPositive(EPS eps,PetscBool* is)
802: {
803: PetscFunctionBegin;
805: PetscAssertPointer(is,2);
806: *is = eps->ispositive;
807: PetscFunctionReturn(PETSC_SUCCESS);
808: }
810: /*@
811: EPSIsStructured - Ask if the `EPS` object corresponds to a structured
812: eigenvalue problem.
814: Not Collective
816: Input Parameter:
817: . eps - the linear eigensolver context
819: Output Parameter:
820: . is - `PETSC_TRUE` if the problem is structured
822: Note:
823: The result will be true if the problem type has been set to some
824: structured type such as `EPS_BSE`. This is independent of whether the
825: input matrix has been built with a certain structure with a helper function.
827: Level: intermediate
829: .seealso: [](ch:eps), `EPSIsGeneralized()`, `EPSIsHermitian()`, `EPSIsPositive()`, `EPSSetProblemType()`
830: @*/
831: PetscErrorCode EPSIsStructured(EPS eps,PetscBool* is)
832: {
833: PetscFunctionBegin;
835: PetscAssertPointer(is,2);
836: *is = eps->isstructured;
837: PetscFunctionReturn(PETSC_SUCCESS);
838: }