1: /*
2: The basic EPS 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/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_EIG;
71: eps->which = (EPSWhich)0;
72: eps->inta = 0.0;
73: eps->intb = 0.0;
74: eps->problem_type = (EPSProblemType)0;
75: eps->extraction = EPS_RITZ;
76: eps->balance = EPS_BALANCE_NONE;
77: eps->balance_its = 5;
78: eps->balance_cutoff = 1e-8;
79: eps->trueres = PETSC_FALSE;
80: eps->trackall = PETSC_FALSE;
81: eps->purify = PETSC_TRUE;
83: eps->converged = EPSConvergedEigRelative;
84: eps->convergeddestroy= NULL;
85: eps->arbitrary = NULL;
86: eps->convergedctx = NULL;
87: eps->arbitraryctx = NULL;
88: eps->numbermonitors = 0;
90: eps->st = NULL;
91: eps->ds = NULL;
92: eps->V = NULL;
93: eps->rg = NULL;
94: eps->rand = NULL;
95: eps->D = NULL;
96: eps->IS = NULL;
97: eps->defl = NULL;
98: eps->eigr = NULL;
99: eps->eigi = NULL;
100: eps->errest = NULL;
101: eps->rr = NULL;
102: eps->ri = NULL;
103: eps->perm = NULL;
104: eps->nwork = 0;
105: eps->work = NULL;
106: eps->data = NULL;
108: eps->state = EPS_STATE_INITIAL;
109: eps->nconv = 0;
110: eps->its = 0;
111: eps->nloc = 0;
112: eps->nrma = 0.0;
113: eps->nrmb = 0.0;
114: eps->isgeneralized = PETSC_FALSE;
115: eps->ispositive = PETSC_FALSE;
116: eps->ishermitian = PETSC_FALSE;
117: eps->reason = EPS_CONVERGED_ITERATING;
119: PetscNewLog(eps,&eps->sc);
120: PetscRandomCreate(comm,&eps->rand);
121: PetscRandomSetSeed(eps->rand,0x12345678);
122: PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->rand);
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 EPS134: 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(), EPSType157: @*/
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 EPS253: 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 EPS290: 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: PetscRandomDestroy(&(*eps)->rand);
311: PetscFree((*eps)->sc);
312: /* just in case the initial vectors have not been used */
313: SlepcBasisDestroy_Private(&(*eps)->nds,&(*eps)->defl);
314: SlepcBasisDestroy_Private(&(*eps)->nini,&(*eps)->IS);
315: if ((*eps)->convergeddestroy) {
316: (*(*eps)->convergeddestroy)((*eps)->convergedctx);
317: }
318: EPSMonitorCancel(*eps);
319: PetscHeaderDestroy(eps);
320: return(0);
321: }
325: /*@
326: EPSSetTarget - Sets the value of the target.
328: Logically Collective on EPS330: Input Parameters:
331: + eps - eigensolver context
332: - target - the value of the target
334: Options Database Key:
335: . -eps_target <scalar> - the value of the target
337: Notes:
338: The target is a scalar value used to determine the portion of the spectrum
339: of interest. It is used in combination with EPSSetWhichEigenpairs().
341: In the case of complex scalars, a complex value can be provided in the
342: command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
343: -eps_target 1.0+2.0i
345: Level: intermediate
347: .seealso: EPSGetTarget(), EPSSetWhichEigenpairs()
348: @*/
349: PetscErrorCode EPSSetTarget(EPS eps,PetscScalar target)350: {
356: eps->target = target;
357: if (!eps->st) { EPSGetST(eps,&eps->st); }
358: STSetDefaultShift(eps->st,target);
359: return(0);
360: }
364: /*@
365: EPSGetTarget - Gets the value of the target.
367: Not Collective
369: Input Parameter:
370: . eps - eigensolver context
372: Output Parameter:
373: . target - the value of the target
375: Note:
376: If the target was not set by the user, then zero is returned.
378: Level: intermediate
380: .seealso: EPSSetTarget()
381: @*/
382: PetscErrorCode EPSGetTarget(EPS eps,PetscScalar* target)383: {
387: *target = eps->target;
388: return(0);
389: }
393: /*@
394: EPSSetInterval - Defines the computational interval for spectrum slicing.
396: Logically Collective on EPS398: Input Parameters:
399: + eps - eigensolver context
400: . inta - left end of the interval
401: - intb - right end of the interval
403: Options Database Key:
404: . -eps_interval <a,b> - set [a,b] as the interval of interest
406: Notes:
407: Spectrum slicing is a technique employed for computing all eigenvalues of
408: symmetric eigenproblems in a given interval. This function provides the
409: interval to be considered. It must be used in combination with EPS_ALL, see
410: EPSSetWhichEigenpairs().
412: In the command-line option, two values must be provided. For an open interval,
413: one can give an infinite, e.g., -eps_interval 1.0,inf or -eps_interval -inf,1.0.
414: An open interval in the programmatic interface can be specified with
415: PETSC_MAX_REAL and -PETSC_MAX_REAL.
417: Level: intermediate
419: .seealso: EPSGetInterval(), EPSSetWhichEigenpairs()
420: @*/
421: PetscErrorCode EPSSetInterval(EPS eps,PetscReal inta,PetscReal intb)422: {
427: if (inta>=intb) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be inta<intb");
428: eps->inta = inta;
429: eps->intb = intb;
430: eps->state = EPS_STATE_INITIAL;
431: return(0);
432: }
436: /*@
437: EPSGetInterval - Gets the computational interval for spectrum slicing.
439: Not Collective
441: Input Parameter:
442: . eps - eigensolver context
444: Output Parameters:
445: + inta - left end of the interval
446: - intb - right end of the interval
448: Level: intermediate
450: Note:
451: If the interval was not set by the user, then zeros are returned.
453: .seealso: EPSSetInterval()
454: @*/
455: PetscErrorCode EPSGetInterval(EPS eps,PetscReal* inta,PetscReal* intb)456: {
461: if (inta) *inta = eps->inta;
462: if (intb) *intb = eps->intb;
463: return(0);
464: }
468: /*@
469: EPSSetST - Associates a spectral transformation object to the eigensolver.
471: Collective on EPS473: Input Parameters:
474: + eps - eigensolver context obtained from EPSCreate()
475: - st - the spectral transformation object
477: Note:
478: Use EPSGetST() to retrieve the spectral transformation context (for example,
479: to free it at the end of the computations).
481: Level: advanced
483: .seealso: EPSGetST()
484: @*/
485: PetscErrorCode EPSSetST(EPS eps,ST st)486: {
493: PetscObjectReference((PetscObject)st);
494: STDestroy(&eps->st);
495: eps->st = st;
496: PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->st);
497: return(0);
498: }
502: /*@
503: EPSGetST - Obtain the spectral transformation (ST) object associated
504: to the eigensolver object.
506: Not Collective
508: Input Parameters:
509: . eps - eigensolver context obtained from EPSCreate()
511: Output Parameter:
512: . st - spectral transformation context
514: Level: intermediate
516: .seealso: EPSSetST()
517: @*/
518: PetscErrorCode EPSGetST(EPS eps,ST *st)519: {
525: if (!eps->st) {
526: STCreate(PetscObjectComm((PetscObject)eps),&eps->st);
527: PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->st);
528: }
529: *st = eps->st;
530: return(0);
531: }
535: /*@
536: EPSSetBV - Associates a basis vectors object to the eigensolver.
538: Collective on EPS540: Input Parameters:
541: + eps - eigensolver context obtained from EPSCreate()
542: - V - the basis vectors object
544: Note:
545: Use EPSGetBV() to retrieve the basis vectors context (for example,
546: to free them at the end of the computations).
548: Level: advanced
550: .seealso: EPSGetBV()
551: @*/
552: PetscErrorCode EPSSetBV(EPS eps,BV V)553: {
560: PetscObjectReference((PetscObject)V);
561: BVDestroy(&eps->V);
562: eps->V = V;
563: PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->V);
564: return(0);
565: }
569: /*@
570: EPSGetBV - Obtain the basis vectors object associated to the eigensolver object.
572: Not Collective
574: Input Parameters:
575: . eps - eigensolver context obtained from EPSCreate()
577: Output Parameter:
578: . V - basis vectors context
580: Level: advanced
582: .seealso: EPSSetBV()
583: @*/
584: PetscErrorCode EPSGetBV(EPS eps,BV *V)585: {
591: if (!eps->V) {
592: BVCreate(PetscObjectComm((PetscObject)eps),&eps->V);
593: PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->V);
594: }
595: *V = eps->V;
596: return(0);
597: }
601: /*@
602: EPSSetRG - Associates a region object to the eigensolver.
604: Collective on EPS606: Input Parameters:
607: + eps - eigensolver context obtained from EPSCreate()
608: - rg - the region object
610: Note:
611: Use EPSGetRG() to retrieve the region context (for example,
612: to free it at the end of the computations).
614: Level: advanced
616: .seealso: EPSGetRG()
617: @*/
618: PetscErrorCode EPSSetRG(EPS eps,RG rg)619: {
626: PetscObjectReference((PetscObject)rg);
627: RGDestroy(&eps->rg);
628: eps->rg = rg;
629: PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->rg);
630: return(0);
631: }
635: /*@
636: EPSGetRG - Obtain the region object associated to the eigensolver.
638: Not Collective
640: Input Parameters:
641: . eps - eigensolver context obtained from EPSCreate()
643: Output Parameter:
644: . rg - region context
646: Level: advanced
648: .seealso: EPSSetRG()
649: @*/
650: PetscErrorCode EPSGetRG(EPS eps,RG *rg)651: {
657: if (!eps->rg) {
658: RGCreate(PetscObjectComm((PetscObject)eps),&eps->rg);
659: PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->rg);
660: }
661: *rg = eps->rg;
662: return(0);
663: }
667: /*@
668: EPSSetDS - Associates a direct solver object to the eigensolver.
670: Collective on EPS672: Input Parameters:
673: + eps - eigensolver context obtained from EPSCreate()
674: - ds - the direct solver object
676: Note:
677: Use EPSGetDS() to retrieve the direct solver context (for example,
678: to free it at the end of the computations).
680: Level: advanced
682: .seealso: EPSGetDS()
683: @*/
684: PetscErrorCode EPSSetDS(EPS eps,DS ds)685: {
692: PetscObjectReference((PetscObject)ds);
693: DSDestroy(&eps->ds);
694: eps->ds = ds;
695: PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->ds);
696: return(0);
697: }
701: /*@
702: EPSGetDS - Obtain the direct solver object associated to the eigensolver object.
704: Not Collective
706: Input Parameters:
707: . eps - eigensolver context obtained from EPSCreate()
709: Output Parameter:
710: . ds - direct solver context
712: Level: advanced
714: .seealso: EPSSetDS()
715: @*/
716: PetscErrorCode EPSGetDS(EPS eps,DS *ds)717: {
723: if (!eps->ds) {
724: DSCreate(PetscObjectComm((PetscObject)eps),&eps->ds);
725: PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->ds);
726: }
727: *ds = eps->ds;
728: return(0);
729: }
733: /*@
734: EPSIsGeneralized - Ask if the EPS object corresponds to a generalized
735: eigenvalue problem.
737: Not collective
739: Input Parameter:
740: . eps - the eigenproblem solver context
742: Output Parameter:
743: . is - the answer
745: Level: intermediate
747: .seealso: EPSIsHermitian(), EPSIsPositive()
748: @*/
749: PetscErrorCode EPSIsGeneralized(EPS eps,PetscBool* is)750: {
754: *is = eps->isgeneralized;
755: return(0);
756: }
760: /*@
761: EPSIsHermitian - Ask if the EPS object corresponds to a Hermitian
762: eigenvalue problem.
764: Not collective
766: Input Parameter:
767: . eps - the eigenproblem solver context
769: Output Parameter:
770: . is - the answer
772: Level: intermediate
774: .seealso: EPSIsGeneralized(), EPSIsPositive()
775: @*/
776: PetscErrorCode EPSIsHermitian(EPS eps,PetscBool* is)777: {
781: *is = eps->ishermitian;
782: return(0);
783: }
787: /*@
788: EPSIsPositive - Ask if the EPS object corresponds to an eigenvalue
789: problem type that requires a positive (semi-) definite matrix B.
791: Not collective
793: Input Parameter:
794: . eps - the eigenproblem solver context
796: Output Parameter:
797: . is - the answer
799: Level: intermediate
801: .seealso: EPSIsGeneralized(), EPSIsHermitian()
802: @*/
803: PetscErrorCode EPSIsPositive(EPS eps,PetscBool* is)804: {
808: *is = eps->ispositive;
809: return(0);
810: }