Actual source code: epsbasic.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*
  2:    The basic EPS routines, Create, Destroy, etc. are here.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc/private/epsimpl.h>      /*I "slepceps.h" I*/

 26: PetscFunctionList EPSList = 0;
 27: PetscBool         EPSRegisterAllCalled = PETSC_FALSE;
 28: PetscClassId      EPS_CLASSID = 0;
 29: PetscLogEvent     EPS_SetUp = 0,EPS_Solve = 0;

 33: /*@
 34:    EPSCreate - Creates the default EPS context.

 36:    Collective on MPI_Comm

 38:    Input Parameter:
 39: .  comm - MPI communicator

 41:    Output Parameter:
 42: .  eps - location to put the EPS context

 44:    Note:
 45:    The default EPS type is EPSKRYLOVSCHUR

 47:    Level: beginner

 49: .seealso: EPSSetUp(), EPSSolve(), EPSDestroy(), EPS
 50: @*/
 51: PetscErrorCode EPSCreate(MPI_Comm comm,EPS *outeps)
 52: {
 54:   EPS            eps;

 58:   *outeps = 0;
 59:   EPSInitializePackage();
 60:   SlepcHeaderCreate(eps,EPS_CLASSID,"EPS","Eigenvalue Problem Solver","EPS",comm,EPSDestroy,EPSView);

 62:   eps->max_it          = 0;
 63:   eps->nev             = 1;
 64:   eps->ncv             = 0;
 65:   eps->mpd             = 0;
 66:   eps->nini            = 0;
 67:   eps->nds             = 0;
 68:   eps->target          = 0.0;
 69:   eps->tol             = PETSC_DEFAULT;
 70:   eps->conv            = EPS_CONV_REL;
 71:   eps->stop            = EPS_STOP_BASIC;
 72:   eps->which           = (EPSWhich)0;
 73:   eps->inta            = 0.0;
 74:   eps->intb            = 0.0;
 75:   eps->problem_type    = (EPSProblemType)0;
 76:   eps->extraction      = EPS_RITZ;
 77:   eps->balance         = EPS_BALANCE_NONE;
 78:   eps->balance_its     = 5;
 79:   eps->balance_cutoff  = 1e-8;
 80:   eps->trueres         = PETSC_FALSE;
 81:   eps->trackall        = PETSC_FALSE;
 82:   eps->purify          = PETSC_TRUE;

 84:   eps->converged       = EPSConvergedRelative;
 85:   eps->convergeddestroy= NULL;
 86:   eps->stopping        = EPSStoppingBasic;
 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->rg              = NULL;
 98:   eps->D               = NULL;
 99:   eps->IS              = NULL;
100:   eps->defl            = NULL;
101:   eps->eigr            = NULL;
102:   eps->eigi            = NULL;
103:   eps->errest          = NULL;
104:   eps->rr              = NULL;
105:   eps->ri              = NULL;
106:   eps->perm            = NULL;
107:   eps->nwork           = 0;
108:   eps->work            = NULL;
109:   eps->data            = NULL;

111:   eps->state           = EPS_STATE_INITIAL;
112:   eps->nconv           = 0;
113:   eps->its             = 0;
114:   eps->nloc            = 0;
115:   eps->nrma            = 0.0;
116:   eps->nrmb            = 0.0;
117:   eps->isgeneralized   = PETSC_FALSE;
118:   eps->ispositive      = PETSC_FALSE;
119:   eps->ishermitian     = PETSC_FALSE;
120:   eps->reason          = EPS_CONVERGED_ITERATING;

122:   PetscNewLog(eps,&eps->sc);
123:   *outeps = eps;
124:   return(0);
125: }

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

132:    Logically Collective on EPS

134:    Input Parameters:
135: +  eps  - the eigensolver context
136: -  type - a known method

138:    Options Database Key:
139: .  -eps_type <method> - Sets the method; use -help for a list
140:     of available methods

142:    Notes:
143:    See "slepc/include/slepceps.h" for available methods. The default
144:    is EPSKRYLOVSCHUR.

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

154:    Level: intermediate

156: .seealso: STSetType(), EPSType
157: @*/
158: PetscErrorCode EPSSetType(EPS eps,EPSType type)
159: {
160:   PetscErrorCode ierr,(*r)(EPS);
161:   PetscBool      match;


167:   PetscObjectTypeCompare((PetscObject)eps,type,&match);
168:   if (match) return(0);

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

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

176:   eps->state = EPS_STATE_INITIAL;
177:   PetscObjectChangeTypeName((PetscObject)eps,type);
178:   (*r)(eps);
179:   return(0);
180: }

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

187:    Not Collective

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

192:    Output Parameter:
193: .  name - name of EPS method

195:    Level: intermediate

197: .seealso: EPSSetType()
198: @*/
199: PetscErrorCode EPSGetType(EPS eps,EPSType *type)
200: {
204:   *type = ((PetscObject)eps)->type_name;
205:   return(0);
206: }

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

213:    Not Collective

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

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

222:    Sample usage:
223: .vb
224:    EPSRegister("my_solver",MySolverCreate);
225: .ve

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

232:    Level: advanced

234: .seealso: EPSRegisterAll()
235: @*/
236: PetscErrorCode EPSRegister(const char *name,PetscErrorCode (*function)(EPS))
237: {

241:   PetscFunctionListAdd(&EPSList,name,function);
242:   return(0);
243: }

247: /*@
248:    EPSReset - Resets the EPS context to the initial state and removes any
249:    allocated objects.

251:    Collective on EPS

253:    Input Parameter:
254: .  eps - eigensolver context obtained from EPSCreate()

256:    Level: advanced

258: .seealso: EPSDestroy()
259: @*/
260: PetscErrorCode EPSReset(EPS eps)
261: {
263:   PetscInt       ncols;

267:   if (eps->ops->reset) { (eps->ops->reset)(eps); }
268:   if (eps->st) { STReset(eps->st); }
269:   if (eps->ds) { DSReset(eps->ds); }
270:   VecDestroy(&eps->D);
271:   BVGetSizes(eps->V,NULL,NULL,&ncols);
272:   if (ncols) {
273:     PetscFree4(eps->eigr,eps->eigi,eps->errest,eps->perm);
274:     PetscFree2(eps->rr,eps->ri);
275:   }
276:   BVDestroy(&eps->V);
277:   VecDestroyVecs(eps->nwork,&eps->work);
278:   eps->nwork = 0;
279:   eps->state = EPS_STATE_INITIAL;
280:   return(0);
281: }

285: /*@
286:    EPSDestroy - Destroys the EPS context.

288:    Collective on EPS

290:    Input Parameter:
291: .  eps - eigensolver context obtained from EPSCreate()

293:    Level: beginner

295: .seealso: EPSCreate(), EPSSetUp(), EPSSolve()
296: @*/
297: PetscErrorCode EPSDestroy(EPS *eps)
298: {

302:   if (!*eps) return(0);
304:   if (--((PetscObject)(*eps))->refct > 0) { *eps = 0; return(0); }
305:   EPSReset(*eps);
306:   if ((*eps)->ops->destroy) { (*(*eps)->ops->destroy)(*eps); }
307:   STDestroy(&(*eps)->st);
308:   RGDestroy(&(*eps)->rg);
309:   DSDestroy(&(*eps)->ds);
310:   PetscFree((*eps)->sc);
311:   /* just in case the initial vectors have not been used */
312:   SlepcBasisDestroy_Private(&(*eps)->nds,&(*eps)->defl);
313:   SlepcBasisDestroy_Private(&(*eps)->nini,&(*eps)->IS);
314:   if ((*eps)->convergeddestroy) {
315:     (*(*eps)->convergeddestroy)((*eps)->convergedctx);
316:   }
317:   EPSMonitorCancel(*eps);
318:   PetscHeaderDestroy(eps);
319:   return(0);
320: }

324: /*@
325:    EPSSetTarget - Sets the value of the target.

327:    Logically Collective on EPS

329:    Input Parameters:
330: +  eps    - eigensolver context
331: -  target - the value of the target

333:    Options Database Key:
334: .  -eps_target <scalar> - the value of the target

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

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

344:    Level: intermediate

346: .seealso: EPSGetTarget(), EPSSetWhichEigenpairs()
347: @*/
348: PetscErrorCode EPSSetTarget(EPS eps,PetscScalar target)
349: {

355:   eps->target = target;
356:   if (!eps->st) { EPSGetST(eps,&eps->st); }
357:   STSetDefaultShift(eps->st,target);
358:   return(0);
359: }

363: /*@
364:    EPSGetTarget - Gets the value of the target.

366:    Not Collective

368:    Input Parameter:
369: .  eps - eigensolver context

371:    Output Parameter:
372: .  target - the value of the target

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

377:    Level: intermediate

379: .seealso: EPSSetTarget()
380: @*/
381: PetscErrorCode EPSGetTarget(EPS eps,PetscScalar* target)
382: {
386:   *target = eps->target;
387:   return(0);
388: }

392: /*@
393:    EPSSetInterval - Defines the computational interval for spectrum slicing.

395:    Logically Collective on EPS

397:    Input Parameters:
398: +  eps  - eigensolver context
399: .  inta - left end of the interval
400: -  intb - right end of the interval

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

405:    Notes:
406:    Spectrum slicing is a technique employed for computing all eigenvalues of
407:    symmetric eigenproblems in a given interval. This function provides the
408:    interval to be considered. It must be used in combination with EPS_ALL, see
409:    EPSSetWhichEigenpairs().

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

416:    Level: intermediate

418: .seealso: EPSGetInterval(), EPSSetWhichEigenpairs()
419: @*/
420: PetscErrorCode EPSSetInterval(EPS eps,PetscReal inta,PetscReal intb)
421: {
426:   if (inta>=intb) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be inta<intb");
427:   eps->inta = inta;
428:   eps->intb = intb;
429:   eps->state = EPS_STATE_INITIAL;
430:   return(0);
431: }

435: /*@
436:    EPSGetInterval - Gets the computational interval for spectrum slicing.

438:    Not Collective

440:    Input Parameter:
441: .  eps - eigensolver context

443:    Output Parameters:
444: +  inta - left end of the interval
445: -  intb - right end of the interval

447:    Level: intermediate

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

452: .seealso: EPSSetInterval()
453: @*/
454: PetscErrorCode EPSGetInterval(EPS eps,PetscReal* inta,PetscReal* intb)
455: {
460:   if (inta) *inta = eps->inta;
461:   if (intb) *intb = eps->intb;
462:   return(0);
463: }

467: /*@
468:    EPSSetST - Associates a spectral transformation object to the eigensolver.

470:    Collective on EPS

472:    Input Parameters:
473: +  eps - eigensolver context obtained from EPSCreate()
474: -  st   - the spectral transformation object

476:    Note:
477:    Use EPSGetST() to retrieve the spectral transformation context (for example,
478:    to free it at the end of the computations).

480:    Level: advanced

482: .seealso: EPSGetST()
483: @*/
484: PetscErrorCode EPSSetST(EPS eps,ST st)
485: {

492:   PetscObjectReference((PetscObject)st);
493:   STDestroy(&eps->st);
494:   eps->st = st;
495:   PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->st);
496:   return(0);
497: }

501: /*@
502:    EPSGetST - Obtain the spectral transformation (ST) object associated
503:    to the eigensolver object.

505:    Not Collective

507:    Input Parameters:
508: .  eps - eigensolver context obtained from EPSCreate()

510:    Output Parameter:
511: .  st - spectral transformation context

513:    Level: intermediate

515: .seealso: EPSSetST()
516: @*/
517: PetscErrorCode EPSGetST(EPS eps,ST *st)
518: {

524:   if (!eps->st) {
525:     STCreate(PetscObjectComm((PetscObject)eps),&eps->st);
526:     PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->st);
527:   }
528:   *st = eps->st;
529:   return(0);
530: }

534: /*@
535:    EPSSetBV - Associates a basis vectors object to the eigensolver.

537:    Collective on EPS

539:    Input Parameters:
540: +  eps - eigensolver context obtained from EPSCreate()
541: -  V   - the basis vectors object

543:    Note:
544:    Use EPSGetBV() to retrieve the basis vectors context (for example,
545:    to free them at the end of the computations).

547:    Level: advanced

549: .seealso: EPSGetBV()
550: @*/
551: PetscErrorCode EPSSetBV(EPS eps,BV V)
552: {

559:   PetscObjectReference((PetscObject)V);
560:   BVDestroy(&eps->V);
561:   eps->V = V;
562:   PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->V);
563:   return(0);
564: }

568: /*@
569:    EPSGetBV - Obtain the basis vectors object associated to the eigensolver object.

571:    Not Collective

573:    Input Parameters:
574: .  eps - eigensolver context obtained from EPSCreate()

576:    Output Parameter:
577: .  V - basis vectors context

579:    Level: advanced

581: .seealso: EPSSetBV()
582: @*/
583: PetscErrorCode EPSGetBV(EPS eps,BV *V)
584: {

590:   if (!eps->V) {
591:     BVCreate(PetscObjectComm((PetscObject)eps),&eps->V);
592:     PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->V);
593:   }
594:   *V = eps->V;
595:   return(0);
596: }

600: /*@
601:    EPSSetRG - Associates a region object to the eigensolver.

603:    Collective on EPS

605:    Input Parameters:
606: +  eps - eigensolver context obtained from EPSCreate()
607: -  rg  - the region object

609:    Note:
610:    Use EPSGetRG() to retrieve the region context (for example,
611:    to free it at the end of the computations).

613:    Level: advanced

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

625:   PetscObjectReference((PetscObject)rg);
626:   RGDestroy(&eps->rg);
627:   eps->rg = rg;
628:   PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->rg);
629:   return(0);
630: }

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

637:    Not Collective

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

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

645:    Level: advanced

647: .seealso: EPSSetRG()
648: @*/
649: PetscErrorCode EPSGetRG(EPS eps,RG *rg)
650: {

656:   if (!eps->rg) {
657:     RGCreate(PetscObjectComm((PetscObject)eps),&eps->rg);
658:     PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->rg);
659:   }
660:   *rg = eps->rg;
661:   return(0);
662: }

666: /*@
667:    EPSSetDS - Associates a direct solver object to the eigensolver.

669:    Collective on EPS

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

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

679:    Level: advanced

681: .seealso: EPSGetDS()
682: @*/
683: PetscErrorCode EPSSetDS(EPS eps,DS ds)
684: {

691:   PetscObjectReference((PetscObject)ds);
692:   DSDestroy(&eps->ds);
693:   eps->ds = ds;
694:   PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->ds);
695:   return(0);
696: }

700: /*@
701:    EPSGetDS - Obtain the direct solver object associated to the eigensolver object.

703:    Not Collective

705:    Input Parameters:
706: .  eps - eigensolver context obtained from EPSCreate()

708:    Output Parameter:
709: .  ds - direct solver context

711:    Level: advanced

713: .seealso: EPSSetDS()
714: @*/
715: PetscErrorCode EPSGetDS(EPS eps,DS *ds)
716: {

722:   if (!eps->ds) {
723:     DSCreate(PetscObjectComm((PetscObject)eps),&eps->ds);
724:     PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->ds);
725:   }
726:   *ds = eps->ds;
727:   return(0);
728: }

732: /*@
733:    EPSIsGeneralized - Ask if the EPS object corresponds to a generalized
734:    eigenvalue problem.

736:    Not collective

738:    Input Parameter:
739: .  eps - the eigenproblem solver context

741:    Output Parameter:
742: .  is - the answer

744:    Level: intermediate

746: .seealso: EPSIsHermitian(), EPSIsPositive()
747: @*/
748: PetscErrorCode EPSIsGeneralized(EPS eps,PetscBool* is)
749: {
753:   *is = eps->isgeneralized;
754:   return(0);
755: }

759: /*@
760:    EPSIsHermitian - Ask if the EPS object corresponds to a Hermitian
761:    eigenvalue problem.

763:    Not collective

765:    Input Parameter:
766: .  eps - the eigenproblem solver context

768:    Output Parameter:
769: .  is - the answer

771:    Level: intermediate

773: .seealso: EPSIsGeneralized(), EPSIsPositive()
774: @*/
775: PetscErrorCode EPSIsHermitian(EPS eps,PetscBool* is)
776: {
780:   *is = eps->ishermitian;
781:   return(0);
782: }

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

790:    Not collective

792:    Input Parameter:
793: .  eps - the eigenproblem solver context

795:    Output Parameter:
796: .  is - the answer

798:    Level: intermediate

800: .seealso: EPSIsGeneralized(), EPSIsHermitian()
801: @*/
802: PetscErrorCode EPSIsPositive(EPS eps,PetscBool* is)
803: {
807:   *is = eps->ispositive;
808:   return(0);
809: }