Actual source code: pepbasic.c
slepc-main 2024-11-22
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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>
16: /* Logging support */
17: PetscClassId PEP_CLASSID = 0;
18: PetscLogEvent PEP_SetUp = 0,PEP_Solve = 0,PEP_Refine = 0,PEP_CISS_SVD = 0;
20: /* List of registered PEP routines */
21: PetscFunctionList PEPList = NULL;
22: PetscBool PEPRegisterAllCalled = PETSC_FALSE;
24: /* List of registered PEP monitors */
25: PetscFunctionList PEPMonitorList = NULL;
26: PetscFunctionList PEPMonitorCreateList = NULL;
27: PetscFunctionList PEPMonitorDestroyList = NULL;
28: PetscBool PEPMonitorRegisterAllCalled = PETSC_FALSE;
30: /*@
31: PEPCreate - Creates the default PEP context.
33: Collective
35: Input Parameter:
36: . comm - MPI communicator
38: Output Parameter:
39: . outpep - location to put the PEP context
41: Note:
42: The default PEP type is PEPTOAR
44: Level: beginner
46: .seealso: PEPSetUp(), PEPSolve(), PEPDestroy(), PEP
47: @*/
48: PetscErrorCode PEPCreate(MPI_Comm comm,PEP *outpep)
49: {
50: PEP pep;
52: PetscFunctionBegin;
53: PetscAssertPointer(outpep,2);
54: PetscCall(PEPInitializePackage());
55: PetscCall(SlepcHeaderCreate(pep,PEP_CLASSID,"PEP","Polynomial Eigenvalue Problem","PEP",comm,PEPDestroy,PEPView));
57: pep->max_it = PETSC_DETERMINE;
58: pep->nev = 1;
59: pep->ncv = PETSC_DETERMINE;
60: pep->mpd = PETSC_DETERMINE;
61: pep->nini = 0;
62: pep->target = 0.0;
63: pep->tol = PETSC_DETERMINE;
64: pep->conv = PEP_CONV_REL;
65: pep->stop = PEP_STOP_BASIC;
66: pep->which = (PEPWhich)0;
67: pep->basis = PEP_BASIS_MONOMIAL;
68: pep->problem_type = (PEPProblemType)0;
69: pep->scale = PEP_SCALE_NONE;
70: pep->sfactor = 1.0;
71: pep->dsfactor = 1.0;
72: pep->sits = 5;
73: pep->slambda = 1.0;
74: pep->refine = PEP_REFINE_NONE;
75: pep->npart = 1;
76: pep->rtol = PETSC_DETERMINE;
77: pep->rits = PETSC_DETERMINE;
78: pep->scheme = (PEPRefineScheme)0;
79: pep->extract = (PEPExtract)0;
80: pep->trackall = PETSC_FALSE;
82: pep->converged = PEPConvergedRelative;
83: pep->convergeduser = NULL;
84: pep->convergeddestroy= NULL;
85: pep->stopping = PEPStoppingBasic;
86: pep->stoppinguser = NULL;
87: pep->stoppingdestroy = NULL;
88: pep->convergedctx = NULL;
89: pep->stoppingctx = NULL;
90: pep->numbermonitors = 0;
92: pep->st = NULL;
93: pep->ds = NULL;
94: pep->V = NULL;
95: pep->rg = 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: PetscCall(PetscNew(&pep->sc));
124: *outpep = pep;
125: PetscFunctionReturn(PETSC_SUCCESS);
126: }
128: /*@
129: PEPSetType - Selects the particular solver to be used in the PEP object.
131: Logically Collective
133: Input Parameters:
134: + pep - the polynomial eigensolver context
135: - type - a known method
137: Options Database Key:
138: . -pep_type <method> - Sets the method; use -help for a list
139: of available methods
141: Notes:
142: See "slepc/include/slepcpep.h" for available methods. The default
143: is PEPTOAR.
145: Normally, it is best to use the PEPSetFromOptions() command and
146: then set the PEP type from the options database rather than by using
147: this routine. Using the options database provides the user with
148: maximum flexibility in evaluating the different available methods.
149: The PEPSetType() routine is provided for those situations where it
150: is necessary to set the iterative solver independently of the command
151: line or options database.
153: Level: intermediate
155: .seealso: PEPType
156: @*/
157: PetscErrorCode PEPSetType(PEP pep,PEPType type)
158: {
159: PetscErrorCode (*r)(PEP);
160: PetscBool match;
162: PetscFunctionBegin;
164: PetscAssertPointer(type,2);
166: PetscCall(PetscObjectTypeCompare((PetscObject)pep,type,&match));
167: if (match) PetscFunctionReturn(PETSC_SUCCESS);
169: PetscCall(PetscFunctionListFind(PEPList,type,&r));
170: PetscCheck(r,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown PEP type given: %s",type);
172: PetscTryTypeMethod(pep,destroy);
173: PetscCall(PetscMemzero(pep->ops,sizeof(struct _PEPOps)));
175: pep->state = PEP_STATE_INITIAL;
176: PetscCall(PetscObjectChangeTypeName((PetscObject)pep,type));
177: PetscCall((*r)(pep));
178: PetscFunctionReturn(PETSC_SUCCESS);
179: }
181: /*@
182: PEPGetType - Gets the PEP type as a string from the PEP object.
184: Not Collective
186: Input Parameter:
187: . pep - the eigensolver context
189: Output Parameter:
190: . type - name of PEP method
192: Level: intermediate
194: .seealso: PEPSetType()
195: @*/
196: PetscErrorCode PEPGetType(PEP pep,PEPType *type)
197: {
198: PetscFunctionBegin;
200: PetscAssertPointer(type,2);
201: *type = ((PetscObject)pep)->type_name;
202: PetscFunctionReturn(PETSC_SUCCESS);
203: }
205: /*@C
206: PEPRegister - Adds a method to the polynomial eigenproblem solver package.
208: Not Collective
210: Input Parameters:
211: + name - name of a new user-defined solver
212: - function - routine to create the solver context
214: Notes:
215: PEPRegister() may be called multiple times to add several user-defined solvers.
217: Example Usage:
218: .vb
219: PEPRegister("my_solver",MySolverCreate);
220: .ve
222: Then, your solver can be chosen with the procedural interface via
223: $ PEPSetType(pep,"my_solver")
224: or at runtime via the option
225: $ -pep_type my_solver
227: Level: advanced
229: .seealso: PEPRegisterAll()
230: @*/
231: PetscErrorCode PEPRegister(const char *name,PetscErrorCode (*function)(PEP))
232: {
233: PetscFunctionBegin;
234: PetscCall(PEPInitializePackage());
235: PetscCall(PetscFunctionListAdd(&PEPList,name,function));
236: PetscFunctionReturn(PETSC_SUCCESS);
237: }
239: /*@C
240: PEPMonitorRegister - Adds PEP monitor routine.
242: Not Collective
244: Input Parameters:
245: + name - name of a new monitor routine
246: . vtype - a PetscViewerType for the output
247: . format - a PetscViewerFormat for the output
248: . monitor - monitor routine
249: . create - creation routine, or NULL
250: - destroy - destruction routine, or NULL
252: Notes:
253: PEPMonitorRegister() may be called multiple times to add several user-defined monitors.
255: Example Usage:
256: .vb
257: PEPMonitorRegister("my_monitor",PETSCVIEWERASCII,PETSC_VIEWER_ASCII_INFO_DETAIL,MyMonitor,NULL,NULL);
258: .ve
260: Then, your monitor can be chosen with the procedural interface via
261: $ PEPMonitorSetFromOptions(pep,"-pep_monitor_my_monitor","my_monitor",NULL)
262: or at runtime via the option
263: $ -pep_monitor_my_monitor
265: Level: advanced
267: .seealso: PEPMonitorRegisterAll()
268: @*/
269: PetscErrorCode PEPMonitorRegister(const char name[],PetscViewerType vtype,PetscViewerFormat format,PetscErrorCode (*monitor)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode (*create)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode (*destroy)(PetscViewerAndFormat**))
270: {
271: char key[PETSC_MAX_PATH_LEN];
273: PetscFunctionBegin;
274: PetscCall(PEPInitializePackage());
275: PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
276: PetscCall(PetscFunctionListAdd(&PEPMonitorList,key,monitor));
277: if (create) PetscCall(PetscFunctionListAdd(&PEPMonitorCreateList,key,create));
278: if (destroy) PetscCall(PetscFunctionListAdd(&PEPMonitorDestroyList,key,destroy));
279: PetscFunctionReturn(PETSC_SUCCESS);
280: }
282: /*@
283: PEPReset - Resets the PEP context to the initial state (prior to setup)
284: and destroys any allocated Vecs and Mats.
286: Collective
288: Input Parameter:
289: . pep - eigensolver context obtained from PEPCreate()
291: Level: advanced
293: .seealso: PEPDestroy()
294: @*/
295: PetscErrorCode PEPReset(PEP pep)
296: {
297: PetscFunctionBegin;
299: if (!pep) PetscFunctionReturn(PETSC_SUCCESS);
300: PetscTryTypeMethod(pep,reset);
301: if (pep->st) PetscCall(STReset(pep->st));
302: if (pep->refineksp) PetscCall(KSPReset(pep->refineksp));
303: if (pep->nmat) {
304: PetscCall(MatDestroyMatrices(pep->nmat,&pep->A));
305: PetscCall(PetscFree2(pep->pbc,pep->nrma));
306: PetscCall(PetscFree(pep->solvematcoeffs));
307: pep->nmat = 0;
308: }
309: PetscCall(VecDestroy(&pep->Dl));
310: PetscCall(VecDestroy(&pep->Dr));
311: PetscCall(BVDestroy(&pep->V));
312: PetscCall(VecDestroyVecs(pep->nwork,&pep->work));
313: pep->nwork = 0;
314: pep->state = PEP_STATE_INITIAL;
315: PetscFunctionReturn(PETSC_SUCCESS);
316: }
318: /*@
319: PEPDestroy - Destroys the PEP context.
321: Collective
323: Input Parameter:
324: . pep - eigensolver context obtained from PEPCreate()
326: Level: beginner
328: .seealso: PEPCreate(), PEPSetUp(), PEPSolve()
329: @*/
330: PetscErrorCode PEPDestroy(PEP *pep)
331: {
332: PetscFunctionBegin;
333: if (!*pep) PetscFunctionReturn(PETSC_SUCCESS);
335: if (--((PetscObject)*pep)->refct > 0) { *pep = NULL; PetscFunctionReturn(PETSC_SUCCESS); }
336: PetscCall(PEPReset(*pep));
337: PetscTryTypeMethod(*pep,destroy);
338: if ((*pep)->eigr) PetscCall(PetscFree4((*pep)->eigr,(*pep)->eigi,(*pep)->errest,(*pep)->perm));
339: PetscCall(STDestroy(&(*pep)->st));
340: PetscCall(RGDestroy(&(*pep)->rg));
341: PetscCall(DSDestroy(&(*pep)->ds));
342: PetscCall(KSPDestroy(&(*pep)->refineksp));
343: PetscCall(PetscSubcommDestroy(&(*pep)->refinesubc));
344: PetscCall(PetscFree((*pep)->sc));
345: /* just in case the initial vectors have not been used */
346: PetscCall(SlepcBasisDestroy_Private(&(*pep)->nini,&(*pep)->IS));
347: if ((*pep)->convergeddestroy) PetscCall((*(*pep)->convergeddestroy)(&(*pep)->convergedctx));
348: if ((*pep)->stoppingdestroy) PetscCall((*(*pep)->stoppingdestroy)(&(*pep)->stoppingctx));
349: PetscCall(PEPMonitorCancel(*pep));
350: PetscCall(PetscHeaderDestroy(pep));
351: PetscFunctionReturn(PETSC_SUCCESS);
352: }
354: /*@
355: PEPSetBV - Associates a basis vectors object to the polynomial eigensolver.
357: Collective
359: Input Parameters:
360: + pep - eigensolver context obtained from PEPCreate()
361: - bv - the basis vectors object
363: Note:
364: Use PEPGetBV() to retrieve the basis vectors context (for example,
365: to free it at the end of the computations).
367: Level: advanced
369: .seealso: PEPGetBV()
370: @*/
371: PetscErrorCode PEPSetBV(PEP pep,BV bv)
372: {
373: PetscFunctionBegin;
376: PetscCheckSameComm(pep,1,bv,2);
377: PetscCall(PetscObjectReference((PetscObject)bv));
378: PetscCall(BVDestroy(&pep->V));
379: pep->V = bv;
380: PetscFunctionReturn(PETSC_SUCCESS);
381: }
383: /*@
384: PEPGetBV - Obtain the basis vectors object associated to the polynomial
385: eigensolver object.
387: Not Collective
389: Input Parameters:
390: . pep - eigensolver context obtained from PEPCreate()
392: Output Parameter:
393: . bv - basis vectors context
395: Level: advanced
397: .seealso: PEPSetBV()
398: @*/
399: PetscErrorCode PEPGetBV(PEP pep,BV *bv)
400: {
401: PetscFunctionBegin;
403: PetscAssertPointer(bv,2);
404: if (!pep->V) {
405: PetscCall(BVCreate(PetscObjectComm((PetscObject)pep),&pep->V));
406: PetscCall(PetscObjectIncrementTabLevel((PetscObject)pep->V,(PetscObject)pep,0));
407: PetscCall(PetscObjectSetOptions((PetscObject)pep->V,((PetscObject)pep)->options));
408: }
409: *bv = pep->V;
410: PetscFunctionReturn(PETSC_SUCCESS);
411: }
413: /*@
414: PEPSetRG - Associates a region object to the polynomial eigensolver.
416: Collective
418: Input Parameters:
419: + pep - eigensolver context obtained from PEPCreate()
420: - rg - the region object
422: Note:
423: Use PEPGetRG() to retrieve the region context (for example,
424: to free it at the end of the computations).
426: Level: advanced
428: .seealso: PEPGetRG()
429: @*/
430: PetscErrorCode PEPSetRG(PEP pep,RG rg)
431: {
432: PetscFunctionBegin;
434: if (rg) {
436: PetscCheckSameComm(pep,1,rg,2);
437: }
438: PetscCall(PetscObjectReference((PetscObject)rg));
439: PetscCall(RGDestroy(&pep->rg));
440: pep->rg = rg;
441: PetscFunctionReturn(PETSC_SUCCESS);
442: }
444: /*@
445: PEPGetRG - Obtain the region object associated to the
446: polynomial eigensolver object.
448: Not Collective
450: Input Parameters:
451: . pep - eigensolver context obtained from PEPCreate()
453: Output Parameter:
454: . rg - region context
456: Level: advanced
458: .seealso: PEPSetRG()
459: @*/
460: PetscErrorCode PEPGetRG(PEP pep,RG *rg)
461: {
462: PetscFunctionBegin;
464: PetscAssertPointer(rg,2);
465: if (!pep->rg) {
466: PetscCall(RGCreate(PetscObjectComm((PetscObject)pep),&pep->rg));
467: PetscCall(PetscObjectIncrementTabLevel((PetscObject)pep->rg,(PetscObject)pep,0));
468: PetscCall(PetscObjectSetOptions((PetscObject)pep->rg,((PetscObject)pep)->options));
469: }
470: *rg = pep->rg;
471: PetscFunctionReturn(PETSC_SUCCESS);
472: }
474: /*@
475: PEPSetDS - Associates a direct solver object to the polynomial eigensolver.
477: Collective
479: Input Parameters:
480: + pep - eigensolver context obtained from PEPCreate()
481: - ds - the direct solver object
483: Note:
484: Use PEPGetDS() to retrieve the direct solver context (for example,
485: to free it at the end of the computations).
487: Level: advanced
489: .seealso: PEPGetDS()
490: @*/
491: PetscErrorCode PEPSetDS(PEP pep,DS ds)
492: {
493: PetscFunctionBegin;
496: PetscCheckSameComm(pep,1,ds,2);
497: PetscCall(PetscObjectReference((PetscObject)ds));
498: PetscCall(DSDestroy(&pep->ds));
499: pep->ds = ds;
500: PetscFunctionReturn(PETSC_SUCCESS);
501: }
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: {
521: PetscFunctionBegin;
523: PetscAssertPointer(ds,2);
524: if (!pep->ds) {
525: PetscCall(DSCreate(PetscObjectComm((PetscObject)pep),&pep->ds));
526: PetscCall(PetscObjectIncrementTabLevel((PetscObject)pep->ds,(PetscObject)pep,0));
527: PetscCall(PetscObjectSetOptions((PetscObject)pep->ds,((PetscObject)pep)->options));
528: }
529: *ds = pep->ds;
530: PetscFunctionReturn(PETSC_SUCCESS);
531: }
533: /*@
534: PEPSetST - Associates a spectral transformation object to the eigensolver.
536: Collective
538: Input Parameters:
539: + pep - eigensolver context obtained from PEPCreate()
540: - st - the spectral transformation object
542: Note:
543: Use PEPGetST() to retrieve the spectral transformation context (for example,
544: to free it at the end of the computations).
546: Level: advanced
548: .seealso: PEPGetST()
549: @*/
550: PetscErrorCode PEPSetST(PEP pep,ST st)
551: {
552: PetscFunctionBegin;
555: PetscCheckSameComm(pep,1,st,2);
556: PetscCall(PetscObjectReference((PetscObject)st));
557: PetscCall(STDestroy(&pep->st));
558: pep->st = st;
559: PetscFunctionReturn(PETSC_SUCCESS);
560: }
562: /*@
563: PEPGetST - Obtain the spectral transformation (ST) object associated
564: to the eigensolver object.
566: Not Collective
568: Input Parameters:
569: . pep - eigensolver context obtained from PEPCreate()
571: Output Parameter:
572: . st - spectral transformation context
574: Level: intermediate
576: .seealso: PEPSetST()
577: @*/
578: PetscErrorCode PEPGetST(PEP pep,ST *st)
579: {
580: PetscFunctionBegin;
582: PetscAssertPointer(st,2);
583: if (!pep->st) {
584: PetscCall(STCreate(PetscObjectComm((PetscObject)pep),&pep->st));
585: PetscCall(PetscObjectIncrementTabLevel((PetscObject)pep->st,(PetscObject)pep,0));
586: PetscCall(PetscObjectSetOptions((PetscObject)pep->st,((PetscObject)pep)->options));
587: }
588: *st = pep->st;
589: PetscFunctionReturn(PETSC_SUCCESS);
590: }
592: /*@
593: PEPRefineGetKSP - Obtain the ksp object used by the eigensolver
594: object in the refinement phase.
596: Collective
598: Input Parameters:
599: . pep - eigensolver context obtained from PEPCreate()
601: Output Parameter:
602: . ksp - ksp context
604: Level: advanced
606: .seealso: PEPSetRefine()
607: @*/
608: PetscErrorCode PEPRefineGetKSP(PEP pep,KSP *ksp)
609: {
610: MPI_Comm comm;
612: PetscFunctionBegin;
614: PetscAssertPointer(ksp,2);
615: if (!pep->refineksp) {
616: if (pep->npart>1) {
617: /* Split in subcomunicators */
618: PetscCall(PetscSubcommCreate(PetscObjectComm((PetscObject)pep),&pep->refinesubc));
619: PetscCall(PetscSubcommSetNumber(pep->refinesubc,pep->npart));
620: PetscCall(PetscSubcommSetType(pep->refinesubc,PETSC_SUBCOMM_CONTIGUOUS));
621: PetscCall(PetscSubcommGetChild(pep->refinesubc,&comm));
622: } else PetscCall(PetscObjectGetComm((PetscObject)pep,&comm));
623: PetscCall(KSPCreate(comm,&pep->refineksp));
624: PetscCall(PetscObjectIncrementTabLevel((PetscObject)pep->refineksp,(PetscObject)pep,0));
625: PetscCall(PetscObjectSetOptions((PetscObject)pep->refineksp,((PetscObject)pep)->options));
626: PetscCall(KSPSetOptionsPrefix(*ksp,((PetscObject)pep)->prefix));
627: PetscCall(KSPAppendOptionsPrefix(*ksp,"pep_refine_"));
628: PetscCall(KSPSetTolerances(pep->refineksp,SlepcDefaultTol(pep->rtol),PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
629: }
630: *ksp = pep->refineksp;
631: PetscFunctionReturn(PETSC_SUCCESS);
632: }
634: /*@
635: PEPSetTarget - Sets the value of the target.
637: Logically Collective
639: Input Parameters:
640: + pep - eigensolver context
641: - target - the value of the target
643: Options Database Key:
644: . -pep_target <scalar> - the value of the target
646: Notes:
647: The target is a scalar value used to determine the portion of the spectrum
648: of interest. It is used in combination with PEPSetWhichEigenpairs().
650: In the case of complex scalars, a complex value can be provided in the
651: command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
652: -pep_target 1.0+2.0i
654: Level: intermediate
656: .seealso: PEPGetTarget(), PEPSetWhichEigenpairs()
657: @*/
658: PetscErrorCode PEPSetTarget(PEP pep,PetscScalar target)
659: {
660: PetscFunctionBegin;
663: pep->target = target;
664: if (!pep->st) PetscCall(PEPGetST(pep,&pep->st));
665: PetscCall(STSetDefaultShift(pep->st,target));
666: PetscFunctionReturn(PETSC_SUCCESS);
667: }
669: /*@
670: PEPGetTarget - Gets the value of the target.
672: Not Collective
674: Input Parameter:
675: . pep - eigensolver context
677: Output Parameter:
678: . target - the value of the target
680: Note:
681: If the target was not set by the user, then zero is returned.
683: Level: intermediate
685: .seealso: PEPSetTarget()
686: @*/
687: PetscErrorCode PEPGetTarget(PEP pep,PetscScalar* target)
688: {
689: PetscFunctionBegin;
691: PetscAssertPointer(target,2);
692: *target = pep->target;
693: PetscFunctionReturn(PETSC_SUCCESS);
694: }
696: /*@
697: PEPSetInterval - Defines the computational interval for spectrum slicing.
699: Logically Collective
701: Input Parameters:
702: + pep - eigensolver context
703: . inta - left end of the interval
704: - intb - right end of the interval
706: Options Database Key:
707: . -pep_interval <a,b> - set [a,b] as the interval of interest
709: Notes:
710: Spectrum slicing is a technique employed for computing all eigenvalues of
711: symmetric eigenproblems in a given interval. This function provides the
712: interval to be considered. It must be used in combination with PEP_ALL, see
713: PEPSetWhichEigenpairs().
715: In the command-line option, two values must be provided. For an open interval,
716: one can give an infinite, e.g., -pep_interval 1.0,inf or -pep_interval -inf,1.0.
717: An open interval in the programmatic interface can be specified with
718: PETSC_MAX_REAL and -PETSC_MAX_REAL.
720: Level: intermediate
722: .seealso: PEPGetInterval(), PEPSetWhichEigenpairs()
723: @*/
724: PetscErrorCode PEPSetInterval(PEP pep,PetscReal inta,PetscReal intb)
725: {
726: PetscFunctionBegin;
730: PetscCheck(inta<intb,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be inta<intb");
731: if (pep->inta != inta || pep->intb != intb) {
732: pep->inta = inta;
733: pep->intb = intb;
734: pep->state = PEP_STATE_INITIAL;
735: }
736: PetscFunctionReturn(PETSC_SUCCESS);
737: }
739: /*@
740: PEPGetInterval - Gets the computational interval for spectrum slicing.
742: Not Collective
744: Input Parameter:
745: . pep - eigensolver context
747: Output Parameters:
748: + inta - left end of the interval
749: - intb - right end of the interval
751: Level: intermediate
753: Note:
754: If the interval was not set by the user, then zeros are returned.
756: .seealso: PEPSetInterval()
757: @*/
758: PetscErrorCode PEPGetInterval(PEP pep,PetscReal* inta,PetscReal* intb)
759: {
760: PetscFunctionBegin;
762: if (inta) *inta = pep->inta;
763: if (intb) *intb = pep->intb;
764: PetscFunctionReturn(PETSC_SUCCESS);
765: }