Actual source code: pepbasic.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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>      /*I "slepcpep.h" I*/

 16: PetscFunctionList PEPList = 0;
 17: PetscBool         PEPRegisterAllCalled = PETSC_FALSE;
 18: PetscClassId      PEP_CLASSID = 0;
 19: PetscLogEvent     PEP_SetUp = 0,PEP_Solve = 0,PEP_Refine = 0;

 21: /*@
 22:    PEPCreate - Creates the default PEP context.

 24:    Collective on MPI_Comm

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

 29:    Output Parameter:
 30: .  pep - location to put the PEP context

 32:    Note:
 33:    The default PEP type is PEPTOAR

 35:    Level: beginner

 37: .seealso: PEPSetUp(), PEPSolve(), PEPDestroy(), PEP
 38: @*/
 39: PetscErrorCode PEPCreate(MPI_Comm comm,PEP *outpep)
 40: {
 42:   PEP            pep;

 46:   *outpep = 0;
 47:   PEPInitializePackage();
 48:   SlepcHeaderCreate(pep,PEP_CLASSID,"PEP","Polynomial Eigenvalue Problem","PEP",comm,PEPDestroy,PEPView);

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

 75:   pep->converged       = PEPConvergedRelative;
 76:   pep->convergeduser   = NULL;
 77:   pep->convergeddestroy= NULL;
 78:   pep->stopping        = PEPStoppingBasic;
 79:   pep->stoppinguser    = NULL;
 80:   pep->stoppingdestroy = NULL;
 81:   pep->convergedctx    = NULL;
 82:   pep->stoppingctx     = NULL;
 83:   pep->numbermonitors  = 0;

 85:   pep->st              = NULL;
 86:   pep->ds              = NULL;
 87:   pep->V               = NULL;
 88:   pep->rg              = NULL;
 89:   pep->A               = NULL;
 90:   pep->nmat            = 0;
 91:   pep->Dl              = NULL;
 92:   pep->Dr              = NULL;
 93:   pep->IS              = NULL;
 94:   pep->eigr            = NULL;
 95:   pep->eigi            = NULL;
 96:   pep->errest          = NULL;
 97:   pep->perm            = NULL;
 98:   pep->pbc             = NULL;
 99:   pep->solvematcoeffs  = NULL;
100:   pep->nwork           = 0;
101:   pep->work            = NULL;
102:   pep->refineksp       = NULL;
103:   pep->refinesubc      = NULL;
104:   pep->data            = NULL;

106:   pep->state           = PEP_STATE_INITIAL;
107:   pep->nconv           = 0;
108:   pep->its             = 0;
109:   pep->n               = 0;
110:   pep->nloc            = 0;
111:   pep->nrma            = NULL;
112:   pep->sfactor_set     = PETSC_FALSE;
113:   pep->lineariz        = PETSC_FALSE;
114:   pep->reason          = PEP_CONVERGED_ITERATING;

116:   PetscNewLog(pep,&pep->sc);
117:   *outpep = pep;
118:   return(0);
119: }

121: /*@C
122:    PEPSetType - Selects the particular solver to be used in the PEP object.

124:    Logically Collective on PEP

126:    Input Parameters:
127: +  pep      - the polynomial eigensolver context
128: -  type     - a known method

130:    Options Database Key:
131: .  -pep_type <method> - Sets the method; use -help for a list
132:     of available methods

134:    Notes:
135:    See "slepc/include/slepcpep.h" for available methods. The default
136:    is PEPTOAR.

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

146:    Level: intermediate

148: .seealso: PEPType
149: @*/
150: PetscErrorCode PEPSetType(PEP pep,PEPType type)
151: {
152:   PetscErrorCode ierr,(*r)(PEP);
153:   PetscBool      match;


159:   PetscObjectTypeCompare((PetscObject)pep,type,&match);
160:   if (match) return(0);

162:   PetscFunctionListFind(PEPList,type,&r);
163:   if (!r) SETERRQ1(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown PEP type given: %s",type);

165:   if (pep->ops->destroy) { (*pep->ops->destroy)(pep); }
166:   PetscMemzero(pep->ops,sizeof(struct _PEPOps));

168:   pep->state = PEP_STATE_INITIAL;
169:   PetscObjectChangeTypeName((PetscObject)pep,type);
170:   (*r)(pep);
171:   return(0);
172: }

174: /*@C
175:    PEPGetType - Gets the PEP type as a string from the PEP object.

177:    Not Collective

179:    Input Parameter:
180: .  pep - the eigensolver context

182:    Output Parameter:
183: .  name - name of PEP method

185:    Level: intermediate

187: .seealso: PEPSetType()
188: @*/
189: PetscErrorCode PEPGetType(PEP pep,PEPType *type)
190: {
194:   *type = ((PetscObject)pep)->type_name;
195:   return(0);
196: }

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

201:    Not Collective

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

207:    Notes:
208:    PEPRegister() may be called multiple times to add several user-defined solvers.

210:    Sample usage:
211: .vb
212:     PEPRegister("my_solver",MySolverCreate);
213: .ve

215:    Then, your solver can be chosen with the procedural interface via
216: $     PEPSetType(pep,"my_solver")
217:    or at runtime via the option
218: $     -pep_type my_solver

220:    Level: advanced

222: .seealso: PEPRegisterAll()
223: @*/
224: PetscErrorCode PEPRegister(const char *name,PetscErrorCode (*function)(PEP))
225: {

229:   PetscFunctionListAdd(&PEPList,name,function);
230:   return(0);
231: }

233: /*@
234:    PEPReset - Resets the PEP context to the initial state (prior to setup)
235:    and destroys any allocated Vecs and Mats.

237:    Collective on PEP

239:    Input Parameter:
240: .  pep - eigensolver context obtained from PEPCreate()

242:    Level: advanced

244: .seealso: PEPDestroy()
245: @*/
246: PetscErrorCode PEPReset(PEP pep)
247: {

252:   if (!pep) return(0);
253:   if (pep->ops->reset) { (pep->ops->reset)(pep); }
254:   if (pep->st) { STReset(pep->st); }
255:   if (pep->refineksp) { KSPReset(pep->refineksp); }
256:   if (pep->nmat) {
257:     MatDestroyMatrices(pep->nmat,&pep->A);
258:     PetscFree2(pep->pbc,pep->nrma);
259:     PetscFree(pep->solvematcoeffs);
260:     pep->nmat = 0;
261:   }
262:   VecDestroy(&pep->Dl);
263:   VecDestroy(&pep->Dr);
264:   BVDestroy(&pep->V);
265:   VecDestroyVecs(pep->nwork,&pep->work);
266:   pep->nwork = 0;
267:   pep->state = PEP_STATE_INITIAL;
268:   return(0);
269: }

271: /*@
272:    PEPDestroy - Destroys the PEP context.

274:    Collective on PEP

276:    Input Parameter:
277: .  pep - eigensolver context obtained from PEPCreate()

279:    Level: beginner

281: .seealso: PEPCreate(), PEPSetUp(), PEPSolve()
282: @*/
283: PetscErrorCode PEPDestroy(PEP *pep)
284: {

288:   if (!*pep) return(0);
290:   if (--((PetscObject)(*pep))->refct > 0) { *pep = 0; return(0); }
291:   PEPReset(*pep);
292:   if ((*pep)->ops->destroy) { (*(*pep)->ops->destroy)(*pep); }
293:   if ((*pep)->eigr) {
294:     PetscFree4((*pep)->eigr,(*pep)->eigi,(*pep)->errest,(*pep)->perm);
295:   }
296:   STDestroy(&(*pep)->st);
297:   RGDestroy(&(*pep)->rg);
298:   DSDestroy(&(*pep)->ds);
299:   KSPDestroy(&(*pep)->refineksp);
300:   PetscSubcommDestroy(&(*pep)->refinesubc);
301:   PetscFree((*pep)->sc);
302:   /* just in case the initial vectors have not been used */
303:   SlepcBasisDestroy_Private(&(*pep)->nini,&(*pep)->IS);
304:   if ((*pep)->convergeddestroy) {
305:     (*(*pep)->convergeddestroy)((*pep)->convergedctx);
306:   }
307:   PEPMonitorCancel(*pep);
308:   PetscHeaderDestroy(pep);
309:   return(0);
310: }

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

315:    Collective on PEP

317:    Input Parameters:
318: +  pep - eigensolver context obtained from PEPCreate()
319: -  bv  - the basis vectors object

321:    Note:
322:    Use PEPGetBV() to retrieve the basis vectors context (for example,
323:    to free it at the end of the computations).

325:    Level: advanced

327: .seealso: PEPGetBV()
328: @*/
329: PetscErrorCode PEPSetBV(PEP pep,BV bv)
330: {

337:   PetscObjectReference((PetscObject)bv);
338:   BVDestroy(&pep->V);
339:   pep->V = bv;
340:   PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->V);
341:   return(0);
342: }

344: /*@
345:    PEPGetBV - Obtain the basis vectors object associated to the polynomial
346:    eigensolver object.

348:    Not Collective

350:    Input Parameters:
351: .  pep - eigensolver context obtained from PEPCreate()

353:    Output Parameter:
354: .  bv - basis vectors context

356:    Level: advanced

358: .seealso: PEPSetBV()
359: @*/
360: PetscErrorCode PEPGetBV(PEP pep,BV *bv)
361: {

367:   if (!pep->V) {
368:     BVCreate(PetscObjectComm((PetscObject)pep),&pep->V);
369:     PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->V);
370:   }
371:   *bv = pep->V;
372:   return(0);
373: }

375: /*@
376:    PEPSetRG - Associates a region object to the polynomial eigensolver.

378:    Collective on PEP

380:    Input Parameters:
381: +  pep - eigensolver context obtained from PEPCreate()
382: -  rg  - the region object

384:    Note:
385:    Use PEPGetRG() to retrieve the region context (for example,
386:    to free it at the end of the computations).

388:    Level: advanced

390: .seealso: PEPGetRG()
391: @*/
392: PetscErrorCode PEPSetRG(PEP pep,RG rg)
393: {

400:   PetscObjectReference((PetscObject)rg);
401:   RGDestroy(&pep->rg);
402:   pep->rg = rg;
403:   PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->rg);
404:   return(0);
405: }

407: /*@
408:    PEPGetRG - Obtain the region object associated to the
409:    polynomial eigensolver object.

411:    Not Collective

413:    Input Parameters:
414: .  pep - eigensolver context obtained from PEPCreate()

416:    Output Parameter:
417: .  rg - region context

419:    Level: advanced

421: .seealso: PEPSetRG()
422: @*/
423: PetscErrorCode PEPGetRG(PEP pep,RG *rg)
424: {

430:   if (!pep->rg) {
431:     RGCreate(PetscObjectComm((PetscObject)pep),&pep->rg);
432:     PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->rg);
433:   }
434:   *rg = pep->rg;
435:   return(0);
436: }

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

441:    Collective on PEP

443:    Input Parameters:
444: +  pep - eigensolver context obtained from PEPCreate()
445: -  ds  - the direct solver object

447:    Note:
448:    Use PEPGetDS() to retrieve the direct solver context (for example,
449:    to free it at the end of the computations).

451:    Level: advanced

453: .seealso: PEPGetDS()
454: @*/
455: PetscErrorCode PEPSetDS(PEP pep,DS ds)
456: {

463:   PetscObjectReference((PetscObject)ds);
464:   DSDestroy(&pep->ds);
465:   pep->ds = ds;
466:   PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->ds);
467:   return(0);
468: }

470: /*@
471:    PEPGetDS - Obtain the direct solver object associated to the
472:    polynomial eigensolver object.

474:    Not Collective

476:    Input Parameters:
477: .  pep - eigensolver context obtained from PEPCreate()

479:    Output Parameter:
480: .  ds - direct solver context

482:    Level: advanced

484: .seealso: PEPSetDS()
485: @*/
486: PetscErrorCode PEPGetDS(PEP pep,DS *ds)
487: {

493:   if (!pep->ds) {
494:     DSCreate(PetscObjectComm((PetscObject)pep),&pep->ds);
495:     PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->ds);
496:   }
497:   *ds = pep->ds;
498:   return(0);
499: }

501: /*@
502:    PEPSetST - Associates a spectral transformation object to the eigensolver.

504:    Collective on PEP

506:    Input Parameters:
507: +  pep - eigensolver context obtained from PEPCreate()
508: -  st   - the spectral transformation object

510:    Note:
511:    Use PEPGetST() to retrieve the spectral transformation context (for example,
512:    to free it at the end of the computations).

514:    Level: advanced

516: .seealso: PEPGetST()
517: @*/
518: PetscErrorCode PEPSetST(PEP pep,ST st)
519: {

526:   PetscObjectReference((PetscObject)st);
527:   STDestroy(&pep->st);
528:   pep->st = st;
529:   PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->st);
530:   return(0);
531: }

533: /*@
534:    PEPGetST - Obtain the spectral transformation (ST) object associated
535:    to the eigensolver object.

537:    Not Collective

539:    Input Parameters:
540: .  pep - eigensolver context obtained from PEPCreate()

542:    Output Parameter:
543: .  st - spectral transformation context

545:    Level: intermediate

547: .seealso: PEPSetST()
548: @*/
549: PetscErrorCode PEPGetST(PEP pep,ST *st)
550: {

556:   if (!pep->st) {
557:     STCreate(PetscObjectComm((PetscObject)pep),&pep->st);
558:     PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->st);
559:   }
560:   *st = pep->st;
561:   return(0);
562: }

564: /*@
565:    PEPRefineGetKSP - Obtain the ksp object used by the eigensolver
566:    object in the refinement phase.

568:    Not Collective

570:    Input Parameters:
571: .  pep - eigensolver context obtained from PEPCreate()

573:    Output Parameter:
574: .  ksp - ksp context

576:    Level: advanced

578: .seealso: PEPSetRefine()
579: @*/
580: PetscErrorCode PEPRefineGetKSP(PEP pep,KSP *ksp)
581: {

587:   if (!pep->refineksp) {
588:     if (pep->npart>1) {
589:       /* Split in subcomunicators */
590:       PetscSubcommCreate(PetscObjectComm((PetscObject)pep),&pep->refinesubc);
591:       PetscSubcommSetNumber(pep->refinesubc,pep->npart);
592:       PetscSubcommSetType(pep->refinesubc,PETSC_SUBCOMM_CONTIGUOUS);
593:       PetscLogObjectMemory((PetscObject)pep,sizeof(PetscSubcomm));
594:     }
595:     KSPCreate((pep->npart==1)?PetscObjectComm((PetscObject)pep):PetscSubcommChild(pep->refinesubc),&pep->refineksp);
596:     PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->refineksp);
597:     KSPSetOptionsPrefix(*ksp,((PetscObject)pep)->prefix);
598:     KSPAppendOptionsPrefix(*ksp,"pep_refine_");
599:     KSPSetTolerances(pep->refineksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
600:   }
601:   *ksp = pep->refineksp;
602:   return(0);
603: }

605: /*@
606:    PEPSetTarget - Sets the value of the target.

608:    Logically Collective on PEP

610:    Input Parameters:
611: +  pep    - eigensolver context
612: -  target - the value of the target

614:    Options Database Key:
615: .  -pep_target <scalar> - the value of the target

617:    Notes:
618:    The target is a scalar value used to determine the portion of the spectrum
619:    of interest. It is used in combination with PEPSetWhichEigenpairs().

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

625:    Level: intermediate

627: .seealso: PEPGetTarget(), PEPSetWhichEigenpairs()
628: @*/
629: PetscErrorCode PEPSetTarget(PEP pep,PetscScalar target)
630: {

636:   pep->target = target;
637:   if (!pep->st) { PEPGetST(pep,&pep->st); }
638:   STSetDefaultShift(pep->st,target);
639:   return(0);
640: }

642: /*@
643:    PEPGetTarget - Gets the value of the target.

645:    Not Collective

647:    Input Parameter:
648: .  pep - eigensolver context

650:    Output Parameter:
651: .  target - the value of the target

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

656:    Level: intermediate

658: .seealso: PEPSetTarget()
659: @*/
660: PetscErrorCode PEPGetTarget(PEP pep,PetscScalar* target)
661: {
665:   *target = pep->target;
666:   return(0);
667: }