Actual source code: pepbasic.c
slepc-3.7.2 2016-07-19
1: /*
2: The basic PEP 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/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_REL;
70: pep->stop = PEP_STOP_BASIC;
71: pep->which = (PEPWhich)0;
72: pep->basis = PEP_BASIS_MONOMIAL;
73: pep->problem_type = (PEPProblemType)0;
74: pep->scale = PEP_SCALE_NONE;
75: pep->sfactor = 1.0;
76: pep->dsfactor = 1.0;
77: pep->sits = 5;
78: pep->slambda = 1.0;
79: pep->refine = PEP_REFINE_NONE;
80: pep->npart = 1;
81: pep->rtol = PETSC_DEFAULT;
82: pep->rits = PETSC_DEFAULT;
83: pep->scheme = (PEPRefineScheme)0;
84: pep->extract = (PEPExtract)0;
85: pep->trackall = PETSC_FALSE;
87: pep->converged = PEPConvergedRelative;
88: pep->convergeddestroy= NULL;
89: pep->stopping = PEPStoppingBasic;
90: pep->stoppingdestroy = NULL;
91: pep->convergedctx = NULL;
92: pep->stoppingctx = NULL;
93: pep->numbermonitors = 0;
95: pep->st = NULL;
96: pep->ds = NULL;
97: pep->V = NULL;
98: pep->rg = NULL;
99: pep->A = NULL;
100: pep->nmat = 0;
101: pep->Dl = NULL;
102: pep->Dr = NULL;
103: pep->IS = NULL;
104: pep->eigr = NULL;
105: pep->eigi = NULL;
106: pep->errest = NULL;
107: pep->perm = NULL;
108: pep->pbc = NULL;
109: pep->solvematcoeffs = NULL;
110: pep->nwork = 0;
111: pep->work = NULL;
112: pep->refineksp = NULL;
113: pep->refinesubc = NULL;
114: pep->data = NULL;
116: pep->state = PEP_STATE_INITIAL;
117: pep->nconv = 0;
118: pep->its = 0;
119: pep->n = 0;
120: pep->nloc = 0;
121: pep->nrma = NULL;
122: pep->sfactor_set = PETSC_FALSE;
123: pep->lineariz = PETSC_FALSE;
124: pep->reason = PEP_CONVERGED_ITERATING;
126: PetscNewLog(pep,&pep->sc);
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: PetscFree((*pep)->sc);
323: /* just in case the initial vectors have not been used */
324: SlepcBasisDestroy_Private(&(*pep)->nini,&(*pep)->IS);
325: if ((*pep)->convergeddestroy) {
326: (*(*pep)->convergeddestroy)((*pep)->convergedctx);
327: }
328: PEPMonitorCancel(*pep);
329: PetscHeaderDestroy(pep);
330: return(0);
331: }
335: /*@
336: PEPSetBV - Associates a basis vectors object to the polynomial eigensolver.
338: Collective on PEP
340: Input Parameters:
341: + pep - eigensolver context obtained from PEPCreate()
342: - bv - the basis vectors object
344: Note:
345: Use PEPGetBV() to retrieve the basis vectors context (for example,
346: to free it at the end of the computations).
348: Level: advanced
350: .seealso: PEPGetBV()
351: @*/
352: PetscErrorCode PEPSetBV(PEP pep,BV bv)
353: {
360: PetscObjectReference((PetscObject)bv);
361: BVDestroy(&pep->V);
362: pep->V = bv;
363: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->V);
364: return(0);
365: }
369: /*@
370: PEPGetBV - Obtain the basis vectors object associated to the polynomial
371: eigensolver object.
373: Not Collective
375: Input Parameters:
376: . pep - eigensolver context obtained from PEPCreate()
378: Output Parameter:
379: . bv - basis vectors context
381: Level: advanced
383: .seealso: PEPSetBV()
384: @*/
385: PetscErrorCode PEPGetBV(PEP pep,BV *bv)
386: {
392: if (!pep->V) {
393: BVCreate(PetscObjectComm((PetscObject)pep),&pep->V);
394: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->V);
395: }
396: *bv = pep->V;
397: return(0);
398: }
402: /*@
403: PEPSetRG - Associates a region object to the polynomial eigensolver.
405: Collective on PEP
407: Input Parameters:
408: + pep - eigensolver context obtained from PEPCreate()
409: - rg - the region object
411: Note:
412: Use PEPGetRG() to retrieve the region context (for example,
413: to free it at the end of the computations).
415: Level: advanced
417: .seealso: PEPGetRG()
418: @*/
419: PetscErrorCode PEPSetRG(PEP pep,RG rg)
420: {
427: PetscObjectReference((PetscObject)rg);
428: RGDestroy(&pep->rg);
429: pep->rg = rg;
430: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->rg);
431: return(0);
432: }
436: /*@
437: PEPGetRG - Obtain the region object associated to the
438: polynomial eigensolver object.
440: Not Collective
442: Input Parameters:
443: . pep - eigensolver context obtained from PEPCreate()
445: Output Parameter:
446: . rg - region context
448: Level: advanced
450: .seealso: PEPSetRG()
451: @*/
452: PetscErrorCode PEPGetRG(PEP pep,RG *rg)
453: {
459: if (!pep->rg) {
460: RGCreate(PetscObjectComm((PetscObject)pep),&pep->rg);
461: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->rg);
462: }
463: *rg = pep->rg;
464: return(0);
465: }
469: /*@
470: PEPSetDS - Associates a direct solver object to the polynomial eigensolver.
472: Collective on PEP
474: Input Parameters:
475: + pep - eigensolver context obtained from PEPCreate()
476: - ds - the direct solver object
478: Note:
479: Use PEPGetDS() to retrieve the direct solver context (for example,
480: to free it at the end of the computations).
482: Level: advanced
484: .seealso: PEPGetDS()
485: @*/
486: PetscErrorCode PEPSetDS(PEP pep,DS ds)
487: {
494: PetscObjectReference((PetscObject)ds);
495: DSDestroy(&pep->ds);
496: pep->ds = ds;
497: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->ds);
498: return(0);
499: }
503: /*@
504: PEPGetDS - Obtain the direct solver object associated to the
505: polynomial eigensolver object.
507: Not Collective
509: Input Parameters:
510: . pep - eigensolver context obtained from PEPCreate()
512: Output Parameter:
513: . ds - direct solver context
515: Level: advanced
517: .seealso: PEPSetDS()
518: @*/
519: PetscErrorCode PEPGetDS(PEP pep,DS *ds)
520: {
526: if (!pep->ds) {
527: DSCreate(PetscObjectComm((PetscObject)pep),&pep->ds);
528: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->ds);
529: }
530: *ds = pep->ds;
531: return(0);
532: }
536: /*@
537: PEPSetST - Associates a spectral transformation object to the eigensolver.
539: Collective on PEP
541: Input Parameters:
542: + pep - eigensolver context obtained from PEPCreate()
543: - st - the spectral transformation object
545: Note:
546: Use PEPGetST() to retrieve the spectral transformation context (for example,
547: to free it at the end of the computations).
549: Level: advanced
551: .seealso: PEPGetST()
552: @*/
553: PetscErrorCode PEPSetST(PEP pep,ST st)
554: {
561: PetscObjectReference((PetscObject)st);
562: STDestroy(&pep->st);
563: pep->st = st;
564: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->st);
565: return(0);
566: }
570: /*@
571: PEPGetST - Obtain the spectral transformation (ST) object associated
572: to the eigensolver object.
574: Not Collective
576: Input Parameters:
577: . pep - eigensolver context obtained from PEPCreate()
579: Output Parameter:
580: . st - spectral transformation context
582: Level: intermediate
584: .seealso: PEPSetST()
585: @*/
586: PetscErrorCode PEPGetST(PEP pep,ST *st)
587: {
593: if (!pep->st) {
594: STCreate(PetscObjectComm((PetscObject)pep),&pep->st);
595: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->st);
596: }
597: *st = pep->st;
598: return(0);
599: }
603: /*@
604: PEPRefineGetKSP - Obtain the ksp object used by the eigensolver
605: object in the refinement phase.
607: Not Collective
609: Input Parameters:
610: . pep - eigensolver context obtained from PEPCreate()
612: Output Parameter:
613: . ksp - ksp context
615: Level: advanced
617: .seealso: PEPSetRefine()
618: @*/
619: PetscErrorCode PEPRefineGetKSP(PEP pep,KSP *ksp)
620: {
626: if (!pep->refineksp) {
627: if (pep->npart>1) {
628: /* Split in subcomunicators */
629: PetscSubcommCreate(PetscObjectComm((PetscObject)pep),&pep->refinesubc);
630: PetscSubcommSetNumber(pep->refinesubc,pep->npart);
631: PetscSubcommSetType(pep->refinesubc,PETSC_SUBCOMM_CONTIGUOUS);
632: PetscLogObjectMemory((PetscObject)pep,sizeof(PetscSubcomm));
633: }
634: KSPCreate((pep->npart==1)?PetscObjectComm((PetscObject)pep):PetscSubcommChild(pep->refinesubc),&pep->refineksp);
635: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->refineksp);
636: KSPSetOptionsPrefix(*ksp,((PetscObject)pep)->prefix);
637: KSPAppendOptionsPrefix(*ksp,"pep_refine_");
638: }
639: *ksp = pep->refineksp;
640: return(0);
641: }
645: /*@
646: PEPSetTarget - Sets the value of the target.
648: Logically Collective on PEP
650: Input Parameters:
651: + pep - eigensolver context
652: - target - the value of the target
654: Options Database Key:
655: . -pep_target <scalar> - the value of the target
657: Notes:
658: The target is a scalar value used to determine the portion of the spectrum
659: of interest. It is used in combination with PEPSetWhichEigenpairs().
661: In the case of complex scalars, a complex value can be provided in the
662: command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
663: -pep_target 1.0+2.0i
665: Level: intermediate
667: .seealso: PEPGetTarget(), PEPSetWhichEigenpairs()
668: @*/
669: PetscErrorCode PEPSetTarget(PEP pep,PetscScalar target)
670: {
676: pep->target = target;
677: if (!pep->st) { PEPGetST(pep,&pep->st); }
678: STSetDefaultShift(pep->st,target);
679: return(0);
680: }
684: /*@
685: PEPGetTarget - Gets the value of the target.
687: Not Collective
689: Input Parameter:
690: . pep - eigensolver context
692: Output Parameter:
693: . target - the value of the target
695: Note:
696: If the target was not set by the user, then zero is returned.
698: Level: intermediate
700: .seealso: PEPSetTarget()
701: @*/
702: PetscErrorCode PEPGetTarget(PEP pep,PetscScalar* target)
703: {
707: *target = pep->target;
708: return(0);
709: }