Actual source code: epsbasic.c

slepc-3.9.2 2018-07-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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>      /*I "slepceps.h" I*/

 16: PetscFunctionList EPSList = 0;
 17: PetscBool         EPSRegisterAllCalled = PETSC_FALSE;
 18: PetscClassId      EPS_CLASSID = 0;
 19: PetscLogEvent     EPS_SetUp = 0,EPS_Solve = 0;

 21: /*@
 22:    EPSCreate - Creates the default EPS context.

 24:    Collective on MPI_Comm

 26:    Input Parameter:
 27: .  comm - MPI communicator

 29:    Output Parameter:
 30: .  eps - location to put the EPS context

 32:    Note:
 33:    The default EPS type is EPSKRYLOVSCHUR

 35:    Level: beginner

 37: .seealso: EPSSetUp(), EPSSolve(), EPSDestroy(), EPS
 38: @*/
 39: PetscErrorCode EPSCreate(MPI_Comm comm,EPS *outeps)
 40: {
 42:   EPS            eps;

 46:   *outeps = 0;
 47:   EPSInitializePackage();
 48:   SlepcHeaderCreate(eps,EPS_CLASSID,"EPS","Eigenvalue Problem Solver","EPS",comm,EPSDestroy,EPSView);

 50:   eps->max_it          = 0;
 51:   eps->nev             = 1;
 52:   eps->ncv             = 0;
 53:   eps->mpd             = 0;
 54:   eps->nini            = 0;
 55:   eps->nds             = 0;
 56:   eps->target          = 0.0;
 57:   eps->tol             = PETSC_DEFAULT;
 58:   eps->conv            = EPS_CONV_REL;
 59:   eps->stop            = EPS_STOP_BASIC;
 60:   eps->which           = (EPSWhich)0;
 61:   eps->inta            = 0.0;
 62:   eps->intb            = 0.0;
 63:   eps->problem_type    = (EPSProblemType)0;
 64:   eps->extraction      = EPS_RITZ;
 65:   eps->balance         = EPS_BALANCE_NONE;
 66:   eps->balance_its     = 5;
 67:   eps->balance_cutoff  = 1e-8;
 68:   eps->trueres         = PETSC_FALSE;
 69:   eps->trackall        = PETSC_FALSE;
 70:   eps->purify          = PETSC_TRUE;

 72:   eps->converged       = EPSConvergedRelative;
 73:   eps->convergeduser   = NULL;
 74:   eps->convergeddestroy= NULL;
 75:   eps->stopping        = EPSStoppingBasic;
 76:   eps->stoppinguser    = NULL;
 77:   eps->stoppingdestroy = NULL;
 78:   eps->arbitrary       = NULL;
 79:   eps->convergedctx    = NULL;
 80:   eps->stoppingctx     = NULL;
 81:   eps->arbitraryctx    = NULL;
 82:   eps->numbermonitors  = 0;

 84:   eps->st              = NULL;
 85:   eps->ds              = NULL;
 86:   eps->V               = NULL;
 87:   eps->rg              = NULL;
 88:   eps->D               = NULL;
 89:   eps->IS              = NULL;
 90:   eps->defl            = NULL;
 91:   eps->eigr            = NULL;
 92:   eps->eigi            = NULL;
 93:   eps->errest          = NULL;
 94:   eps->rr              = NULL;
 95:   eps->ri              = NULL;
 96:   eps->perm            = NULL;
 97:   eps->nwork           = 0;
 98:   eps->work            = NULL;
 99:   eps->data            = NULL;

101:   eps->state           = EPS_STATE_INITIAL;
102:   eps->categ           = EPS_CATEGORY_KRYLOV;
103:   eps->nconv           = 0;
104:   eps->its             = 0;
105:   eps->nloc            = 0;
106:   eps->nrma            = 0.0;
107:   eps->nrmb            = 0.0;
108:   eps->useds           = PETSC_FALSE;
109:   eps->isgeneralized   = PETSC_FALSE;
110:   eps->ispositive      = PETSC_FALSE;
111:   eps->ishermitian     = PETSC_FALSE;
112:   eps->reason          = EPS_CONVERGED_ITERATING;

114:   PetscNewLog(eps,&eps->sc);
115:   *outeps = eps;
116:   return(0);
117: }

119: /*@C
120:    EPSSetType - Selects the particular solver to be used in the EPS object.

122:    Logically Collective on EPS

124:    Input Parameters:
125: +  eps  - the eigensolver context
126: -  type - a known method

128:    Options Database Key:
129: .  -eps_type <method> - Sets the method; use -help for a list
130:     of available methods

132:    Notes:
133:    See "slepc/include/slepceps.h" for available methods. The default
134:    is EPSKRYLOVSCHUR.

136:    Normally, it is best to use the EPSSetFromOptions() command and
137:    then set the EPS type from the options database rather than by using
138:    this routine.  Using the options database provides the user with
139:    maximum flexibility in evaluating the different available methods.
140:    The EPSSetType() routine is provided for those situations where it
141:    is necessary to set the iterative solver independently of the command
142:    line or options database.

144:    Level: intermediate

146: .seealso: STSetType(), EPSType
147: @*/
148: PetscErrorCode EPSSetType(EPS eps,EPSType type)
149: {
150:   PetscErrorCode ierr,(*r)(EPS);
151:   PetscBool      match;


157:   PetscObjectTypeCompare((PetscObject)eps,type,&match);
158:   if (match) return(0);

160:   PetscFunctionListFind(EPSList,type,&r);
161:   if (!r) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown EPS type given: %s",type);

163:   if (eps->ops->destroy) { (*eps->ops->destroy)(eps); }
164:   PetscMemzero(eps->ops,sizeof(struct _EPSOps));

166:   eps->state = EPS_STATE_INITIAL;
167:   PetscObjectChangeTypeName((PetscObject)eps,type);
168:   (*r)(eps);
169:   return(0);
170: }

172: /*@C
173:    EPSGetType - Gets the EPS type as a string from the EPS object.

175:    Not Collective

177:    Input Parameter:
178: .  eps - the eigensolver context

180:    Output Parameter:
181: .  name - name of EPS method

183:    Level: intermediate

185: .seealso: EPSSetType()
186: @*/
187: PetscErrorCode EPSGetType(EPS eps,EPSType *type)
188: {
192:   *type = ((PetscObject)eps)->type_name;
193:   return(0);
194: }

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

199:    Not Collective

201:    Input Parameters:
202: +  name - name of a new user-defined solver
203: -  function - routine to create the solver context

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

208:    Sample usage:
209: .vb
210:     EPSRegister("my_solver",MySolverCreate);
211: .ve

213:    Then, your solver can be chosen with the procedural interface via
214: $     EPSSetType(eps,"my_solver")
215:    or at runtime via the option
216: $     -eps_type my_solver

218:    Level: advanced

220: .seealso: EPSRegisterAll()
221: @*/
222: PetscErrorCode EPSRegister(const char *name,PetscErrorCode (*function)(EPS))
223: {

227:   PetscFunctionListAdd(&EPSList,name,function);
228:   return(0);
229: }

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

235:    Collective on EPS

237:    Input Parameter:
238: .  eps - eigensolver context obtained from EPSCreate()

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

246:    Level: advanced

248: .seealso: EPSDestroy()
249: @*/
250: PetscErrorCode EPSReset(EPS eps)
251: {

256:   if (!eps) return(0);
257:   if (eps->ops->reset) { (eps->ops->reset)(eps); }
258:   if (eps->st) { STReset(eps->st); }
259:   VecDestroy(&eps->D);
260:   BVDestroy(&eps->V);
261:   VecDestroyVecs(eps->nwork,&eps->work);
262:   eps->nwork = 0;
263:   eps->state = EPS_STATE_INITIAL;
264:   return(0);
265: }

267: /*@
268:    EPSDestroy - Destroys the EPS context.

270:    Collective on EPS

272:    Input Parameter:
273: .  eps - eigensolver context obtained from EPSCreate()

275:    Level: beginner

277: .seealso: EPSCreate(), EPSSetUp(), EPSSolve()
278: @*/
279: PetscErrorCode EPSDestroy(EPS *eps)
280: {

284:   if (!*eps) return(0);
286:   if (--((PetscObject)(*eps))->refct > 0) { *eps = 0; return(0); }
287:   EPSReset(*eps);
288:   if ((*eps)->ops->destroy) { (*(*eps)->ops->destroy)(*eps); }
289:   if ((*eps)->eigr) {
290:     PetscFree4((*eps)->eigr,(*eps)->eigi,(*eps)->errest,(*eps)->perm);
291:   }
292:   if ((*eps)->rr) {
293:     PetscFree2((*eps)->rr,(*eps)->ri);
294:   }
295:   STDestroy(&(*eps)->st);
296:   RGDestroy(&(*eps)->rg);
297:   DSDestroy(&(*eps)->ds);
298:   PetscFree((*eps)->sc);
299:   /* just in case the initial vectors have not been used */
300:   SlepcBasisDestroy_Private(&(*eps)->nds,&(*eps)->defl);
301:   SlepcBasisDestroy_Private(&(*eps)->nini,&(*eps)->IS);
302:   if ((*eps)->convergeddestroy) {
303:     (*(*eps)->convergeddestroy)((*eps)->convergedctx);
304:   }
305:   EPSMonitorCancel(*eps);
306:   PetscHeaderDestroy(eps);
307:   return(0);
308: }

310: /*@
311:    EPSSetTarget - Sets the value of the target.

313:    Logically Collective on EPS

315:    Input Parameters:
316: +  eps    - eigensolver context
317: -  target - the value of the target

319:    Options Database Key:
320: .  -eps_target <scalar> - the value of the target

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

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

330:    Level: intermediate

332: .seealso: EPSGetTarget(), EPSSetWhichEigenpairs()
333: @*/
334: PetscErrorCode EPSSetTarget(EPS eps,PetscScalar target)
335: {

341:   eps->target = target;
342:   if (!eps->st) { EPSGetST(eps,&eps->st); }
343:   STSetDefaultShift(eps->st,target);
344:   return(0);
345: }

347: /*@
348:    EPSGetTarget - Gets the value of the target.

350:    Not Collective

352:    Input Parameter:
353: .  eps - eigensolver context

355:    Output Parameter:
356: .  target - the value of the target

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

361:    Level: intermediate

363: .seealso: EPSSetTarget()
364: @*/
365: PetscErrorCode EPSGetTarget(EPS eps,PetscScalar* target)
366: {
370:   *target = eps->target;
371:   return(0);
372: }

374: /*@
375:    EPSSetInterval - Defines the computational interval for spectrum slicing.

377:    Logically Collective on EPS

379:    Input Parameters:
380: +  eps  - eigensolver context
381: .  inta - left end of the interval
382: -  intb - right end of the interval

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

387:    Notes:
388:    Spectrum slicing is a technique employed for computing all eigenvalues of
389:    symmetric eigenproblems in a given interval. This function provides the
390:    interval to be considered. It must be used in combination with EPS_ALL, see
391:    EPSSetWhichEigenpairs().

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

398:    Level: intermediate

400: .seealso: EPSGetInterval(), EPSSetWhichEigenpairs()
401: @*/
402: PetscErrorCode EPSSetInterval(EPS eps,PetscReal inta,PetscReal intb)
403: {
408:   if (inta>intb) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be inta<intb");
409:   if (eps->inta != inta || eps->intb != intb) {
410:     eps->inta = inta;
411:     eps->intb = intb;
412:     eps->state = EPS_STATE_INITIAL;
413:   }
414:   return(0);
415: }

417: /*@
418:    EPSGetInterval - Gets the computational interval for spectrum slicing.

420:    Not Collective

422:    Input Parameter:
423: .  eps - eigensolver context

425:    Output Parameters:
426: +  inta - left end of the interval
427: -  intb - right end of the interval

429:    Level: intermediate

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

434: .seealso: EPSSetInterval()
435: @*/
436: PetscErrorCode EPSGetInterval(EPS eps,PetscReal* inta,PetscReal* intb)
437: {
440:   if (inta) *inta = eps->inta;
441:   if (intb) *intb = eps->intb;
442:   return(0);
443: }

445: /*@
446:    EPSSetST - Associates a spectral transformation object to the eigensolver.

448:    Collective on EPS

450:    Input Parameters:
451: +  eps - eigensolver context obtained from EPSCreate()
452: -  st   - the spectral transformation object

454:    Note:
455:    Use EPSGetST() to retrieve the spectral transformation context (for example,
456:    to free it at the end of the computations).

458:    Level: advanced

460: .seealso: EPSGetST()
461: @*/
462: PetscErrorCode EPSSetST(EPS eps,ST st)
463: {

470:   PetscObjectReference((PetscObject)st);
471:   STDestroy(&eps->st);
472:   eps->st = st;
473:   PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->st);
474:   return(0);
475: }

477: /*@
478:    EPSGetST - Obtain the spectral transformation (ST) object associated
479:    to the eigensolver object.

481:    Not Collective

483:    Input Parameters:
484: .  eps - eigensolver context obtained from EPSCreate()

486:    Output Parameter:
487: .  st - spectral transformation context

489:    Level: intermediate

491: .seealso: EPSSetST()
492: @*/
493: PetscErrorCode EPSGetST(EPS eps,ST *st)
494: {

500:   if (!eps->st) {
501:     STCreate(PetscObjectComm((PetscObject)eps),&eps->st);
502:     PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->st);
503:   }
504:   *st = eps->st;
505:   return(0);
506: }

508: /*@
509:    EPSSetBV - Associates a basis vectors object to the eigensolver.

511:    Collective on EPS

513:    Input Parameters:
514: +  eps - eigensolver context obtained from EPSCreate()
515: -  V   - the basis vectors object

517:    Note:
518:    Use EPSGetBV() to retrieve the basis vectors context (for example,
519:    to free them at the end of the computations).

521:    Level: advanced

523: .seealso: EPSGetBV()
524: @*/
525: PetscErrorCode EPSSetBV(EPS eps,BV V)
526: {

533:   PetscObjectReference((PetscObject)V);
534:   BVDestroy(&eps->V);
535:   eps->V = V;
536:   PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->V);
537:   return(0);
538: }

540: /*@
541:    EPSGetBV - Obtain the basis vectors object associated to the eigensolver object.

543:    Not Collective

545:    Input Parameters:
546: .  eps - eigensolver context obtained from EPSCreate()

548:    Output Parameter:
549: .  V - basis vectors context

551:    Level: advanced

553: .seealso: EPSSetBV()
554: @*/
555: PetscErrorCode EPSGetBV(EPS eps,BV *V)
556: {

562:   if (!eps->V) {
563:     BVCreate(PetscObjectComm((PetscObject)eps),&eps->V);
564:     PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->V);
565:   }
566:   *V = eps->V;
567:   return(0);
568: }

570: /*@
571:    EPSSetRG - Associates a region object to the eigensolver.

573:    Collective on EPS

575:    Input Parameters:
576: +  eps - eigensolver context obtained from EPSCreate()
577: -  rg  - the region object

579:    Note:
580:    Use EPSGetRG() to retrieve the region context (for example,
581:    to free it at the end of the computations).

583:    Level: advanced

585: .seealso: EPSGetRG()
586: @*/
587: PetscErrorCode EPSSetRG(EPS eps,RG rg)
588: {

595:   PetscObjectReference((PetscObject)rg);
596:   RGDestroy(&eps->rg);
597:   eps->rg = rg;
598:   PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->rg);
599:   return(0);
600: }

602: /*@
603:    EPSGetRG - Obtain the region object associated to the eigensolver.

605:    Not Collective

607:    Input Parameters:
608: .  eps - eigensolver context obtained from EPSCreate()

610:    Output Parameter:
611: .  rg - region context

613:    Level: advanced

615: .seealso: EPSSetRG()
616: @*/
617: PetscErrorCode EPSGetRG(EPS eps,RG *rg)
618: {

624:   if (!eps->rg) {
625:     RGCreate(PetscObjectComm((PetscObject)eps),&eps->rg);
626:     PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->rg);
627:   }
628:   *rg = eps->rg;
629:   return(0);
630: }

632: /*@
633:    EPSSetDS - Associates a direct solver object to the eigensolver.

635:    Collective on EPS

637:    Input Parameters:
638: +  eps - eigensolver context obtained from EPSCreate()
639: -  ds  - the direct solver object

641:    Note:
642:    Use EPSGetDS() to retrieve the direct solver context (for example,
643:    to free it at the end of the computations).

645:    Level: advanced

647: .seealso: EPSGetDS()
648: @*/
649: PetscErrorCode EPSSetDS(EPS eps,DS ds)
650: {

657:   PetscObjectReference((PetscObject)ds);
658:   DSDestroy(&eps->ds);
659:   eps->ds = ds;
660:   PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->ds);
661:   return(0);
662: }

664: /*@
665:    EPSGetDS - Obtain the direct solver object associated to the eigensolver object.

667:    Not Collective

669:    Input Parameters:
670: .  eps - eigensolver context obtained from EPSCreate()

672:    Output Parameter:
673: .  ds - direct solver context

675:    Level: advanced

677: .seealso: EPSSetDS()
678: @*/
679: PetscErrorCode EPSGetDS(EPS eps,DS *ds)
680: {

686:   if (!eps->ds) {
687:     DSCreate(PetscObjectComm((PetscObject)eps),&eps->ds);
688:     PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->ds);
689:   }
690:   *ds = eps->ds;
691:   return(0);
692: }

694: /*@
695:    EPSIsGeneralized - Ask if the EPS object corresponds to a generalized
696:    eigenvalue problem.

698:    Not collective

700:    Input Parameter:
701: .  eps - the eigenproblem solver context

703:    Output Parameter:
704: .  is - the answer

706:    Level: intermediate

708: .seealso: EPSIsHermitian(), EPSIsPositive()
709: @*/
710: PetscErrorCode EPSIsGeneralized(EPS eps,PetscBool* is)
711: {
715:   *is = eps->isgeneralized;
716:   return(0);
717: }

719: /*@
720:    EPSIsHermitian - Ask if the EPS object corresponds to a Hermitian
721:    eigenvalue problem.

723:    Not collective

725:    Input Parameter:
726: .  eps - the eigenproblem solver context

728:    Output Parameter:
729: .  is - the answer

731:    Level: intermediate

733: .seealso: EPSIsGeneralized(), EPSIsPositive()
734: @*/
735: PetscErrorCode EPSIsHermitian(EPS eps,PetscBool* is)
736: {
740:   *is = eps->ishermitian;
741:   return(0);
742: }

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

748:    Not collective

750:    Input Parameter:
751: .  eps - the eigenproblem solver context

753:    Output Parameter:
754: .  is - the answer

756:    Level: intermediate

758: .seealso: EPSIsGeneralized(), EPSIsHermitian()
759: @*/
760: PetscErrorCode EPSIsPositive(EPS eps,PetscBool* is)
761: {
765:   *is = eps->ispositive;
766:   return(0);
767: }