Actual source code: epsbasic.c

slepc-main 2024-12-17
Report Typos and Errors
  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             = 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->convergedctx    = NULL;
 90:   eps->stoppingctx     = NULL;
 91:   eps->arbitraryctx    = NULL;
 92:   eps->numbermonitors  = 0;

 94:   eps->st              = NULL;
 95:   eps->ds              = NULL;
 96:   eps->V               = NULL;
 97:   eps->W               = NULL;
 98:   eps->rg              = NULL;
 99:   eps->D               = NULL;
100:   eps->IS              = NULL;
101:   eps->ISL             = NULL;
102:   eps->defl            = NULL;
103:   eps->eigr            = NULL;
104:   eps->eigi            = NULL;
105:   eps->errest          = NULL;
106:   eps->rr              = NULL;
107:   eps->ri              = NULL;
108:   eps->perm            = NULL;
109:   eps->nwork           = 0;
110:   eps->work            = NULL;
111:   eps->data            = NULL;

113:   eps->state           = EPS_STATE_INITIAL;
114:   eps->categ           = EPS_CATEGORY_KRYLOV;
115:   eps->nconv           = 0;
116:   eps->its             = 0;
117:   eps->nloc            = 0;
118:   eps->nrma            = 0.0;
119:   eps->nrmb            = 0.0;
120:   eps->useds           = PETSC_FALSE;
121:   eps->isgeneralized   = PETSC_FALSE;
122:   eps->ispositive      = PETSC_FALSE;
123:   eps->ishermitian     = PETSC_FALSE;
124:   eps->isstructured    = PETSC_FALSE;
125:   eps->reason          = EPS_CONVERGED_ITERATING;

127:   PetscCall(PetscNew(&eps->sc));
128:   *outeps = eps;
129:   PetscFunctionReturn(PETSC_SUCCESS);
130: }

132: /*@
133:    EPSSetType - Selects the particular solver to be used in the EPS object.

135:    Logically Collective

137:    Input Parameters:
138: +  eps  - the eigensolver context
139: -  type - a known method

141:    Options Database Key:
142: .  -eps_type <method> - Sets the method; use -help for a list
143:     of available methods

145:    Notes:
146:    See "slepc/include/slepceps.h" for available methods. The default
147:    is EPSKRYLOVSCHUR.

149:    Normally, it is best to use the EPSSetFromOptions() command and
150:    then set the EPS type from the options database rather than by using
151:    this routine.  Using the options database provides the user with
152:    maximum flexibility in evaluating the different available methods.
153:    The EPSSetType() routine is provided for those situations where it
154:    is necessary to set the iterative solver independently of the command
155:    line or options database.

157:    Level: intermediate

159: .seealso: STSetType(), EPSType
160: @*/
161: PetscErrorCode EPSSetType(EPS eps,EPSType type)
162: {
163:   PetscErrorCode (*r)(EPS);
164:   PetscBool      match;

166:   PetscFunctionBegin;
168:   PetscAssertPointer(type,2);

170:   PetscCall(PetscObjectTypeCompare((PetscObject)eps,type,&match));
171:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

173:   PetscCall(PetscFunctionListFind(EPSList,type,&r));
174:   PetscCheck(r,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown EPS type given: %s",type);

176:   PetscTryTypeMethod(eps,destroy);
177:   PetscCall(PetscMemzero(eps->ops,sizeof(struct _EPSOps)));

179:   eps->state = EPS_STATE_INITIAL;
180:   PetscCall(PetscObjectChangeTypeName((PetscObject)eps,type));
181:   PetscCall((*r)(eps));
182:   PetscFunctionReturn(PETSC_SUCCESS);
183: }

185: /*@
186:    EPSGetType - Gets the EPS type as a string from the EPS object.

188:    Not Collective

190:    Input Parameter:
191: .  eps - the eigensolver context

193:    Output Parameter:
194: .  type - name of EPS method

196:    Level: intermediate

198: .seealso: EPSSetType()
199: @*/
200: PetscErrorCode EPSGetType(EPS eps,EPSType *type)
201: {
202:   PetscFunctionBegin;
204:   PetscAssertPointer(type,2);
205:   *type = ((PetscObject)eps)->type_name;
206:   PetscFunctionReturn(PETSC_SUCCESS);
207: }

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

212:    Not Collective

214:    Input Parameters:
215: +  name - name of a new user-defined solver
216: -  function - routine to create the solver context

218:    Notes:
219:    EPSRegister() may be called multiple times to add several user-defined solvers.

221:    Example Usage:
222: .vb
223:     EPSRegister("my_solver",MySolverCreate);
224: .ve

226:    Then, your solver can be chosen with the procedural interface via
227: $     EPSSetType(eps,"my_solver")
228:    or at runtime via the option
229: $     -eps_type my_solver

231:    Level: advanced

233: .seealso: 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 - Adds EPS monitor routine.

246:    Not Collective

248:    Input Parameters:
249: +  name    - name of a new monitor routine
250: .  vtype   - a PetscViewerType for the output
251: .  format  - a PetscViewerFormat for the output
252: .  monitor - monitor routine
253: .  create  - creation routine, or NULL
254: -  destroy - destruction routine, or NULL

256:    Notes:
257:    EPSMonitorRegister() may be called multiple times to add several user-defined monitors.

259:    Example Usage:
260: .vb
261:    EPSMonitorRegister("my_monitor",PETSCVIEWERASCII,PETSC_VIEWER_ASCII_INFO_DETAIL,MyMonitor,NULL,NULL);
262: .ve

264:    Then, your monitor can be chosen with the procedural interface via
265: $      EPSMonitorSetFromOptions(eps,"-eps_monitor_my_monitor","my_monitor",NULL)
266:    or at runtime via the option
267: $      -eps_monitor_my_monitor

269:    Level: advanced

271: .seealso: EPSMonitorRegisterAll()
272: @*/
273: 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**))
274: {
275:   char           key[PETSC_MAX_PATH_LEN];

277:   PetscFunctionBegin;
278:   PetscCall(EPSInitializePackage());
279:   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
280:   PetscCall(PetscFunctionListAdd(&EPSMonitorList,key,monitor));
281:   if (create)  PetscCall(PetscFunctionListAdd(&EPSMonitorCreateList,key,create));
282:   if (destroy) PetscCall(PetscFunctionListAdd(&EPSMonitorDestroyList,key,destroy));
283:   PetscFunctionReturn(PETSC_SUCCESS);
284: }

286: /*@
287:    EPSReset - Resets the EPS context to the initial state (prior to setup)
288:    and destroys any allocated Vecs and Mats.

290:    Collective

292:    Input Parameter:
293: .  eps - eigensolver context obtained from EPSCreate()

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

301:    Level: advanced

303: .seealso: EPSDestroy()
304: @*/
305: PetscErrorCode EPSReset(EPS eps)
306: {
307:   PetscFunctionBegin;
309:   if (!eps) PetscFunctionReturn(PETSC_SUCCESS);
310:   PetscTryTypeMethod(eps,reset);
311:   if (eps->st) PetscCall(STReset(eps->st));
312:   PetscCall(VecDestroy(&eps->D));
313:   PetscCall(BVDestroy(&eps->V));
314:   PetscCall(BVDestroy(&eps->W));
315:   PetscCall(VecDestroyVecs(eps->nwork,&eps->work));
316:   eps->nwork = 0;
317:   eps->state = EPS_STATE_INITIAL;
318:   PetscFunctionReturn(PETSC_SUCCESS);
319: }

321: /*@
322:    EPSDestroy - Destroys the EPS context.

324:    Collective

326:    Input Parameter:
327: .  eps - eigensolver context obtained from EPSCreate()

329:    Level: beginner

331: .seealso: EPSCreate(), EPSSetUp(), EPSSolve()
332: @*/
333: PetscErrorCode EPSDestroy(EPS *eps)
334: {
335:   PetscFunctionBegin;
336:   if (!*eps) PetscFunctionReturn(PETSC_SUCCESS);
338:   if (--((PetscObject)*eps)->refct > 0) { *eps = NULL; PetscFunctionReturn(PETSC_SUCCESS); }
339:   PetscCall(EPSReset(*eps));
340:   PetscTryTypeMethod(*eps,destroy);
341:   if ((*eps)->eigr) PetscCall(PetscFree4((*eps)->eigr,(*eps)->eigi,(*eps)->errest,(*eps)->perm));
342:   if ((*eps)->rr) PetscCall(PetscFree2((*eps)->rr,(*eps)->ri));
343:   PetscCall(STDestroy(&(*eps)->st));
344:   PetscCall(RGDestroy(&(*eps)->rg));
345:   PetscCall(DSDestroy(&(*eps)->ds));
346:   PetscCall(PetscFree((*eps)->sc));
347:   /* just in case the initial vectors have not been used */
348:   PetscCall(SlepcBasisDestroy_Private(&(*eps)->nds,&(*eps)->defl));
349:   PetscCall(SlepcBasisDestroy_Private(&(*eps)->nini,&(*eps)->IS));
350:   PetscCall(SlepcBasisDestroy_Private(&(*eps)->ninil,&(*eps)->ISL));
351:   if ((*eps)->convergeddestroy) PetscCall((*(*eps)->convergeddestroy)(&(*eps)->convergedctx));
352:   if ((*eps)->stoppingdestroy) PetscCall((*(*eps)->stoppingdestroy)(&(*eps)->stoppingctx));
353:   PetscCall(EPSMonitorCancel(*eps));
354:   PetscCall(PetscHeaderDestroy(eps));
355:   PetscFunctionReturn(PETSC_SUCCESS);
356: }

358: /*@
359:    EPSSetTarget - Sets the value of the target.

361:    Logically Collective

363:    Input Parameters:
364: +  eps    - eigensolver context
365: -  target - the value of the target

367:    Options Database Key:
368: .  -eps_target <scalar> - the value of the target

370:    Notes:
371:    The target is a scalar value used to determine the portion of the spectrum
372:    of interest. It is used in combination with EPSSetWhichEigenpairs().

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

378:    Level: intermediate

380: .seealso: EPSGetTarget(), EPSSetWhichEigenpairs()
381: @*/
382: PetscErrorCode EPSSetTarget(EPS eps,PetscScalar target)
383: {
384:   PetscFunctionBegin;
387:   eps->target = target;
388:   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
389:   PetscCall(STSetDefaultShift(eps->st,target));
390:   PetscFunctionReturn(PETSC_SUCCESS);
391: }

393: /*@
394:    EPSGetTarget - Gets the value of the target.

396:    Not Collective

398:    Input Parameter:
399: .  eps - eigensolver context

401:    Output Parameter:
402: .  target - the value of the target

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

407:    Level: intermediate

409: .seealso: EPSSetTarget()
410: @*/
411: PetscErrorCode EPSGetTarget(EPS eps,PetscScalar* target)
412: {
413:   PetscFunctionBegin;
415:   PetscAssertPointer(target,2);
416:   *target = eps->target;
417:   PetscFunctionReturn(PETSC_SUCCESS);
418: }

420: /*@
421:    EPSSetInterval - Defines the computational interval for spectrum slicing.

423:    Logically Collective

425:    Input Parameters:
426: +  eps  - eigensolver context
427: .  inta - left end of the interval
428: -  intb - right end of the interval

430:    Options Database Key:
431: .  -eps_interval <a,b> - set [a,b] as the interval of interest

433:    Notes:
434:    Spectrum slicing is a technique employed for computing all eigenvalues of
435:    symmetric eigenproblems in a given interval. This function provides the
436:    interval to be considered. It must be used in combination with EPS_ALL, see
437:    EPSSetWhichEigenpairs().

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

444:    Level: intermediate

446: .seealso: EPSGetInterval(), EPSSetWhichEigenpairs()
447: @*/
448: PetscErrorCode EPSSetInterval(EPS eps,PetscReal inta,PetscReal intb)
449: {
450:   PetscFunctionBegin;
454:   PetscCheck(inta<intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be inta<intb");
455:   if (eps->inta != inta || eps->intb != intb) {
456:     eps->inta = inta;
457:     eps->intb = intb;
458:     eps->state = EPS_STATE_INITIAL;
459:   }
460:   PetscFunctionReturn(PETSC_SUCCESS);
461: }

463: /*@
464:    EPSGetInterval - Gets the computational interval for spectrum slicing.

466:    Not Collective

468:    Input Parameter:
469: .  eps - eigensolver context

471:    Output Parameters:
472: +  inta - left end of the interval
473: -  intb - right end of the interval

475:    Level: intermediate

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

480: .seealso: EPSSetInterval()
481: @*/
482: PetscErrorCode EPSGetInterval(EPS eps,PetscReal* inta,PetscReal* intb)
483: {
484:   PetscFunctionBegin;
486:   if (inta) *inta = eps->inta;
487:   if (intb) *intb = eps->intb;
488:   PetscFunctionReturn(PETSC_SUCCESS);
489: }

491: /*@
492:    EPSSetST - Associates a spectral transformation object to the eigensolver.

494:    Collective

496:    Input Parameters:
497: +  eps - eigensolver context obtained from EPSCreate()
498: -  st   - the spectral transformation object

500:    Note:
501:    Use EPSGetST() to retrieve the spectral transformation context (for example,
502:    to free it at the end of the computations).

504:    Level: advanced

506: .seealso: EPSGetST()
507: @*/
508: PetscErrorCode EPSSetST(EPS eps,ST st)
509: {
510:   PetscFunctionBegin;
513:   PetscCheckSameComm(eps,1,st,2);
514:   PetscCall(PetscObjectReference((PetscObject)st));
515:   PetscCall(STDestroy(&eps->st));
516:   eps->st = st;
517:   PetscFunctionReturn(PETSC_SUCCESS);
518: }

520: /*@
521:    EPSGetST - Obtain the spectral transformation (ST) object associated
522:    to the eigensolver object.

524:    Not Collective

526:    Input Parameters:
527: .  eps - eigensolver context obtained from EPSCreate()

529:    Output Parameter:
530: .  st - spectral transformation context

532:    Level: intermediate

534: .seealso: EPSSetST()
535: @*/
536: PetscErrorCode EPSGetST(EPS eps,ST *st)
537: {
538:   PetscFunctionBegin;
540:   PetscAssertPointer(st,2);
541:   if (!eps->st) {
542:     PetscCall(STCreate(PetscObjectComm((PetscObject)eps),&eps->st));
543:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->st,(PetscObject)eps,0));
544:     PetscCall(PetscObjectSetOptions((PetscObject)eps->st,((PetscObject)eps)->options));
545:   }
546:   *st = eps->st;
547:   PetscFunctionReturn(PETSC_SUCCESS);
548: }

550: /*@
551:    EPSSetBV - Associates a basis vectors object to the eigensolver.

553:    Collective

555:    Input Parameters:
556: +  eps - eigensolver context obtained from EPSCreate()
557: -  V   - the basis vectors object

559:    Level: advanced

561: .seealso: EPSGetBV()
562: @*/
563: PetscErrorCode EPSSetBV(EPS eps,BV V)
564: {
565:   PetscFunctionBegin;
568:   PetscCheckSameComm(eps,1,V,2);
569:   PetscCall(PetscObjectReference((PetscObject)V));
570:   PetscCall(BVDestroy(&eps->V));
571:   eps->V = V;
572:   PetscFunctionReturn(PETSC_SUCCESS);
573: }

575: /*@
576:    EPSGetBV - Obtain the basis vectors object associated to the eigensolver object.

578:    Not Collective

580:    Input Parameters:
581: .  eps - eigensolver context obtained from EPSCreate()

583:    Output Parameter:
584: .  V - basis vectors context

586:    Level: advanced

588: .seealso: EPSSetBV()
589: @*/
590: PetscErrorCode EPSGetBV(EPS eps,BV *V)
591: {
592:   PetscFunctionBegin;
594:   PetscAssertPointer(V,2);
595:   if (!eps->V) {
596:     PetscCall(BVCreate(PetscObjectComm((PetscObject)eps),&eps->V));
597:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->V,(PetscObject)eps,0));
598:     PetscCall(PetscObjectSetOptions((PetscObject)eps->V,((PetscObject)eps)->options));
599:   }
600:   *V = eps->V;
601:   PetscFunctionReturn(PETSC_SUCCESS);
602: }

604: /*@
605:    EPSSetRG - Associates a region object to the eigensolver.

607:    Collective

609:    Input Parameters:
610: +  eps - eigensolver context obtained from EPSCreate()
611: -  rg  - the region object

613:    Note:
614:    Use EPSGetRG() to retrieve the region context (for example,
615:    to free it at the end of the computations).

617:    Level: advanced

619: .seealso: EPSGetRG()
620: @*/
621: PetscErrorCode EPSSetRG(EPS eps,RG rg)
622: {
623:   PetscFunctionBegin;
625:   if (rg) {
627:     PetscCheckSameComm(eps,1,rg,2);
628:   }
629:   PetscCall(PetscObjectReference((PetscObject)rg));
630:   PetscCall(RGDestroy(&eps->rg));
631:   eps->rg = rg;
632:   PetscFunctionReturn(PETSC_SUCCESS);
633: }

635: /*@
636:    EPSGetRG - Obtain the region object associated to the eigensolver.

638:    Not Collective

640:    Input Parameters:
641: .  eps - eigensolver context obtained from EPSCreate()

643:    Output Parameter:
644: .  rg - region context

646:    Level: advanced

648: .seealso: EPSSetRG()
649: @*/
650: PetscErrorCode EPSGetRG(EPS eps,RG *rg)
651: {
652:   PetscFunctionBegin;
654:   PetscAssertPointer(rg,2);
655:   if (!eps->rg) {
656:     PetscCall(RGCreate(PetscObjectComm((PetscObject)eps),&eps->rg));
657:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->rg,(PetscObject)eps,0));
658:     PetscCall(PetscObjectSetOptions((PetscObject)eps->rg,((PetscObject)eps)->options));
659:   }
660:   *rg = eps->rg;
661:   PetscFunctionReturn(PETSC_SUCCESS);
662: }

664: /*@
665:    EPSSetDS - Associates a direct solver object to the eigensolver.

667:    Collective

669:    Input Parameters:
670: +  eps - eigensolver context obtained from EPSCreate()
671: -  ds  - the direct solver object

673:    Note:
674:    Use EPSGetDS() to retrieve the direct solver context (for example,
675:    to free it at the end of the computations).

677:    Level: advanced

679: .seealso: EPSGetDS()
680: @*/
681: PetscErrorCode EPSSetDS(EPS eps,DS ds)
682: {
683:   PetscFunctionBegin;
686:   PetscCheckSameComm(eps,1,ds,2);
687:   PetscCall(PetscObjectReference((PetscObject)ds));
688:   PetscCall(DSDestroy(&eps->ds));
689:   eps->ds = ds;
690:   PetscFunctionReturn(PETSC_SUCCESS);
691: }

693: /*@
694:    EPSGetDS - Obtain the direct solver object associated to the eigensolver object.

696:    Not Collective

698:    Input Parameters:
699: .  eps - eigensolver context obtained from EPSCreate()

701:    Output Parameter:
702: .  ds - direct solver context

704:    Level: advanced

706: .seealso: EPSSetDS()
707: @*/
708: PetscErrorCode EPSGetDS(EPS eps,DS *ds)
709: {
710:   PetscFunctionBegin;
712:   PetscAssertPointer(ds,2);
713:   if (!eps->ds) {
714:     PetscCall(DSCreate(PetscObjectComm((PetscObject)eps),&eps->ds));
715:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)eps->ds,(PetscObject)eps,0));
716:     PetscCall(PetscObjectSetOptions((PetscObject)eps->ds,((PetscObject)eps)->options));
717:   }
718:   *ds = eps->ds;
719:   PetscFunctionReturn(PETSC_SUCCESS);
720: }

722: /*@
723:    EPSIsGeneralized - Ask if the EPS object corresponds to a generalized
724:    eigenvalue problem.

726:    Not Collective

728:    Input Parameter:
729: .  eps - the eigenproblem solver context

731:    Output Parameter:
732: .  is - the answer

734:    Level: intermediate

736: .seealso: EPSIsHermitian(), EPSIsPositive(), EPSIsStructured()
737: @*/
738: PetscErrorCode EPSIsGeneralized(EPS eps,PetscBool* is)
739: {
740:   PetscFunctionBegin;
742:   PetscAssertPointer(is,2);
743:   *is = eps->isgeneralized;
744:   PetscFunctionReturn(PETSC_SUCCESS);
745: }

747: /*@
748:    EPSIsHermitian - Ask if the EPS object corresponds to a Hermitian
749:    eigenvalue problem.

751:    Not Collective

753:    Input Parameter:
754: .  eps - the eigenproblem solver context

756:    Output Parameter:
757: .  is - the answer

759:    Level: intermediate

761: .seealso: EPSIsGeneralized(), EPSIsPositive(), EPSIsStructured()
762: @*/
763: PetscErrorCode EPSIsHermitian(EPS eps,PetscBool* is)
764: {
765:   PetscFunctionBegin;
767:   PetscAssertPointer(is,2);
768:   *is = eps->ishermitian;
769:   PetscFunctionReturn(PETSC_SUCCESS);
770: }

772: /*@
773:    EPSIsPositive - Ask if the EPS object corresponds to an eigenvalue
774:    problem type that requires a positive (semi-) definite matrix B.

776:    Not Collective

778:    Input Parameter:
779: .  eps - the eigenproblem solver context

781:    Output Parameter:
782: .  is - the answer

784:    Level: intermediate

786: .seealso: EPSIsGeneralized(), EPSIsHermitian(), EPSIsStructured()
787: @*/
788: PetscErrorCode EPSIsPositive(EPS eps,PetscBool* is)
789: {
790:   PetscFunctionBegin;
792:   PetscAssertPointer(is,2);
793:   *is = eps->ispositive;
794:   PetscFunctionReturn(PETSC_SUCCESS);
795: }

797: /*@
798:    EPSIsStructured - Ask if the EPS object corresponds to a structured
799:    eigenvalue problem.

801:    Not Collective

803:    Input Parameter:
804: .  eps - the eigenproblem solver context

806:    Output Parameter:
807: .  is - the answer

809:    Note:
810:    The result will be true if the problem type has been set to some
811:    structured type such as EPS_BSE. This is independent of whether the
812:    input matrix has been built with a certain structure with a helper function.

814:    Level: intermediate

816: .seealso: EPSIsGeneralized(), EPSIsHermitian(), EPSIsPositive(), EPSSetProblemType()
817: @*/
818: PetscErrorCode EPSIsStructured(EPS eps,PetscBool* is)
819: {
820:   PetscFunctionBegin;
822:   PetscAssertPointer(is,2);
823:   *is = eps->isstructured;
824:   PetscFunctionReturn(PETSC_SUCCESS);
825: }