Actual source code: pepbasic.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 PEP routines
 12: */

 14: #include <slepc/private/pepimpl.h>

 16: /* Logging support */
 17: PetscClassId      PEP_CLASSID = 0;
 18: PetscLogEvent     PEP_SetUp = 0,PEP_Solve = 0,PEP_Refine = 0,PEP_CISS_SVD = 0;

 20: /* List of registered PEP routines */
 21: PetscFunctionList PEPList = NULL;
 22: PetscBool         PEPRegisterAllCalled = PETSC_FALSE;

 24: /* List of registered PEP monitors */
 25: PetscFunctionList PEPMonitorList              = NULL;
 26: PetscFunctionList PEPMonitorCreateList        = NULL;
 27: PetscFunctionList PEPMonitorDestroyList       = NULL;
 28: PetscBool         PEPMonitorRegisterAllCalled = PETSC_FALSE;

 30: /*@
 31:    PEPCreate - Creates the `PEP` context.

 33:    Collective

 35:    Input Parameter:
 36: .  comm - MPI communicator

 38:    Output Parameter:
 39: .  outpep - location to put the `PEP` context

 41:    Note:
 42:    The default `PEP` type is `PEPTOAR`.

 44:    Level: beginner

 46: .seealso: [](ch:pep), `PEPSetUp()`, `PEPSolve()`, `PEPDestroy()`, `PEP`
 47: @*/
 48: PetscErrorCode PEPCreate(MPI_Comm comm,PEP *outpep)
 49: {
 50:   PEP            pep;

 52:   PetscFunctionBegin;
 53:   PetscAssertPointer(outpep,2);
 54:   PetscCall(PEPInitializePackage());
 55:   PetscCall(SlepcHeaderCreate(pep,PEP_CLASSID,"PEP","Polynomial Eigenvalue Problem","PEP",comm,PEPDestroy,PEPView));

 57:   pep->max_it          = PETSC_DETERMINE;
 58:   pep->nev             = 1;
 59:   pep->ncv             = PETSC_DETERMINE;
 60:   pep->mpd             = PETSC_DETERMINE;
 61:   pep->nini            = 0;
 62:   pep->target          = 0.0;
 63:   pep->tol             = PETSC_DETERMINE;
 64:   pep->conv            = PEP_CONV_REL;
 65:   pep->stop            = PEP_STOP_BASIC;
 66:   pep->which           = (PEPWhich)0;
 67:   pep->basis           = PEP_BASIS_MONOMIAL;
 68:   pep->problem_type    = (PEPProblemType)0;
 69:   pep->scale           = PEP_SCALE_NONE;
 70:   pep->sfactor         = 1.0;
 71:   pep->dsfactor        = 1.0;
 72:   pep->sits            = 5;
 73:   pep->slambda         = 1.0;
 74:   pep->refine          = PEP_REFINE_NONE;
 75:   pep->npart           = 1;
 76:   pep->rtol            = PETSC_DETERMINE;
 77:   pep->rits            = PETSC_DETERMINE;
 78:   pep->scheme          = (PEPRefineScheme)0;
 79:   pep->extract         = (PEPExtract)0;
 80:   pep->trackall        = PETSC_FALSE;

 82:   pep->converged       = PEPConvergedRelative;
 83:   pep->convergeduser   = NULL;
 84:   pep->convergeddestroy= NULL;
 85:   pep->stopping        = PEPStoppingBasic;
 86:   pep->stoppinguser    = NULL;
 87:   pep->stoppingdestroy = NULL;
 88:   pep->convergedctx    = NULL;
 89:   pep->stoppingctx     = NULL;
 90:   pep->numbermonitors  = 0;

 92:   pep->st              = NULL;
 93:   pep->ds              = NULL;
 94:   pep->V               = NULL;
 95:   pep->rg              = NULL;
 96:   pep->A               = NULL;
 97:   pep->nmat            = 0;
 98:   pep->Dl              = NULL;
 99:   pep->Dr              = NULL;
100:   pep->IS              = NULL;
101:   pep->eigr            = NULL;
102:   pep->eigi            = NULL;
103:   pep->errest          = NULL;
104:   pep->perm            = NULL;
105:   pep->pbc             = NULL;
106:   pep->solvematcoeffs  = NULL;
107:   pep->nwork           = 0;
108:   pep->work            = NULL;
109:   pep->refineksp       = NULL;
110:   pep->refinesubc      = NULL;
111:   pep->data            = NULL;

113:   pep->state           = PEP_STATE_INITIAL;
114:   pep->nconv           = 0;
115:   pep->its             = 0;
116:   pep->n               = 0;
117:   pep->nloc            = 0;
118:   pep->nrma            = NULL;
119:   pep->sfactor_set     = PETSC_FALSE;
120:   pep->lineariz        = PETSC_FALSE;
121:   pep->reason          = PEP_CONVERGED_ITERATING;

123:   PetscCall(PetscNew(&pep->sc));
124:   *outpep = pep;
125:   PetscFunctionReturn(PETSC_SUCCESS);
126: }

128: /*@
129:    PEPSetType - Selects the particular solver to be used in the `PEP` object.

131:    Logically Collective

133:    Input Parameters:
134: +  pep  - the polynomial eigensolver context
135: -  type - a known method

137:    Options Database Key:
138: .  -pep_type \<type\> - sets the method; use `-help` for a list of available methods

140:    Notes:
141:    See `PEPType` for available methods. The default is `PEPTOAR`.

143:    Normally, it is best to use the `PEPSetFromOptions()` command and
144:    then set the `PEP` type from the options database rather than by using
145:    this routine.  Using the options database provides the user with
146:    maximum flexibility in evaluating the different available methods.
147:    The `PEPSetType()` routine is provided for those situations where it
148:    is necessary to set the iterative solver independently of the command
149:    line or options database.

151:    Level: intermediate

153: .seealso: [](ch:pep), `PEPType`
154: @*/
155: PetscErrorCode PEPSetType(PEP pep,PEPType type)
156: {
157:   PetscErrorCode (*r)(PEP);
158:   PetscBool      match;

160:   PetscFunctionBegin;
162:   PetscAssertPointer(type,2);

164:   PetscCall(PetscObjectTypeCompare((PetscObject)pep,type,&match));
165:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

167:   PetscCall(PetscFunctionListFind(PEPList,type,&r));
168:   PetscCheck(r,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown PEP type given: %s",type);

170:   PetscTryTypeMethod(pep,destroy);
171:   PetscCall(PetscMemzero(pep->ops,sizeof(struct _PEPOps)));

173:   pep->state = PEP_STATE_INITIAL;
174:   PetscCall(PetscObjectChangeTypeName((PetscObject)pep,type));
175:   PetscCall((*r)(pep));
176:   PetscFunctionReturn(PETSC_SUCCESS);
177: }

179: /*@
180:    PEPGetType - Gets the `PEP` type as a string from the `PEP` object.

182:    Not Collective

184:    Input Parameter:
185: .  pep - the polynomial eigensolver context

187:    Output Parameter:
188: .  type - name of `PEP` method

190:    Level: intermediate

192: .seealso: [](ch:pep), `PEPSetType()`
193: @*/
194: PetscErrorCode PEPGetType(PEP pep,PEPType *type)
195: {
196:   PetscFunctionBegin;
198:   PetscAssertPointer(type,2);
199:   *type = ((PetscObject)pep)->type_name;
200:   PetscFunctionReturn(PETSC_SUCCESS);
201: }

203: /*@C
204:    PEPRegister - Adds a method to the polynomial eigenproblem solver package.

206:    Not Collective

208:    Input Parameters:
209: +  name     - name of a new user-defined solver
210: -  function - routine to create the solver context

212:    Note:
213:    `PEPRegister()` may be called multiple times to add several user-defined solvers.

215:    Example Usage:
216: .vb
217:    PEPRegister("my_solver",MySolverCreate);
218: .ve

220:    Then, your solver can be chosen with the procedural interface via
221: .vb
222:    PEPSetType(pep,"my_solver")
223: .ve
224:    or at runtime via the option `-pep_type my_solver`.

226:    Level: advanced

228: .seealso: [](ch:pep), `PEPRegisterAll()`
229: @*/
230: PetscErrorCode PEPRegister(const char *name,PetscErrorCode (*function)(PEP))
231: {
232:   PetscFunctionBegin;
233:   PetscCall(PEPInitializePackage());
234:   PetscCall(PetscFunctionListAdd(&PEPList,name,function));
235:   PetscFunctionReturn(PETSC_SUCCESS);
236: }

238: /*@C
239:    PEPMonitorRegister - Registers a `PEP` monitor routine that may be accessed with
240:    `PEPMonitorSetFromOptions()`.

242:    Not Collective

244:    Input Parameters:
245: +  name    - name of a new monitor routine
246: .  vtype   - a `PetscViewerType` for the output
247: .  format  - a `PetscViewerFormat` for the output
248: .  monitor - monitor routine, see `PEPMonitorRegisterFn`
249: .  create  - creation routine, or `NULL`
250: -  destroy - destruction routine, or `NULL`

252:    Notes:
253:    `PEPMonitorRegister()` may be called multiple times to add several user-defined monitors.

255:    The calling sequence for the given function matches the calling sequence of `PEPMonitorFn`
256:    functions passed to `PEPMonitorSet()` with the additional requirement that its final argument
257:    be a `PetscViewerAndFormat`.

259:    Example Usage:
260: .vb
261:    PEPMonitorRegister("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: .vb
266:    PEPMonitorSetFromOptions(pep,"-pep_monitor_my_monitor","my_monitor",NULL);
267: .ve
268:    or at runtime via the option `-pep_monitor_my_monitor`.

270:    Level: advanced

272: .seealso: [](ch:pep), `PEPMonitorSet()`, `PEPMonitorRegisterAll()`, `PEPMonitorSetFromOptions()`
273: @*/
274: PetscErrorCode PEPMonitorRegister(const char name[],PetscViewerType vtype,PetscViewerFormat format,PEPMonitorRegisterFn *monitor,PEPMonitorRegisterCreateFn *create,PEPMonitorRegisterDestroyFn *destroy)
275: {
276:   char           key[PETSC_MAX_PATH_LEN];

278:   PetscFunctionBegin;
279:   PetscCall(PEPInitializePackage());
280:   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
281:   PetscCall(PetscFunctionListAdd(&PEPMonitorList,key,monitor));
282:   if (create)  PetscCall(PetscFunctionListAdd(&PEPMonitorCreateList,key,create));
283:   if (destroy) PetscCall(PetscFunctionListAdd(&PEPMonitorDestroyList,key,destroy));
284:   PetscFunctionReturn(PETSC_SUCCESS);
285: }

287: /*@
288:    PEPReset - Resets the `PEP` context to the initial state (prior to setup)
289:    and destroys any allocated `Vec`s and `Mat`s.

291:    Collective

293:    Input Parameter:
294: .  pep - the polynomial eigensolver context

296:    Level: advanced

298: .seealso: [](ch:pep), `PEPDestroy()`
299: @*/
300: PetscErrorCode PEPReset(PEP pep)
301: {
302:   PetscFunctionBegin;
304:   if (!pep) PetscFunctionReturn(PETSC_SUCCESS);
305:   PetscTryTypeMethod(pep,reset);
306:   if (pep->st) PetscCall(STReset(pep->st));
307:   if (pep->refineksp) PetscCall(KSPReset(pep->refineksp));
308:   if (pep->nmat) {
309:     PetscCall(MatDestroyMatrices(pep->nmat,&pep->A));
310:     PetscCall(PetscFree2(pep->pbc,pep->nrma));
311:     PetscCall(PetscFree(pep->solvematcoeffs));
312:     pep->nmat = 0;
313:   }
314:   PetscCall(VecDestroy(&pep->Dl));
315:   PetscCall(VecDestroy(&pep->Dr));
316:   PetscCall(BVDestroy(&pep->V));
317:   PetscCall(VecDestroyVecs(pep->nwork,&pep->work));
318:   pep->nwork = 0;
319:   pep->state = PEP_STATE_INITIAL;
320:   PetscFunctionReturn(PETSC_SUCCESS);
321: }

323: /*@
324:    PEPDestroy - Destroys the `PEP` context.

326:    Collective

328:    Input Parameter:
329: .  pep - the polynomial eigensolver context

331:    Level: beginner

333: .seealso: [](ch:pep), `PEPCreate()`, `PEPSetUp()`, `PEPSolve()`
334: @*/
335: PetscErrorCode PEPDestroy(PEP *pep)
336: {
337:   PetscFunctionBegin;
338:   if (!*pep) PetscFunctionReturn(PETSC_SUCCESS);
340:   if (--((PetscObject)*pep)->refct > 0) { *pep = NULL; PetscFunctionReturn(PETSC_SUCCESS); }
341:   PetscCall(PEPReset(*pep));
342:   PetscTryTypeMethod(*pep,destroy);
343:   if ((*pep)->eigr) PetscCall(PetscFree4((*pep)->eigr,(*pep)->eigi,(*pep)->errest,(*pep)->perm));
344:   PetscCall(STDestroy(&(*pep)->st));
345:   PetscCall(RGDestroy(&(*pep)->rg));
346:   PetscCall(DSDestroy(&(*pep)->ds));
347:   PetscCall(KSPDestroy(&(*pep)->refineksp));
348:   PetscCall(PetscSubcommDestroy(&(*pep)->refinesubc));
349:   PetscCall(PetscFree((*pep)->sc));
350:   /* just in case the initial vectors have not been used */
351:   PetscCall(SlepcBasisDestroy_Private(&(*pep)->nini,&(*pep)->IS));
352:   if ((*pep)->convergeddestroy) PetscCall((*(*pep)->convergeddestroy)(&(*pep)->convergedctx));
353:   if ((*pep)->stoppingdestroy) PetscCall((*(*pep)->stoppingdestroy)(&(*pep)->stoppingctx));
354:   PetscCall(PEPMonitorCancel(*pep));
355:   PetscCall(PetscHeaderDestroy(pep));
356:   PetscFunctionReturn(PETSC_SUCCESS);
357: }

359: /*@
360:    PEPSetBV - Associates a basis vectors object to the polynomial eigensolver.

362:    Collective

364:    Input Parameters:
365: +  pep - the polynomial eigensolver context
366: -  bv  - the basis vectors object

368:    Note:
369:    Use `PEPGetBV()` to retrieve the basis vectors context (for example,
370:    to free it at the end of the computations).

372:    Level: advanced

374: .seealso: [](ch:pep), `PEPGetBV()`
375: @*/
376: PetscErrorCode PEPSetBV(PEP pep,BV bv)
377: {
378:   PetscFunctionBegin;
381:   PetscCheckSameComm(pep,1,bv,2);
382:   PetscCall(PetscObjectReference((PetscObject)bv));
383:   PetscCall(BVDestroy(&pep->V));
384:   pep->V = bv;
385:   PetscFunctionReturn(PETSC_SUCCESS);
386: }

388: /*@
389:    PEPGetBV - Obtain the basis vectors object associated to the polynomial
390:    eigensolver object.

392:    Not Collective

394:    Input Parameter:
395: .  pep - the polynomial eigensolver context

397:    Output Parameter:
398: .  bv - basis vectors context

400:    Level: advanced

402: .seealso: [](ch:pep), `PEPSetBV()`
403: @*/
404: PetscErrorCode PEPGetBV(PEP pep,BV *bv)
405: {
406:   PetscFunctionBegin;
408:   PetscAssertPointer(bv,2);
409:   if (!pep->V) {
410:     PetscCall(BVCreate(PetscObjectComm((PetscObject)pep),&pep->V));
411:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pep->V,(PetscObject)pep,0));
412:     PetscCall(PetscObjectSetOptions((PetscObject)pep->V,((PetscObject)pep)->options));
413:   }
414:   *bv = pep->V;
415:   PetscFunctionReturn(PETSC_SUCCESS);
416: }

418: /*@
419:    PEPSetRG - Associates a region object to the polynomial eigensolver.

421:    Collective

423:    Input Parameters:
424: +  pep - the polynomial eigensolver context
425: -  rg  - the region object

427:    Note:
428:    Use `PEPGetRG()` to retrieve the region context (for example,
429:    to free it at the end of the computations).

431:    Level: advanced

433: .seealso: [](ch:pep), `PEPGetRG()`
434: @*/
435: PetscErrorCode PEPSetRG(PEP pep,RG rg)
436: {
437:   PetscFunctionBegin;
439:   if (rg) {
441:     PetscCheckSameComm(pep,1,rg,2);
442:   }
443:   PetscCall(PetscObjectReference((PetscObject)rg));
444:   PetscCall(RGDestroy(&pep->rg));
445:   pep->rg = rg;
446:   PetscFunctionReturn(PETSC_SUCCESS);
447: }

449: /*@
450:    PEPGetRG - Obtain the region object associated to the
451:    polynomial eigensolver object.

453:    Not Collective

455:    Input Parameter:
456: .  pep - the polynomial eigensolver context

458:    Output Parameter:
459: .  rg - region context

461:    Level: advanced

463: .seealso: [](ch:pep), `PEPSetRG()`
464: @*/
465: PetscErrorCode PEPGetRG(PEP pep,RG *rg)
466: {
467:   PetscFunctionBegin;
469:   PetscAssertPointer(rg,2);
470:   if (!pep->rg) {
471:     PetscCall(RGCreate(PetscObjectComm((PetscObject)pep),&pep->rg));
472:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pep->rg,(PetscObject)pep,0));
473:     PetscCall(PetscObjectSetOptions((PetscObject)pep->rg,((PetscObject)pep)->options));
474:   }
475:   *rg = pep->rg;
476:   PetscFunctionReturn(PETSC_SUCCESS);
477: }

479: /*@
480:    PEPSetDS - Associates a direct solver object to the polynomial eigensolver.

482:    Collective

484:    Input Parameters:
485: +  pep - the polynomial eigensolver context
486: -  ds  - the direct solver object

488:    Note:
489:    Use `PEPGetDS()` to retrieve the direct solver context (for example,
490:    to free it at the end of the computations).

492:    Level: advanced

494: .seealso: [](ch:pep), `PEPGetDS()`
495: @*/
496: PetscErrorCode PEPSetDS(PEP pep,DS ds)
497: {
498:   PetscFunctionBegin;
501:   PetscCheckSameComm(pep,1,ds,2);
502:   PetscCall(PetscObjectReference((PetscObject)ds));
503:   PetscCall(DSDestroy(&pep->ds));
504:   pep->ds = ds;
505:   PetscFunctionReturn(PETSC_SUCCESS);
506: }

508: /*@
509:    PEPGetDS - Obtain the direct solver object associated to the
510:    polynomial eigensolver object.

512:    Not Collective

514:    Input Parameter:
515: .  pep - the polynomial eigensolver context

517:    Output Parameter:
518: .  ds - direct solver context

520:    Level: advanced

522: .seealso: [](ch:pep), `PEPSetDS()`
523: @*/
524: PetscErrorCode PEPGetDS(PEP pep,DS *ds)
525: {
526:   PetscFunctionBegin;
528:   PetscAssertPointer(ds,2);
529:   if (!pep->ds) {
530:     PetscCall(DSCreate(PetscObjectComm((PetscObject)pep),&pep->ds));
531:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pep->ds,(PetscObject)pep,0));
532:     PetscCall(PetscObjectSetOptions((PetscObject)pep->ds,((PetscObject)pep)->options));
533:   }
534:   *ds = pep->ds;
535:   PetscFunctionReturn(PETSC_SUCCESS);
536: }

538: /*@
539:    PEPSetST - Associates a spectral transformation object to the eigensolver.

541:    Collective

543:    Input Parameters:
544: +  pep - the polynomial eigensolver context
545: -  st  - the spectral transformation object

547:    Note:
548:    Use `PEPGetST()` to retrieve the spectral transformation context (for example,
549:    to free it at the end of the computations).

551:    Level: advanced

553: .seealso: [](ch:pep), `PEPGetST()`
554: @*/
555: PetscErrorCode PEPSetST(PEP pep,ST st)
556: {
557:   PetscFunctionBegin;
560:   PetscCheckSameComm(pep,1,st,2);
561:   PetscCall(PetscObjectReference((PetscObject)st));
562:   PetscCall(STDestroy(&pep->st));
563:   pep->st = st;
564:   PetscFunctionReturn(PETSC_SUCCESS);
565: }

567: /*@
568:    PEPGetST - Obtain the spectral transformation (`ST`) object associated
569:    to the eigensolver object.

571:    Not Collective

573:    Input Parameter:
574: .  pep - the polynomial eigensolver context

576:    Output Parameter:
577: .  st - spectral transformation context

579:    Level: intermediate

581: .seealso: [](ch:pep), `PEPSetST()`
582: @*/
583: PetscErrorCode PEPGetST(PEP pep,ST *st)
584: {
585:   PetscFunctionBegin;
587:   PetscAssertPointer(st,2);
588:   if (!pep->st) {
589:     PetscCall(STCreate(PetscObjectComm((PetscObject)pep),&pep->st));
590:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pep->st,(PetscObject)pep,0));
591:     PetscCall(PetscObjectSetOptions((PetscObject)pep->st,((PetscObject)pep)->options));
592:   }
593:   *st = pep->st;
594:   PetscFunctionReturn(PETSC_SUCCESS);
595: }

597: /*@
598:    PEPRefineGetKSP - Obtain the `KSP` object used by the eigensolver
599:    object in the refinement phase.

601:    Collective

603:    Input Parameter:
604: .  pep - the polynomial eigensolver context

606:    Output Parameter:
607: .  ksp - the linear solver context

609:    Level: advanced

611: .seealso: [](ch:pep), `PEPSetRefine()`
612: @*/
613: PetscErrorCode PEPRefineGetKSP(PEP pep,KSP *ksp)
614: {
615:   MPI_Comm       comm;

617:   PetscFunctionBegin;
619:   PetscAssertPointer(ksp,2);
620:   if (!pep->refineksp) {
621:     if (pep->npart>1) {
622:       /* Split in subcomunicators */
623:       PetscCall(PetscSubcommCreate(PetscObjectComm((PetscObject)pep),&pep->refinesubc));
624:       PetscCall(PetscSubcommSetNumber(pep->refinesubc,pep->npart));
625:       PetscCall(PetscSubcommSetType(pep->refinesubc,PETSC_SUBCOMM_CONTIGUOUS));
626:       PetscCall(PetscSubcommGetChild(pep->refinesubc,&comm));
627:     } else PetscCall(PetscObjectGetComm((PetscObject)pep,&comm));
628:     PetscCall(KSPCreate(comm,&pep->refineksp));
629:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pep->refineksp,(PetscObject)pep,0));
630:     PetscCall(PetscObjectSetOptions((PetscObject)pep->refineksp,((PetscObject)pep)->options));
631:     PetscCall(KSPSetOptionsPrefix(*ksp,((PetscObject)pep)->prefix));
632:     PetscCall(KSPAppendOptionsPrefix(*ksp,"pep_refine_"));
633:     PetscCall(KSPSetTolerances(pep->refineksp,SlepcDefaultTol(pep->rtol),PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
634:   }
635:   *ksp = pep->refineksp;
636:   PetscFunctionReturn(PETSC_SUCCESS);
637: }

639: /*@
640:    PEPSetTarget - Sets the value of the target.

642:    Logically Collective

644:    Input Parameters:
645: +  pep    - the polynomial eigensolver context
646: -  target - the value of the target

648:    Options Database Key:
649: .  -pep_target \<target\> - the value of the target

651:    Notes:
652:    The target is a scalar value used to determine the portion of the spectrum
653:    of interest. It is used in combination with `PEPSetWhichEigenpairs()`.

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

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

662:    Level: intermediate

664: .seealso: [](ch:pep), `PEPGetTarget()`, `PEPSetWhichEigenpairs()`
665: @*/
666: PetscErrorCode PEPSetTarget(PEP pep,PetscScalar target)
667: {
668:   PetscFunctionBegin;
671:   pep->target = target;
672:   if (!pep->st) PetscCall(PEPGetST(pep,&pep->st));
673:   PetscCall(STSetDefaultShift(pep->st,target));
674:   PetscFunctionReturn(PETSC_SUCCESS);
675: }

677: /*@
678:    PEPGetTarget - Gets the value of the target.

680:    Not Collective

682:    Input Parameter:
683: .  pep - the polynomial eigensolver context

685:    Output Parameter:
686: .  target - the value of the target

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

691:    Level: intermediate

693: .seealso: [](ch:pep), `PEPSetTarget()`
694: @*/
695: PetscErrorCode PEPGetTarget(PEP pep,PetscScalar* target)
696: {
697:   PetscFunctionBegin;
699:   PetscAssertPointer(target,2);
700:   *target = pep->target;
701:   PetscFunctionReturn(PETSC_SUCCESS);
702: }

704: /*@
705:    PEPSetInterval - Defines the computational interval for spectrum slicing.

707:    Logically Collective

709:    Input Parameters:
710: +  pep  - the polynomial eigensolver context
711: .  inta - left end of the interval
712: -  intb - right end of the interval

714:    Options Database Key:
715: .  -pep_interval <a,b> - set $[a,b]$ as the interval of interest

717:    Notes:
718:    Spectrum slicing is a technique employed for computing all eigenvalues of
719:    symmetric eigenproblems in a given interval, see section [](#sec:qslice).
720:    This function provides the interval to be considered. It must be used in
721:    combination with `PEP_ALL`, see `PEPSetWhichEigenpairs()`. Note that in
722:    polynomial eigenproblems spectrum slicing is implemented in `PEPSTOAR` only.

724:    In the command-line option, two values must be provided. For an open interval,
725:    one can give an infinite, e.g., `-pep_interval 1.0,inf` or `-pep_interval -inf,1.0`.
726:    An open interval in the programmatic interface can be specified with
727:    `PETSC_MAX_REAL` and -`PETSC_MAX_REAL`.

729:    Level: intermediate

731: .seealso: [](ch:pep), [](#sec:qslice), `PEPGetInterval()`, `PEPSetWhichEigenpairs()`, `PEPSTOAR`
732: @*/
733: PetscErrorCode PEPSetInterval(PEP pep,PetscReal inta,PetscReal intb)
734: {
735:   PetscFunctionBegin;
739:   PetscCheck(inta<intb,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be inta<intb");
740:   if (pep->inta != inta || pep->intb != intb) {
741:     pep->inta = inta;
742:     pep->intb = intb;
743:     pep->state = PEP_STATE_INITIAL;
744:   }
745:   PetscFunctionReturn(PETSC_SUCCESS);
746: }

748: /*@
749:    PEPGetInterval - Gets the computational interval for spectrum slicing.

751:    Not Collective

753:    Input Parameter:
754: .  pep - the polynomial eigensolver context

756:    Output Parameters:
757: +  inta - left end of the interval
758: -  intb - right end of the interval

760:    Level: intermediate

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

765: .seealso: [](ch:pep), `PEPSetInterval()`
766: @*/
767: PetscErrorCode PEPGetInterval(PEP pep,PetscReal* inta,PetscReal* intb)
768: {
769:   PetscFunctionBegin;
771:   if (inta) *inta = pep->inta;
772:   if (intb) *intb = pep->intb;
773:   PetscFunctionReturn(PETSC_SUCCESS);
774: }