Actual source code: pepbasic.c
slepc-3.6.1 2015-09-03
1: /*
2: The basic PEP routines, Create, Destroy, etc. are here.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, 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/pepimpl.h> /*I "slepcpep.h" I*/
26: PetscFunctionList PEPList = 0;
27: PetscBool PEPRegisterAllCalled = PETSC_FALSE;
28: PetscClassId PEP_CLASSID = 0;
29: PetscLogEvent PEP_SetUp = 0,PEP_Solve = 0,PEP_Refine = 0;
33: /*@
34: PEPCreate - Creates the default PEP context.
36: Collective on MPI_Comm
38: Input Parameter:
39: . comm - MPI communicator
41: Output Parameter:
42: . pep - location to put the PEP context
44: Note:
45: The default PEP type is PEPTOAR
47: Level: beginner
49: .seealso: PEPSetUp(), PEPSolve(), PEPDestroy(), PEP
50: @*/
51: PetscErrorCode PEPCreate(MPI_Comm comm,PEP *outpep)
52: {
54: PEP pep;
58: *outpep = 0;
59: PEPInitializePackage();
60: SlepcHeaderCreate(pep,PEP_CLASSID,"PEP","Polynomial Eigenvalue Problem","PEP",comm,PEPDestroy,PEPView);
62: pep->max_it = 0;
63: pep->nev = 1;
64: pep->ncv = 0;
65: pep->mpd = 0;
66: pep->nini = 0;
67: pep->target = 0.0;
68: pep->tol = PETSC_DEFAULT;
69: pep->conv = PEP_CONV_EIG;
70: pep->which = (PEPWhich)0;
71: pep->basis = PEP_BASIS_MONOMIAL;
72: pep->problem_type = (PEPProblemType)0;
73: pep->scale = PEP_SCALE_NONE;
74: pep->sfactor = 1.0;
75: pep->dsfactor = 1.0;
76: pep->sits = 5;
77: pep->slambda = 1.0;
78: pep->refine = PEP_REFINE_NONE;
79: pep->npart = 1;
80: pep->rtol = PETSC_DEFAULT;
81: pep->rits = PETSC_DEFAULT;
82: pep->schur = PETSC_FALSE;
83: pep->extract = (PEPExtract)0;
84: pep->trackall = PETSC_FALSE;
86: pep->converged = PEPConvergedEigRelative;
87: pep->convergeddestroy= NULL;
88: pep->convergedctx = NULL;
89: pep->numbermonitors = 0;
91: pep->st = NULL;
92: pep->ds = NULL;
93: pep->V = NULL;
94: pep->rg = NULL;
95: pep->rand = 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: PetscNewLog(pep,&pep->sc);
124: PetscRandomCreate(comm,&pep->rand);
125: PetscRandomSetSeed(pep->rand,0x12345678);
126: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->rand);
127: *outpep = pep;
128: return(0);
129: }
133: /*@C
134: PEPSetType - Selects the particular solver to be used in the PEP object.
136: Logically Collective on PEP
138: Input Parameters:
139: + pep - the polynomial eigensolver context
140: - type - a known method
142: Options Database Key:
143: . -pep_type <method> - Sets the method; use -help for a list
144: of available methods
146: Notes:
147: See "slepc/include/slepcpep.h" for available methods. The default
148: is PEPTOAR.
150: Normally, it is best to use the PEPSetFromOptions() command and
151: then set the PEP type from the options database rather than by using
152: this routine. Using the options database provides the user with
153: maximum flexibility in evaluating the different available methods.
154: The PEPSetType() routine is provided for those situations where it
155: is necessary to set the iterative solver independently of the command
156: line or options database.
158: Level: intermediate
160: .seealso: PEPType
161: @*/
162: PetscErrorCode PEPSetType(PEP pep,PEPType type)
163: {
164: PetscErrorCode ierr,(*r)(PEP);
165: PetscBool match;
171: PetscObjectTypeCompare((PetscObject)pep,type,&match);
172: if (match) return(0);
174: PetscFunctionListFind(PEPList,type,&r);
175: if (!r) SETERRQ1(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown PEP type given: %s",type);
177: if (pep->ops->destroy) { (*pep->ops->destroy)(pep); }
178: PetscMemzero(pep->ops,sizeof(struct _PEPOps));
180: pep->state = PEP_STATE_INITIAL;
181: PetscObjectChangeTypeName((PetscObject)pep,type);
182: (*r)(pep);
183: return(0);
184: }
188: /*@C
189: PEPGetType - Gets the PEP type as a string from the PEP object.
191: Not Collective
193: Input Parameter:
194: . pep - the eigensolver context
196: Output Parameter:
197: . name - name of PEP method
199: Level: intermediate
201: .seealso: PEPSetType()
202: @*/
203: PetscErrorCode PEPGetType(PEP pep,PEPType *type)
204: {
208: *type = ((PetscObject)pep)->type_name;
209: return(0);
210: }
214: /*@C
215: PEPRegister - Adds a method to the polynomial eigenproblem solver package.
217: Not Collective
219: Input Parameters:
220: + name - name of a new user-defined solver
221: - function - routine to create the solver context
223: Notes:
224: PEPRegister() may be called multiple times to add several user-defined solvers.
226: Sample usage:
227: .vb
228: PEPRegister("my_solver",MySolverCreate);
229: .ve
231: Then, your solver can be chosen with the procedural interface via
232: $ PEPSetType(pep,"my_solver")
233: or at runtime via the option
234: $ -pep_type my_solver
236: Level: advanced
238: .seealso: PEPRegisterAll()
239: @*/
240: PetscErrorCode PEPRegister(const char *name,PetscErrorCode (*function)(PEP))
241: {
245: PetscFunctionListAdd(&PEPList,name,function);
246: return(0);
247: }
251: /*@
252: PEPReset - Resets the PEP context to the initial state and removes any
253: allocated objects.
255: Collective on PEP
257: Input Parameter:
258: . pep - eigensolver context obtained from PEPCreate()
260: Level: advanced
262: .seealso: PEPDestroy()
263: @*/
264: PetscErrorCode PEPReset(PEP pep)
265: {
267: PetscInt ncols;
271: if (pep->ops->reset) { (pep->ops->reset)(pep); }
272: if (pep->st) { STReset(pep->st); }
273: if (pep->ds) { DSReset(pep->ds); }
274: if (pep->nmat) {
275: MatDestroyMatrices(pep->nmat,&pep->A);
276: PetscFree2(pep->pbc,pep->nrma);
277: PetscFree(pep->solvematcoeffs);
278: pep->nmat = 0;
279: }
280: VecDestroy(&pep->Dl);
281: VecDestroy(&pep->Dr);
282: BVGetSizes(pep->V,NULL,NULL,&ncols);
283: if (ncols) {
284: PetscFree4(pep->eigr,pep->eigi,pep->errest,pep->perm);
285: }
286: BVDestroy(&pep->V);
287: VecDestroyVecs(pep->nwork,&pep->work);
288: KSPDestroy(&pep->refineksp);
289: PetscSubcommDestroy(&pep->refinesubc);
290: pep->nwork = 0;
291: pep->state = PEP_STATE_INITIAL;
292: return(0);
293: }
297: /*@
298: PEPDestroy - Destroys the PEP context.
300: Collective on PEP
302: Input Parameter:
303: . pep - eigensolver context obtained from PEPCreate()
305: Level: beginner
307: .seealso: PEPCreate(), PEPSetUp(), PEPSolve()
308: @*/
309: PetscErrorCode PEPDestroy(PEP *pep)
310: {
314: if (!*pep) return(0);
316: if (--((PetscObject)(*pep))->refct > 0) { *pep = 0; return(0); }
317: PEPReset(*pep);
318: if ((*pep)->ops->destroy) { (*(*pep)->ops->destroy)(*pep); }
319: STDestroy(&(*pep)->st);
320: RGDestroy(&(*pep)->rg);
321: DSDestroy(&(*pep)->ds);
322: PetscRandomDestroy(&(*pep)->rand);
323: PetscFree((*pep)->sc);
324: /* just in case the initial vectors have not been used */
325: SlepcBasisDestroy_Private(&(*pep)->nini,&(*pep)->IS);
326: if ((*pep)->convergeddestroy) {
327: (*(*pep)->convergeddestroy)((*pep)->convergedctx);
328: }
329: PEPMonitorCancel(*pep);
330: PetscHeaderDestroy(pep);
331: return(0);
332: }
336: /*@
337: PEPSetBV - Associates a basis vectors object to the polynomial eigensolver.
339: Collective on PEP
341: Input Parameters:
342: + pep - eigensolver context obtained from PEPCreate()
343: - bv - the basis vectors object
345: Note:
346: Use PEPGetBV() to retrieve the basis vectors context (for example,
347: to free it at the end of the computations).
349: Level: advanced
351: .seealso: PEPGetBV()
352: @*/
353: PetscErrorCode PEPSetBV(PEP pep,BV bv)
354: {
361: PetscObjectReference((PetscObject)bv);
362: BVDestroy(&pep->V);
363: pep->V = bv;
364: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->V);
365: return(0);
366: }
370: /*@
371: PEPGetBV - Obtain the basis vectors object associated to the polynomial
372: eigensolver object.
374: Not Collective
376: Input Parameters:
377: . pep - eigensolver context obtained from PEPCreate()
379: Output Parameter:
380: . bv - basis vectors context
382: Level: advanced
384: .seealso: PEPSetBV()
385: @*/
386: PetscErrorCode PEPGetBV(PEP pep,BV *bv)
387: {
393: if (!pep->V) {
394: BVCreate(PetscObjectComm((PetscObject)pep),&pep->V);
395: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->V);
396: }
397: *bv = pep->V;
398: return(0);
399: }
403: /*@
404: PEPSetRG - Associates a region object to the polynomial eigensolver.
406: Collective on PEP
408: Input Parameters:
409: + pep - eigensolver context obtained from PEPCreate()
410: - rg - the region object
412: Note:
413: Use PEPGetRG() to retrieve the region context (for example,
414: to free it at the end of the computations).
416: Level: advanced
418: .seealso: PEPGetRG()
419: @*/
420: PetscErrorCode PEPSetRG(PEP pep,RG rg)
421: {
428: PetscObjectReference((PetscObject)rg);
429: RGDestroy(&pep->rg);
430: pep->rg = rg;
431: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->rg);
432: return(0);
433: }
437: /*@
438: PEPGetRG - Obtain the region object associated to the
439: polynomial eigensolver object.
441: Not Collective
443: Input Parameters:
444: . pep - eigensolver context obtained from PEPCreate()
446: Output Parameter:
447: . rg - region context
449: Level: advanced
451: .seealso: PEPSetRG()
452: @*/
453: PetscErrorCode PEPGetRG(PEP pep,RG *rg)
454: {
460: if (!pep->rg) {
461: RGCreate(PetscObjectComm((PetscObject)pep),&pep->rg);
462: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->rg);
463: }
464: *rg = pep->rg;
465: return(0);
466: }
470: /*@
471: PEPSetDS - Associates a direct solver object to the polynomial eigensolver.
473: Collective on PEP
475: Input Parameters:
476: + pep - eigensolver context obtained from PEPCreate()
477: - ds - the direct solver object
479: Note:
480: Use PEPGetDS() to retrieve the direct solver context (for example,
481: to free it at the end of the computations).
483: Level: advanced
485: .seealso: PEPGetDS()
486: @*/
487: PetscErrorCode PEPSetDS(PEP pep,DS ds)
488: {
495: PetscObjectReference((PetscObject)ds);
496: DSDestroy(&pep->ds);
497: pep->ds = ds;
498: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->ds);
499: return(0);
500: }
504: /*@
505: PEPGetDS - Obtain the direct solver object associated to the
506: polynomial eigensolver object.
508: Not Collective
510: Input Parameters:
511: . pep - eigensolver context obtained from PEPCreate()
513: Output Parameter:
514: . ds - direct solver context
516: Level: advanced
518: .seealso: PEPSetDS()
519: @*/
520: PetscErrorCode PEPGetDS(PEP pep,DS *ds)
521: {
527: if (!pep->ds) {
528: DSCreate(PetscObjectComm((PetscObject)pep),&pep->ds);
529: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->ds);
530: }
531: *ds = pep->ds;
532: return(0);
533: }
537: /*@
538: PEPSetST - Associates a spectral transformation object to the eigensolver.
540: Collective on PEP
542: Input Parameters:
543: + pep - eigensolver context obtained from PEPCreate()
544: - st - the spectral transformation object
546: Note:
547: Use PEPGetST() to retrieve the spectral transformation context (for example,
548: to free it at the end of the computations).
550: Level: advanced
552: .seealso: PEPGetST()
553: @*/
554: PetscErrorCode PEPSetST(PEP pep,ST st)
555: {
562: PetscObjectReference((PetscObject)st);
563: STDestroy(&pep->st);
564: pep->st = st;
565: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->st);
566: return(0);
567: }
571: /*@
572: PEPGetST - Obtain the spectral transformation (ST) object associated
573: to the eigensolver object.
575: Not Collective
577: Input Parameters:
578: . pep - eigensolver context obtained from PEPCreate()
580: Output Parameter:
581: . st - spectral transformation context
583: Level: intermediate
585: .seealso: PEPSetST()
586: @*/
587: PetscErrorCode PEPGetST(PEP pep,ST *st)
588: {
594: if (!pep->st) {
595: STCreate(PetscObjectComm((PetscObject)pep),&pep->st);
596: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->st);
597: }
598: *st = pep->st;
599: return(0);
600: }
604: /*@
605: PEPRefineGetKSP - Obtain the ksp object used by the eigensolver
606: object in the refinement phase.
608: Not Collective
610: Input Parameters:
611: . pep - eigensolver context obtained from PEPCreate()
613: Output Parameter:
614: . ksp - ksp context
616: Level: advanced
618: .seealso: PEPSetRefine()
619: @*/
620: PetscErrorCode PEPRefineGetKSP(PEP pep,KSP *ksp)
621: {
627: if (!pep->refineksp) {
628: if (pep->npart>1) {
629: /* Split in subcomunicators */
630: PetscSubcommCreate(PetscObjectComm((PetscObject)pep),&pep->refinesubc);
631: PetscSubcommSetNumber(pep->refinesubc,pep->npart);
632: PetscSubcommSetType(pep->refinesubc,PETSC_SUBCOMM_CONTIGUOUS);
633: PetscLogObjectMemory((PetscObject)pep,sizeof(PetscSubcomm));
634: }
635: KSPCreate((pep->npart==1)?PetscObjectComm((PetscObject)pep):PetscSubcommChild(pep->refinesubc),&pep->refineksp);
636: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->refineksp);
637: KSPSetOptionsPrefix(*ksp,((PetscObject)pep)->prefix);
638: KSPAppendOptionsPrefix(*ksp,"pep_refine_");
639: KSPSetErrorIfNotConverged(*ksp,PETSC_TRUE);
640: }
641: *ksp = pep->refineksp;
642: return(0);
643: }
647: /*@
648: PEPSetTarget - Sets the value of the target.
650: Logically Collective on PEP
652: Input Parameters:
653: + pep - eigensolver context
654: - target - the value of the target
656: Options Database Key:
657: . -pep_target <scalar> - the value of the target
659: Notes:
660: The target is a scalar value used to determine the portion of the spectrum
661: of interest. It is used in combination with PEPSetWhichEigenpairs().
663: In the case of complex scalars, a complex value can be provided in the
664: command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
665: -pep_target 1.0+2.0i
667: Level: intermediate
669: .seealso: PEPGetTarget(), PEPSetWhichEigenpairs()
670: @*/
671: PetscErrorCode PEPSetTarget(PEP pep,PetscScalar target)
672: {
678: pep->target = target;
679: if (!pep->st) { PEPGetST(pep,&pep->st); }
680: STSetDefaultShift(pep->st,target);
681: return(0);
682: }
686: /*@
687: PEPGetTarget - Gets the value of the target.
689: Not Collective
691: Input Parameter:
692: . pep - eigensolver context
694: Output Parameter:
695: . target - the value of the target
697: Note:
698: If the target was not set by the user, then zero is returned.
700: Level: intermediate
702: .seealso: PEPSetTarget()
703: @*/
704: PetscErrorCode PEPGetTarget(PEP pep,PetscScalar* target)
705: {
709: *target = pep->target;
710: return(0);
711: }