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:    Note:
196:    `type` should not be retained for later use as it will be an invalid pointer
197:    if the `EPSType` of `eps` is changed.

199:    Level: intermediate

201: .seealso: [](ch:eps), `EPSSetType()`, `PetscObjectTypeCompare()`, `PetscObjectTypeCompareAny()`
202: @*/
203: PetscErrorCode EPSGetType(EPS eps,EPSType *type)
204: {
205:   PetscFunctionBegin;
207:   PetscAssertPointer(type,2);
208:   *type = ((PetscObject)eps)->type_name;
209:   PetscFunctionReturn(PETSC_SUCCESS);
210: }

212: /*@C
213:    EPSRegister - Adds a method to the eigenproblem solver package.

215:    Not Collective

217:    Input Parameters:
218: +  name - name of a new user-defined solver
219: -  function - routine to create the solver context

221:    Note:
222:    `EPSRegister()` may be called multiple times to add several user-defined solvers.

224:    Example Usage:
225: .vb
226:    EPSRegister("my_solver",MySolverCreate);
227: .ve

229:    Then, your solver can be chosen with the procedural interface via
230: .vb
231:    EPSSetType(eps,"my_solver")
232: .ve
233:    or at runtime via the option `-eps_type my_solver`.

235:    Level: advanced

237: .seealso: [](ch:eps), `EPSRegisterAll()`
238: @*/
239: PetscErrorCode EPSRegister(const char *name,PetscErrorCode (*function)(EPS))
240: {
241:   PetscFunctionBegin;
242:   PetscCall(EPSInitializePackage());
243:   PetscCall(PetscFunctionListAdd(&EPSList,name,function));
244:   PetscFunctionReturn(PETSC_SUCCESS);
245: }

247: /*@C
248:    EPSMonitorRegister - Registers an `EPS` monitor routine that may be accessed with
249:    `EPSMonitorSetFromOptions()`.

251:    Not Collective

253:    Input Parameters:
254: +  name    - name of a new monitor routine
255: .  vtype   - a `PetscViewerType` for the output
256: .  format  - a `PetscViewerFormat` for the output
257: .  monitor - monitor routine, see `EPSMonitorRegisterFn`
258: .  create  - creation routine, or `NULL`
259: -  destroy - destruction routine, or `NULL`

261:    Notes:
262:    `EPSMonitorRegister()` may be called multiple times to add several user-defined monitors.

264:    The calling sequence for the given function matches the calling sequence of `EPSMonitorFn`
265:    functions passed to `EPSMonitorSet()` with the additional requirement that its final argument
266:    be a `PetscViewerAndFormat`.

268:    Example Usage:
269: .vb
270:    EPSMonitorRegister("my_monitor",PETSCVIEWERASCII,PETSC_VIEWER_ASCII_INFO_DETAIL,MyMonitor,NULL,NULL);
271: .ve

273:    Then, your monitor can be chosen with the procedural interface via
274: .vb
275:    EPSMonitorSetFromOptions(eps,"-eps_monitor_my_monitor","my_monitor",NULL)
276: .ve
277:    or at runtime via the option `-eps_monitor_my_monitor`.

279:    Level: advanced

281: .seealso: [](ch:eps), `EPSMonitorSet()`, `EPSMonitorRegisterAll()`, `EPSMonitorSetFromOptions()`
282: @*/
283: PetscErrorCode EPSMonitorRegister(const char name[],PetscViewerType vtype,PetscViewerFormat format,EPSMonitorRegisterFn *monitor,EPSMonitorRegisterCreateFn *create,EPSMonitorRegisterDestroyFn *destroy)
284: {
285:   char           key[PETSC_MAX_PATH_LEN];

287:   PetscFunctionBegin;
288:   PetscCall(EPSInitializePackage());
289:   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
290:   PetscCall(PetscFunctionListAdd(&EPSMonitorList,key,monitor));
291:   if (create)  PetscCall(PetscFunctionListAdd(&EPSMonitorCreateList,key,create));
292:   if (destroy) PetscCall(PetscFunctionListAdd(&EPSMonitorDestroyList,key,destroy));
293:   PetscFunctionReturn(PETSC_SUCCESS);
294: }

296: /*@
297:    EPSReset - Resets the `EPS` context to the initial state (prior to setup)
298:    and destroys any allocated `Vec`s and `Mat`s.

300:    Collective

302:    Input Parameter:
303: .  eps - the linear eigensolver context

305:    Note:
306:    This can be used when a problem of different matrix size wants to be solved.
307:    All options that have previously been set are preserved, so in a next use
308:    the solver configuration is the same, but new sizes for matrices and vectors
309:    are allowed.

311:    Level: advanced

313: .seealso: [](ch:eps), `EPSDestroy()`
314: @*/
315: PetscErrorCode EPSReset(EPS eps)
316: {
317:   PetscFunctionBegin;
319:   if (!eps) PetscFunctionReturn(PETSC_SUCCESS);
320:   PetscTryTypeMethod(eps,reset);
321:   if (eps->st) PetscCall(STReset(eps->st));
322:   PetscCall(VecDestroy(&eps->D));
323:   PetscCall(BVDestroy(&eps->V));
324:   PetscCall(BVDestroy(&eps->W));
325:   PetscCall(VecDestroyVecs(eps->nwork,&eps->work));
326:   eps->nwork = 0;
327:   eps->state = EPS_STATE_INITIAL;
328:   PetscFunctionReturn(PETSC_SUCCESS);
329: }

331: /*@
332:    EPSDestroy - Destroys the `EPS` context.

334:    Collective

336:    Input Parameter:
337: .  eps - the linear eigensolver context

339:    Level: beginner

341: .seealso: [](ch:eps), `EPSCreate()`, `EPSSetUp()`, `EPSSolve()`
342: @*/
343: PetscErrorCode EPSDestroy(EPS *eps)
344: {
345:   PetscFunctionBegin;
346:   if (!*eps) PetscFunctionReturn(PETSC_SUCCESS);
348:   if (--((PetscObject)*eps)->refct > 0) { *eps = NULL; PetscFunctionReturn(PETSC_SUCCESS); }
349:   PetscCall(EPSReset(*eps));
350:   PetscTryTypeMethod(*eps,destroy);
351:   if ((*eps)->eigr) PetscCall(PetscFree4((*eps)->eigr,(*eps)->eigi,(*eps)->errest,(*eps)->perm));
352:   if ((*eps)->rr) PetscCall(PetscFree2((*eps)->rr,(*eps)->ri));
353:   PetscCall(STDestroy(&(*eps)->st));
354:   PetscCall(RGDestroy(&(*eps)->rg));
355:   PetscCall(DSDestroy(&(*eps)->ds));
356:   PetscCall(PetscFree((*eps)->sc));
357:   /* just in case the initial vectors have not been used */
358:   PetscCall(SlepcBasisDestroy_Private(&(*eps)->nds,&(*eps)->defl));
359:   PetscCall(SlepcBasisDestroy_Private(&(*eps)->nini,&(*eps)->IS));
360:   PetscCall(SlepcBasisDestroy_Private(&(*eps)->ninil,&(*eps)->ISL));
361:   if ((*eps)->convergeddestroy) PetscCall((*(*eps)->convergeddestroy)(&(*eps)->convergedctx));
362:   if ((*eps)->stoppingdestroy) PetscCall((*(*eps)->stoppingdestroy)(&(*eps)->stoppingctx));
363:   if ((*eps)->arbitrarydestroy) PetscCall((*(*eps)->arbitrarydestroy)(&(*eps)->arbitraryctx));
364:   PetscCall(EPSMonitorCancel(*eps));
365:   PetscCall(PetscHeaderDestroy(eps));
366:   PetscFunctionReturn(PETSC_SUCCESS);
367: }

369: /*@
370:    EPSSetTarget - Sets the value of the target.

372:    Logically Collective

374:    Input Parameters:
375: +  eps    - the linear eigensolver context
376: -  target - the value of the target

378:    Options Database Key:
379: .  -eps_target target - the value of the target

381:    Notes:
382:    The target is a scalar value used to determine the portion of the spectrum
383:    of interest. It is used in combination with `EPSSetWhichEigenpairs()`.

385:    When PETSc is built with real scalars, it is not possible to specify a
386:    complex target.

388:    In the case of complex scalars, a complex value can be provided in the
389:    command line with `[+/-][realnumber][+/-]realnumberi` with no spaces, e.g.
390:    `-eps_target 1.0+2.0i`.

392:    Level: intermediate

394: .seealso: [](ch:eps), `EPSGetTarget()`, `EPSSetWhichEigenpairs()`
395: @*/
396: PetscErrorCode EPSSetTarget(EPS eps,PetscScalar target)
397: {
398:   PetscFunctionBegin;
401:   eps->target = target;
402:   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
403:   PetscCall(STSetDefaultShift(eps->st,target));
404:   PetscFunctionReturn(PETSC_SUCCESS);
405: }

407: /*@
408:    EPSGetTarget - Gets the value of the target.

410:    Not Collective

412:    Input Parameter:
413: .  eps - the linear eigensolver context

415:    Output Parameter:
416: .  target - the value of the target

418:    Note:
419:    If the target was not set by the user, then zero is returned.

421:    Level: intermediate

423: .seealso: [](ch:eps), `EPSSetTarget()`
424: @*/
425: PetscErrorCode EPSGetTarget(EPS eps,PetscScalar* target)
426: {
427:   PetscFunctionBegin;
429:   PetscAssertPointer(target,2);
430:   *target = eps->target;
431:   PetscFunctionReturn(PETSC_SUCCESS);
432: }

434: /*@
435:    EPSSetInterval - Defines the computational interval for spectrum slicing.

437:    Logically Collective

439:    Input Parameters:
440: +  eps  - the linear eigensolver context
441: .  inta - left end of the interval
442: -  intb - right end of the interval

444:    Options Database Key:
445: .  -eps_interval a,b - set $[a,b]$ as the interval of interest

447:    Notes:
448:    Spectrum slicing is a technique employed for computing all eigenvalues of
449:    symmetric eigenproblems in a given interval, see section [](#sec:slice).
450:    This function provides the interval to be considered. It must be used in
451:    combination with `EPS_ALL`, see `EPSSetWhichEigenpairs()`.

453:    A computational interval is also needed when using polynomial filters,
454:    see `STFILTER` and section [](#sec:filter).

456:    In the command-line option, two values must be provided. For an open interval,
457:    one can give an infinite, e.g., `-eps_interval 1.0,inf` or `-eps_interval -inf,1.0`.
458:    An open interval in the programmatic interface can be specified with
459:    `PETSC_MAX_REAL` and -`PETSC_MAX_REAL`.

461:    Level: intermediate

463: .seealso: [](ch:eps), [](#sec:slice), [](#sec:filter), `EPSGetInterval()`, `EPSSetWhichEigenpairs()`, `STFILTER`
464: @*/
465: PetscErrorCode EPSSetInterval(EPS eps,PetscReal inta,PetscReal intb)
466: {
467:   PetscFunctionBegin;
471:   PetscCheck(inta<intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be inta<intb");
472:   if (eps->inta != inta || eps->intb != intb) {
473:     eps->inta = inta;
474:     eps->intb = intb;
475:     eps->state = EPS_STATE_INITIAL;
476:   }
477:   PetscFunctionReturn(PETSC_SUCCESS);
478: }

480: /*@
481:    EPSGetInterval - Gets the computational interval for spectrum slicing.

483:    Not Collective

485:    Input Parameter:
486: .  eps - the linear eigensolver context

488:    Output Parameters:
489: +  inta - left end of the interval
490: -  intb - right end of the interval

492:    Level: intermediate

494:    Note:
495:    If the interval was not set by the user, then zeros are returned.

497: .seealso: [](ch:eps), `EPSSetInterval()`
498: @*/
499: PetscErrorCode EPSGetInterval(EPS eps,PetscReal* inta,PetscReal* intb)
500: {
501:   PetscFunctionBegin;
503:   if (inta) *inta = eps->inta;
504:   if (intb) *intb = eps->intb;
505:   PetscFunctionReturn(PETSC_SUCCESS);
506: }

508: /*@
509:    EPSSetST - Associates a spectral transformation object to the eigensolver.

511:    Collective

513:    Input Parameters:
514: +  eps - the linear eigensolver context
515: -  st  - the spectral transformation object

517:    Note:
518:    Use `EPSGetST()` to retrieve the spectral transformation context at a later time
519:    (for example, to free it at the end of the computations).

521:    Level: advanced

523: .seealso: [](ch:eps), `EPSGetST()`
524: @*/
525: PetscErrorCode EPSSetST(EPS eps,ST st)
526: {
527:   PetscFunctionBegin;
530:   PetscCheckSameComm(eps,1,st,2);
531:   PetscCall(PetscObjectReference((PetscObject)st));
532:   PetscCall(STDestroy(&eps->st));
533:   eps->st = st;
534:   PetscFunctionReturn(PETSC_SUCCESS);
535: }

537: /*@
538:    EPSGetST - Obtain the spectral transformation (`ST`) object associated
539:    to the eigensolver object.

541:    Not Collective

543:    Input Parameter:
544: .  eps - the linear eigensolver context

546:    Output Parameter:
547: .  st - spectral transformation context

549:    Level: intermediate

551: .seealso: [](ch:eps), `EPSSetST()`
552: @*/
553: PetscErrorCode EPSGetST(EPS eps,ST *st)
554: {
555:   PetscFunctionBegin;
557:   PetscAssertPointer(st,2);
558:   if (!eps->st) {
559:     PetscCall(STCreate(PetscObjectComm((PetscObject)eps),&eps->st));
560:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->st,(PetscObject)eps,0));
561:     PetscCall(PetscObjectSetOptions((PetscObject)eps->st,((PetscObject)eps)->options));
562:   }
563:   *st = eps->st;
564:   PetscFunctionReturn(PETSC_SUCCESS);
565: }

567: /*@
568:    EPSSetBV - Associates a basis vectors object to the eigensolver.

570:    Collective

572:    Input Parameters:
573: +  eps - the linear eigensolver context
574: -  V   - the basis vectors object

576:    Level: advanced

578: .seealso: [](ch:eps), `EPSGetBV()`
579: @*/
580: PetscErrorCode EPSSetBV(EPS eps,BV V)
581: {
582:   PetscFunctionBegin;
585:   PetscCheckSameComm(eps,1,V,2);
586:   PetscCall(PetscObjectReference((PetscObject)V));
587:   PetscCall(BVDestroy(&eps->V));
588:   eps->V = V;
589:   PetscFunctionReturn(PETSC_SUCCESS);
590: }

592: /*@
593:    EPSGetBV - Obtain the basis vectors object associated to the eigensolver object.

595:    Not Collective

597:    Input Parameter:
598: .  eps - the linear eigensolver context

600:    Output Parameter:
601: .  V - basis vectors context

603:    Level: advanced

605: .seealso: [](ch:eps), `EPSSetBV()`
606: @*/
607: PetscErrorCode EPSGetBV(EPS eps,BV *V)
608: {
609:   PetscFunctionBegin;
611:   PetscAssertPointer(V,2);
612:   if (!eps->V) {
613:     PetscCall(BVCreate(PetscObjectComm((PetscObject)eps),&eps->V));
614:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->V,(PetscObject)eps,0));
615:     PetscCall(PetscObjectSetOptions((PetscObject)eps->V,((PetscObject)eps)->options));
616:   }
617:   *V = eps->V;
618:   PetscFunctionReturn(PETSC_SUCCESS);
619: }

621: /*@
622:    EPSSetRG - Associates a region object to the eigensolver.

624:    Collective

626:    Input Parameters:
627: +  eps - the linear eigensolver context
628: -  rg  - the region object

630:    Note:
631:    Use `EPSGetRG()` to retrieve the region context at a later time (for example,
632:    to free it at the end of the computations).

634:    Level: advanced

636: .seealso: [](ch:eps), `EPSGetRG()`
637: @*/
638: PetscErrorCode EPSSetRG(EPS eps,RG rg)
639: {
640:   PetscFunctionBegin;
642:   if (rg) {
644:     PetscCheckSameComm(eps,1,rg,2);
645:   }
646:   PetscCall(PetscObjectReference((PetscObject)rg));
647:   PetscCall(RGDestroy(&eps->rg));
648:   eps->rg = rg;
649:   PetscFunctionReturn(PETSC_SUCCESS);
650: }

652: /*@
653:    EPSGetRG - Obtain the region object associated to the eigensolver.

655:    Not Collective

657:    Input Parameter:
658: .  eps - the linear eigensolver context

660:    Output Parameter:
661: .  rg - region context

663:    Level: advanced

665: .seealso: [](ch:eps), `EPSSetRG()`
666: @*/
667: PetscErrorCode EPSGetRG(EPS eps,RG *rg)
668: {
669:   PetscFunctionBegin;
671:   PetscAssertPointer(rg,2);
672:   if (!eps->rg) {
673:     PetscCall(RGCreate(PetscObjectComm((PetscObject)eps),&eps->rg));
674:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->rg,(PetscObject)eps,0));
675:     PetscCall(PetscObjectSetOptions((PetscObject)eps->rg,((PetscObject)eps)->options));
676:   }
677:   *rg = eps->rg;
678:   PetscFunctionReturn(PETSC_SUCCESS);
679: }

681: /*@
682:    EPSSetDS - Associates a direct solver object to the eigensolver.

684:    Collective

686:    Input Parameters:
687: +  eps - the linear eigensolver context
688: -  ds  - the direct solver object

690:    Note:
691:    Use `EPSGetDS()` to retrieve the direct solver context at a later time (for example,
692:    to free it at the end of the computations).

694:    Level: advanced

696: .seealso: [](ch:eps), `EPSGetDS()`
697: @*/
698: PetscErrorCode EPSSetDS(EPS eps,DS ds)
699: {
700:   PetscFunctionBegin;
703:   PetscCheckSameComm(eps,1,ds,2);
704:   PetscCall(PetscObjectReference((PetscObject)ds));
705:   PetscCall(DSDestroy(&eps->ds));
706:   eps->ds = ds;
707:   PetscFunctionReturn(PETSC_SUCCESS);
708: }

710: /*@
711:    EPSGetDS - Obtain the direct solver object associated to the eigensolver object.

713:    Not Collective

715:    Input Parameter:
716: .  eps - the linear eigensolver context

718:    Output Parameter:
719: .  ds - direct solver context

721:    Level: advanced

723: .seealso: [](ch:eps), `EPSSetDS()`
724: @*/
725: PetscErrorCode EPSGetDS(EPS eps,DS *ds)
726: {
727:   PetscFunctionBegin;
729:   PetscAssertPointer(ds,2);
730:   if (!eps->ds) {
731:     PetscCall(DSCreate(PetscObjectComm((PetscObject)eps),&eps->ds));
732:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->ds,(PetscObject)eps,0));
733:     PetscCall(PetscObjectSetOptions((PetscObject)eps->ds,((PetscObject)eps)->options));
734:   }
735:   *ds = eps->ds;
736:   PetscFunctionReturn(PETSC_SUCCESS);
737: }

739: /*@
740:    EPSIsGeneralized - Ask if the `EPS` object corresponds to a generalized
741:    eigenvalue problem.

743:    Not Collective

745:    Input Parameter:
746: .  eps - the linear eigensolver context

748:    Output Parameter:
749: .  is - `PETSC_TRUE` if the problem is generalized

751:    Level: intermediate

753: .seealso: [](ch:eps), `EPSIsHermitian()`, `EPSIsPositive()`, `EPSIsStructured()`
754: @*/
755: PetscErrorCode EPSIsGeneralized(EPS eps,PetscBool* is)
756: {
757:   PetscFunctionBegin;
759:   PetscAssertPointer(is,2);
760:   *is = eps->isgeneralized;
761:   PetscFunctionReturn(PETSC_SUCCESS);
762: }

764: /*@
765:    EPSIsHermitian - Ask if the `EPS` object corresponds to a Hermitian
766:    eigenvalue problem.

768:    Not Collective

770:    Input Parameter:
771: .  eps - the linear eigensolver context

773:    Output Parameter:
774: .  is - `PETSC_TRUE` if the problem is Hermitian

776:    Level: intermediate

778: .seealso: [](ch:eps), `EPSIsGeneralized()`, `EPSIsPositive()`, `EPSIsStructured()`
779: @*/
780: PetscErrorCode EPSIsHermitian(EPS eps,PetscBool* is)
781: {
782:   PetscFunctionBegin;
784:   PetscAssertPointer(is,2);
785:   *is = eps->ishermitian;
786:   PetscFunctionReturn(PETSC_SUCCESS);
787: }

789: /*@
790:    EPSIsPositive - Ask if the `EPS` object corresponds to an eigenvalue
791:    problem type that requires a positive (semi-) definite matrix $B$.

793:    Not Collective

795:    Input Parameter:
796: .  eps - the linear eigensolver context

798:    Output Parameter:
799: .  is - `PETSC_TRUE` if the problem is positive (semi-) definite

801:    Level: intermediate

803: .seealso: [](ch:eps), `EPSIsGeneralized()`, `EPSIsHermitian()`, `EPSIsStructured()`
804: @*/
805: PetscErrorCode EPSIsPositive(EPS eps,PetscBool* is)
806: {
807:   PetscFunctionBegin;
809:   PetscAssertPointer(is,2);
810:   *is = eps->ispositive;
811:   PetscFunctionReturn(PETSC_SUCCESS);
812: }

814: /*@
815:    EPSIsStructured - Ask if the `EPS` object corresponds to a structured
816:    eigenvalue problem.

818:    Not Collective

820:    Input Parameter:
821: .  eps - the linear eigensolver context

823:    Output Parameter:
824: .  is - `PETSC_TRUE` if the problem is structured

826:    Note:
827:    The result will be true if the problem type has been set to some
828:    structured type such as `EPS_BSE`. This is independent of whether the
829:    input matrix has been built with a certain structure with a helper function.

831:    Level: intermediate

833: .seealso: [](ch:eps), `EPSIsGeneralized()`, `EPSIsHermitian()`, `EPSIsPositive()`, `EPSSetProblemType()`
834: @*/
835: PetscErrorCode EPSIsStructured(EPS eps,PetscBool* is)
836: {
837:   PetscFunctionBegin;
839:   PetscAssertPointer(is,2);
840:   *is = eps->isstructured;
841:   PetscFunctionReturn(PETSC_SUCCESS);
842: }