Actual source code: epsbasic.c
slepc-main 2024-11-15
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 default 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: 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 = 1;
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->conv = EPS_CONV_REL;
66: eps->stop = EPS_STOP_BASIC;
67: eps->which = (EPSWhich)0;
68: eps->inta = 0.0;
69: eps->intb = 0.0;
70: eps->problem_type = (EPSProblemType)0;
71: eps->extraction = EPS_RITZ;
72: eps->balance = EPS_BALANCE_NONE;
73: eps->balance_its = 5;
74: eps->balance_cutoff = 1e-8;
75: eps->trueres = PETSC_FALSE;
76: eps->trackall = PETSC_FALSE;
77: eps->purify = PETSC_TRUE;
78: eps->twosided = PETSC_FALSE;
80: eps->converged = EPSConvergedRelative;
81: eps->convergeduser = NULL;
82: eps->convergeddestroy= NULL;
83: eps->stopping = EPSStoppingBasic;
84: eps->stoppinguser = NULL;
85: eps->stoppingdestroy = NULL;
86: eps->arbitrary = NULL;
87: eps->convergedctx = NULL;
88: eps->stoppingctx = NULL;
89: eps->arbitraryctx = NULL;
90: eps->numbermonitors = 0;
92: eps->st = NULL;
93: eps->ds = NULL;
94: eps->V = NULL;
95: eps->W = NULL;
96: eps->rg = NULL;
97: eps->D = NULL;
98: eps->IS = NULL;
99: eps->ISL = NULL;
100: eps->defl = NULL;
101: eps->eigr = NULL;
102: eps->eigi = NULL;
103: eps->errest = NULL;
104: eps->rr = NULL;
105: eps->ri = NULL;
106: eps->perm = NULL;
107: eps->nwork = 0;
108: eps->work = NULL;
109: eps->data = NULL;
111: eps->state = EPS_STATE_INITIAL;
112: eps->categ = EPS_CATEGORY_KRYLOV;
113: eps->nconv = 0;
114: eps->its = 0;
115: eps->nloc = 0;
116: eps->nrma = 0.0;
117: eps->nrmb = 0.0;
118: eps->useds = PETSC_FALSE;
119: eps->isgeneralized = PETSC_FALSE;
120: eps->ispositive = PETSC_FALSE;
121: eps->ishermitian = PETSC_FALSE;
122: eps->isstructured = PETSC_FALSE;
123: eps->reason = EPS_CONVERGED_ITERATING;
125: PetscCall(PetscNew(&eps->sc));
126: *outeps = eps;
127: PetscFunctionReturn(PETSC_SUCCESS);
128: }
130: /*@
131: EPSSetType - Selects the particular solver to be used in the EPS object.
133: Logically Collective
135: Input Parameters:
136: + eps - the eigensolver context
137: - type - a known method
139: Options Database Key:
140: . -eps_type <method> - Sets the method; use -help for a list
141: of available methods
143: Notes:
144: See "slepc/include/slepceps.h" for available methods. The default
145: is EPSKRYLOVSCHUR.
147: Normally, it is best to use the EPSSetFromOptions() command and
148: then set the EPS type from the options database rather than by using
149: this routine. Using the options database provides the user with
150: maximum flexibility in evaluating the different available methods.
151: The EPSSetType() routine is provided for those situations where it
152: is necessary to set the iterative solver independently of the command
153: line or options database.
155: Level: intermediate
157: .seealso: STSetType(), EPSType
158: @*/
159: PetscErrorCode EPSSetType(EPS eps,EPSType type)
160: {
161: PetscErrorCode (*r)(EPS);
162: PetscBool match;
164: PetscFunctionBegin;
166: PetscAssertPointer(type,2);
168: PetscCall(PetscObjectTypeCompare((PetscObject)eps,type,&match));
169: if (match) PetscFunctionReturn(PETSC_SUCCESS);
171: PetscCall(PetscFunctionListFind(EPSList,type,&r));
172: PetscCheck(r,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown EPS type given: %s",type);
174: PetscTryTypeMethod(eps,destroy);
175: PetscCall(PetscMemzero(eps->ops,sizeof(struct _EPSOps)));
177: eps->state = EPS_STATE_INITIAL;
178: PetscCall(PetscObjectChangeTypeName((PetscObject)eps,type));
179: PetscCall((*r)(eps));
180: PetscFunctionReturn(PETSC_SUCCESS);
181: }
183: /*@
184: EPSGetType - Gets the EPS type as a string from the EPS object.
186: Not Collective
188: Input Parameter:
189: . eps - the eigensolver context
191: Output Parameter:
192: . type - name of EPS method
194: Level: intermediate
196: .seealso: EPSSetType()
197: @*/
198: PetscErrorCode EPSGetType(EPS eps,EPSType *type)
199: {
200: PetscFunctionBegin;
202: PetscAssertPointer(type,2);
203: *type = ((PetscObject)eps)->type_name;
204: PetscFunctionReturn(PETSC_SUCCESS);
205: }
207: /*@C
208: EPSRegister - Adds a method to the eigenproblem solver package.
210: Not Collective
212: Input Parameters:
213: + name - name of a new user-defined solver
214: - function - routine to create the solver context
216: Notes:
217: EPSRegister() may be called multiple times to add several user-defined solvers.
219: Example Usage:
220: .vb
221: EPSRegister("my_solver",MySolverCreate);
222: .ve
224: Then, your solver can be chosen with the procedural interface via
225: $ EPSSetType(eps,"my_solver")
226: or at runtime via the option
227: $ -eps_type my_solver
229: Level: advanced
231: .seealso: EPSRegisterAll()
232: @*/
233: PetscErrorCode EPSRegister(const char *name,PetscErrorCode (*function)(EPS))
234: {
235: PetscFunctionBegin;
236: PetscCall(EPSInitializePackage());
237: PetscCall(PetscFunctionListAdd(&EPSList,name,function));
238: PetscFunctionReturn(PETSC_SUCCESS);
239: }
241: /*@C
242: EPSMonitorRegister - Adds EPS monitor routine.
244: Not Collective
246: Input Parameters:
247: + name - name of a new monitor routine
248: . vtype - a PetscViewerType for the output
249: . format - a PetscViewerFormat for the output
250: . monitor - monitor routine
251: . create - creation routine, or NULL
252: - destroy - destruction routine, or NULL
254: Notes:
255: EPSMonitorRegister() may be called multiple times to add several user-defined monitors.
257: Example Usage:
258: .vb
259: EPSMonitorRegister("my_monitor",PETSCVIEWERASCII,PETSC_VIEWER_ASCII_INFO_DETAIL,MyMonitor,NULL,NULL);
260: .ve
262: Then, your monitor can be chosen with the procedural interface via
263: $ EPSMonitorSetFromOptions(eps,"-eps_monitor_my_monitor","my_monitor",NULL)
264: or at runtime via the option
265: $ -eps_monitor_my_monitor
267: Level: advanced
269: .seealso: EPSMonitorRegisterAll()
270: @*/
271: PetscErrorCode EPSMonitorRegister(const char name[],PetscViewerType vtype,PetscViewerFormat format,PetscErrorCode (*monitor)(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode (*create)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode (*destroy)(PetscViewerAndFormat**))
272: {
273: char key[PETSC_MAX_PATH_LEN];
275: PetscFunctionBegin;
276: PetscCall(EPSInitializePackage());
277: PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
278: PetscCall(PetscFunctionListAdd(&EPSMonitorList,key,monitor));
279: if (create) PetscCall(PetscFunctionListAdd(&EPSMonitorCreateList,key,create));
280: if (destroy) PetscCall(PetscFunctionListAdd(&EPSMonitorDestroyList,key,destroy));
281: PetscFunctionReturn(PETSC_SUCCESS);
282: }
284: /*@
285: EPSReset - Resets the EPS context to the initial state (prior to setup)
286: and destroys any allocated Vecs and Mats.
288: Collective
290: Input Parameter:
291: . eps - eigensolver context obtained from EPSCreate()
293: Note:
294: This can be used when a problem of different matrix size wants to be solved.
295: All options that have previously been set are preserved, so in a next use
296: the solver configuration is the same, but new sizes for matrices and vectors
297: are allowed.
299: Level: advanced
301: .seealso: EPSDestroy()
302: @*/
303: PetscErrorCode EPSReset(EPS eps)
304: {
305: PetscFunctionBegin;
307: if (!eps) PetscFunctionReturn(PETSC_SUCCESS);
308: PetscTryTypeMethod(eps,reset);
309: if (eps->st) PetscCall(STReset(eps->st));
310: PetscCall(VecDestroy(&eps->D));
311: PetscCall(BVDestroy(&eps->V));
312: PetscCall(BVDestroy(&eps->W));
313: PetscCall(VecDestroyVecs(eps->nwork,&eps->work));
314: eps->nwork = 0;
315: eps->state = EPS_STATE_INITIAL;
316: PetscFunctionReturn(PETSC_SUCCESS);
317: }
319: /*@
320: EPSDestroy - Destroys the EPS context.
322: Collective
324: Input Parameter:
325: . eps - eigensolver context obtained from EPSCreate()
327: Level: beginner
329: .seealso: EPSCreate(), EPSSetUp(), EPSSolve()
330: @*/
331: PetscErrorCode EPSDestroy(EPS *eps)
332: {
333: PetscFunctionBegin;
334: if (!*eps) PetscFunctionReturn(PETSC_SUCCESS);
336: if (--((PetscObject)*eps)->refct > 0) { *eps = NULL; PetscFunctionReturn(PETSC_SUCCESS); }
337: PetscCall(EPSReset(*eps));
338: PetscTryTypeMethod(*eps,destroy);
339: if ((*eps)->eigr) PetscCall(PetscFree4((*eps)->eigr,(*eps)->eigi,(*eps)->errest,(*eps)->perm));
340: if ((*eps)->rr) PetscCall(PetscFree2((*eps)->rr,(*eps)->ri));
341: PetscCall(STDestroy(&(*eps)->st));
342: PetscCall(RGDestroy(&(*eps)->rg));
343: PetscCall(DSDestroy(&(*eps)->ds));
344: PetscCall(PetscFree((*eps)->sc));
345: /* just in case the initial vectors have not been used */
346: PetscCall(SlepcBasisDestroy_Private(&(*eps)->nds,&(*eps)->defl));
347: PetscCall(SlepcBasisDestroy_Private(&(*eps)->nini,&(*eps)->IS));
348: PetscCall(SlepcBasisDestroy_Private(&(*eps)->ninil,&(*eps)->ISL));
349: if ((*eps)->convergeddestroy) PetscCall((*(*eps)->convergeddestroy)(&(*eps)->convergedctx));
350: if ((*eps)->stoppingdestroy) PetscCall((*(*eps)->stoppingdestroy)(&(*eps)->stoppingctx));
351: PetscCall(EPSMonitorCancel(*eps));
352: PetscCall(PetscHeaderDestroy(eps));
353: PetscFunctionReturn(PETSC_SUCCESS);
354: }
356: /*@
357: EPSSetTarget - Sets the value of the target.
359: Logically Collective
361: Input Parameters:
362: + eps - eigensolver context
363: - target - the value of the target
365: Options Database Key:
366: . -eps_target <scalar> - the value of the target
368: Notes:
369: The target is a scalar value used to determine the portion of the spectrum
370: of interest. It is used in combination with EPSSetWhichEigenpairs().
372: In the case of complex scalars, a complex value can be provided in the
373: command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
374: -eps_target 1.0+2.0i
376: Level: intermediate
378: .seealso: EPSGetTarget(), EPSSetWhichEigenpairs()
379: @*/
380: PetscErrorCode EPSSetTarget(EPS eps,PetscScalar target)
381: {
382: PetscFunctionBegin;
385: eps->target = target;
386: if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
387: PetscCall(STSetDefaultShift(eps->st,target));
388: PetscFunctionReturn(PETSC_SUCCESS);
389: }
391: /*@
392: EPSGetTarget - Gets the value of the target.
394: Not Collective
396: Input Parameter:
397: . eps - eigensolver context
399: Output Parameter:
400: . target - the value of the target
402: Note:
403: If the target was not set by the user, then zero is returned.
405: Level: intermediate
407: .seealso: EPSSetTarget()
408: @*/
409: PetscErrorCode EPSGetTarget(EPS eps,PetscScalar* target)
410: {
411: PetscFunctionBegin;
413: PetscAssertPointer(target,2);
414: *target = eps->target;
415: PetscFunctionReturn(PETSC_SUCCESS);
416: }
418: /*@
419: EPSSetInterval - Defines the computational interval for spectrum slicing.
421: Logically Collective
423: Input Parameters:
424: + eps - eigensolver context
425: . inta - left end of the interval
426: - intb - right end of the interval
428: Options Database Key:
429: . -eps_interval <a,b> - set [a,b] as the interval of interest
431: Notes:
432: Spectrum slicing is a technique employed for computing all eigenvalues of
433: symmetric eigenproblems in a given interval. This function provides the
434: interval to be considered. It must be used in combination with EPS_ALL, see
435: EPSSetWhichEigenpairs().
437: In the command-line option, two values must be provided. For an open interval,
438: one can give an infinite, e.g., -eps_interval 1.0,inf or -eps_interval -inf,1.0.
439: An open interval in the programmatic interface can be specified with
440: PETSC_MAX_REAL and -PETSC_MAX_REAL.
442: Level: intermediate
444: .seealso: EPSGetInterval(), EPSSetWhichEigenpairs()
445: @*/
446: PetscErrorCode EPSSetInterval(EPS eps,PetscReal inta,PetscReal intb)
447: {
448: PetscFunctionBegin;
452: PetscCheck(inta<intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be inta<intb");
453: if (eps->inta != inta || eps->intb != intb) {
454: eps->inta = inta;
455: eps->intb = intb;
456: eps->state = EPS_STATE_INITIAL;
457: }
458: PetscFunctionReturn(PETSC_SUCCESS);
459: }
461: /*@
462: EPSGetInterval - Gets the computational interval for spectrum slicing.
464: Not Collective
466: Input Parameter:
467: . eps - eigensolver context
469: Output Parameters:
470: + inta - left end of the interval
471: - intb - right end of the interval
473: Level: intermediate
475: Note:
476: If the interval was not set by the user, then zeros are returned.
478: .seealso: EPSSetInterval()
479: @*/
480: PetscErrorCode EPSGetInterval(EPS eps,PetscReal* inta,PetscReal* intb)
481: {
482: PetscFunctionBegin;
484: if (inta) *inta = eps->inta;
485: if (intb) *intb = eps->intb;
486: PetscFunctionReturn(PETSC_SUCCESS);
487: }
489: /*@
490: EPSSetST - Associates a spectral transformation object to the eigensolver.
492: Collective
494: Input Parameters:
495: + eps - eigensolver context obtained from EPSCreate()
496: - st - the spectral transformation object
498: Note:
499: Use EPSGetST() to retrieve the spectral transformation context (for example,
500: to free it at the end of the computations).
502: Level: advanced
504: .seealso: EPSGetST()
505: @*/
506: PetscErrorCode EPSSetST(EPS eps,ST st)
507: {
508: PetscFunctionBegin;
511: PetscCheckSameComm(eps,1,st,2);
512: PetscCall(PetscObjectReference((PetscObject)st));
513: PetscCall(STDestroy(&eps->st));
514: eps->st = st;
515: PetscFunctionReturn(PETSC_SUCCESS);
516: }
518: /*@
519: EPSGetST - Obtain the spectral transformation (ST) object associated
520: to the eigensolver object.
522: Not Collective
524: Input Parameters:
525: . eps - eigensolver context obtained from EPSCreate()
527: Output Parameter:
528: . st - spectral transformation context
530: Level: intermediate
532: .seealso: EPSSetST()
533: @*/
534: PetscErrorCode EPSGetST(EPS eps,ST *st)
535: {
536: PetscFunctionBegin;
538: PetscAssertPointer(st,2);
539: if (!eps->st) {
540: PetscCall(STCreate(PetscObjectComm((PetscObject)eps),&eps->st));
541: PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->st,(PetscObject)eps,0));
542: PetscCall(PetscObjectSetOptions((PetscObject)eps->st,((PetscObject)eps)->options));
543: }
544: *st = eps->st;
545: PetscFunctionReturn(PETSC_SUCCESS);
546: }
548: /*@
549: EPSSetBV - Associates a basis vectors object to the eigensolver.
551: Collective
553: Input Parameters:
554: + eps - eigensolver context obtained from EPSCreate()
555: - V - the basis vectors object
557: Level: advanced
559: .seealso: EPSGetBV()
560: @*/
561: PetscErrorCode EPSSetBV(EPS eps,BV V)
562: {
563: PetscFunctionBegin;
566: PetscCheckSameComm(eps,1,V,2);
567: PetscCall(PetscObjectReference((PetscObject)V));
568: PetscCall(BVDestroy(&eps->V));
569: eps->V = V;
570: PetscFunctionReturn(PETSC_SUCCESS);
571: }
573: /*@
574: EPSGetBV - Obtain the basis vectors object associated to the eigensolver object.
576: Not Collective
578: Input Parameters:
579: . eps - eigensolver context obtained from EPSCreate()
581: Output Parameter:
582: . V - basis vectors context
584: Level: advanced
586: .seealso: EPSSetBV()
587: @*/
588: PetscErrorCode EPSGetBV(EPS eps,BV *V)
589: {
590: PetscFunctionBegin;
592: PetscAssertPointer(V,2);
593: if (!eps->V) {
594: PetscCall(BVCreate(PetscObjectComm((PetscObject)eps),&eps->V));
595: PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->V,(PetscObject)eps,0));
596: PetscCall(PetscObjectSetOptions((PetscObject)eps->V,((PetscObject)eps)->options));
597: }
598: *V = eps->V;
599: PetscFunctionReturn(PETSC_SUCCESS);
600: }
602: /*@
603: EPSSetRG - Associates a region object to the eigensolver.
605: Collective
607: Input Parameters:
608: + eps - eigensolver context obtained from EPSCreate()
609: - rg - the region object
611: Note:
612: Use EPSGetRG() to retrieve the region context (for example,
613: to free it at the end of the computations).
615: Level: advanced
617: .seealso: EPSGetRG()
618: @*/
619: PetscErrorCode EPSSetRG(EPS eps,RG rg)
620: {
621: PetscFunctionBegin;
623: if (rg) {
625: PetscCheckSameComm(eps,1,rg,2);
626: }
627: PetscCall(PetscObjectReference((PetscObject)rg));
628: PetscCall(RGDestroy(&eps->rg));
629: eps->rg = rg;
630: PetscFunctionReturn(PETSC_SUCCESS);
631: }
633: /*@
634: EPSGetRG - Obtain the region object associated to the eigensolver.
636: Not Collective
638: Input Parameters:
639: . eps - eigensolver context obtained from EPSCreate()
641: Output Parameter:
642: . rg - region context
644: Level: advanced
646: .seealso: EPSSetRG()
647: @*/
648: PetscErrorCode EPSGetRG(EPS eps,RG *rg)
649: {
650: PetscFunctionBegin;
652: PetscAssertPointer(rg,2);
653: if (!eps->rg) {
654: PetscCall(RGCreate(PetscObjectComm((PetscObject)eps),&eps->rg));
655: PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->rg,(PetscObject)eps,0));
656: PetscCall(PetscObjectSetOptions((PetscObject)eps->rg,((PetscObject)eps)->options));
657: }
658: *rg = eps->rg;
659: PetscFunctionReturn(PETSC_SUCCESS);
660: }
662: /*@
663: EPSSetDS - Associates a direct solver object to the eigensolver.
665: Collective
667: Input Parameters:
668: + eps - eigensolver context obtained from EPSCreate()
669: - ds - the direct solver object
671: Note:
672: Use EPSGetDS() to retrieve the direct solver context (for example,
673: to free it at the end of the computations).
675: Level: advanced
677: .seealso: EPSGetDS()
678: @*/
679: PetscErrorCode EPSSetDS(EPS eps,DS ds)
680: {
681: PetscFunctionBegin;
684: PetscCheckSameComm(eps,1,ds,2);
685: PetscCall(PetscObjectReference((PetscObject)ds));
686: PetscCall(DSDestroy(&eps->ds));
687: eps->ds = ds;
688: PetscFunctionReturn(PETSC_SUCCESS);
689: }
691: /*@
692: EPSGetDS - Obtain the direct solver object associated to the eigensolver object.
694: Not Collective
696: Input Parameters:
697: . eps - eigensolver context obtained from EPSCreate()
699: Output Parameter:
700: . ds - direct solver context
702: Level: advanced
704: .seealso: EPSSetDS()
705: @*/
706: PetscErrorCode EPSGetDS(EPS eps,DS *ds)
707: {
708: PetscFunctionBegin;
710: PetscAssertPointer(ds,2);
711: if (!eps->ds) {
712: PetscCall(DSCreate(PetscObjectComm((PetscObject)eps),&eps->ds));
713: PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->ds,(PetscObject)eps,0));
714: PetscCall(PetscObjectSetOptions((PetscObject)eps->ds,((PetscObject)eps)->options));
715: }
716: *ds = eps->ds;
717: PetscFunctionReturn(PETSC_SUCCESS);
718: }
720: /*@
721: EPSIsGeneralized - Ask if the EPS object corresponds to a generalized
722: eigenvalue problem.
724: Not Collective
726: Input Parameter:
727: . eps - the eigenproblem solver context
729: Output Parameter:
730: . is - the answer
732: Level: intermediate
734: .seealso: EPSIsHermitian(), EPSIsPositive(), EPSIsStructured()
735: @*/
736: PetscErrorCode EPSIsGeneralized(EPS eps,PetscBool* is)
737: {
738: PetscFunctionBegin;
740: PetscAssertPointer(is,2);
741: *is = eps->isgeneralized;
742: PetscFunctionReturn(PETSC_SUCCESS);
743: }
745: /*@
746: EPSIsHermitian - Ask if the EPS object corresponds to a Hermitian
747: eigenvalue problem.
749: Not Collective
751: Input Parameter:
752: . eps - the eigenproblem solver context
754: Output Parameter:
755: . is - the answer
757: Level: intermediate
759: .seealso: EPSIsGeneralized(), EPSIsPositive(), EPSIsStructured()
760: @*/
761: PetscErrorCode EPSIsHermitian(EPS eps,PetscBool* is)
762: {
763: PetscFunctionBegin;
765: PetscAssertPointer(is,2);
766: *is = eps->ishermitian;
767: PetscFunctionReturn(PETSC_SUCCESS);
768: }
770: /*@
771: EPSIsPositive - Ask if the EPS object corresponds to an eigenvalue
772: problem type that requires a positive (semi-) definite matrix B.
774: Not Collective
776: Input Parameter:
777: . eps - the eigenproblem solver context
779: Output Parameter:
780: . is - the answer
782: Level: intermediate
784: .seealso: EPSIsGeneralized(), EPSIsHermitian(), EPSIsStructured()
785: @*/
786: PetscErrorCode EPSIsPositive(EPS eps,PetscBool* is)
787: {
788: PetscFunctionBegin;
790: PetscAssertPointer(is,2);
791: *is = eps->ispositive;
792: PetscFunctionReturn(PETSC_SUCCESS);
793: }
795: /*@
796: EPSIsStructured - Ask if the EPS object corresponds to a structured
797: eigenvalue problem.
799: Not Collective
801: Input Parameter:
802: . eps - the eigenproblem solver context
804: Output Parameter:
805: . is - the answer
807: Note:
808: The result will be true if the problem type has been set to some
809: structured type such as EPS_BSE. This is independent of whether the
810: input matrix has been built with a certain structure with a helper function.
812: Level: intermediate
814: .seealso: EPSIsGeneralized(), EPSIsHermitian(), EPSIsPositive(), EPSSetProblemType()
815: @*/
816: PetscErrorCode EPSIsStructured(EPS eps,PetscBool* is)
817: {
818: PetscFunctionBegin;
820: PetscAssertPointer(is,2);
821: *is = eps->isstructured;
822: PetscFunctionReturn(PETSC_SUCCESS);
823: }