1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: PEP routines related to options that can be set via the command-line
12: or procedurally
13: */
15: #include <slepc/private/pepimpl.h> 16: #include <petscdraw.h>
18: /*@C
19: PEPMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
20: indicated by the user.
22: Collective on pep
24: Input Parameters:
25: + pep - the polynomial eigensolver context
26: . name - the monitor option name
27: . help - message indicating what monitoring is done
28: . manual - manual page for the monitor
29: . monitor - the monitor function, whose context is a PetscViewerAndFormat
30: - trackall - whether this monitor tracks all eigenvalues or not
32: Level: developer
34: .seealso: PEPMonitorSet(), PEPSetTrackAll(), PEPConvMonitorSetFromOptions()
35: @*/
36: PetscErrorCode PEPMonitorSetFromOptions(PEP pep,const char name[],const char help[],const char manual[],PetscErrorCode (*monitor)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscBool trackall) 37: {
38: PetscErrorCode ierr;
39: PetscBool flg;
40: PetscViewer viewer;
41: PetscViewerFormat format;
42: PetscViewerAndFormat *vf;
45: PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->options,((PetscObject)pep)->prefix,name,&viewer,&format,&flg);
46: if (flg) {
47: PetscViewerAndFormatCreate(viewer,format,&vf);
48: PetscObjectDereference((PetscObject)viewer);
49: PEPMonitorSet(pep,(PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
50: if (trackall) {
51: PEPSetTrackAll(pep,PETSC_TRUE);
52: }
53: }
54: return(0);
55: }
57: /*@C
58: PEPConvMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
59: indicated by the user (for monitors that only show iteration numbers of convergence).
61: Collective on pep
63: Input Parameters:
64: + pep - the polynomial eigensolver context
65: . name - the monitor option name
66: . help - message indicating what monitoring is done
67: . manual - manual page for the monitor
68: - monitor - the monitor function, whose context is a SlepcConvMonitor
70: Level: developer
72: .seealso: PEPMonitorSet(), PEPMonitorSetFromOptions()
73: @*/
74: PetscErrorCode PEPConvMonitorSetFromOptions(PEP pep,const char name[],const char help[],const char manual[],PetscErrorCode (*monitor)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,SlepcConvMonitor)) 75: {
76: PetscErrorCode ierr;
77: PetscBool flg;
78: PetscViewer viewer;
79: PetscViewerFormat format;
80: SlepcConvMonitor ctx;
83: PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->options,((PetscObject)pep)->prefix,name,&viewer,&format,&flg);
84: if (flg) {
85: SlepcConvMonitorCreate(viewer,format,&ctx);
86: PetscObjectDereference((PetscObject)viewer);
87: PEPMonitorSet(pep,(PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))monitor,ctx,(PetscErrorCode (*)(void**))SlepcConvMonitorDestroy);
88: }
89: return(0);
90: }
92: /*@
93: PEPSetFromOptions - Sets PEP options from the options database.
94: This routine must be called before PEPSetUp() if the user is to be
95: allowed to set the solver type.
97: Collective on pep
99: Input Parameters:
100: . pep - the polynomial eigensolver context
102: Notes:
103: To see all options, run your program with the -help option.
105: Level: beginner
106: @*/
107: PetscErrorCode PEPSetFromOptions(PEP pep)108: {
109: PetscErrorCode ierr;
110: char type[256];
111: PetscBool set,flg,flg1,flg2,flg3,flg4,flg5;
112: PetscReal r,t,array[2]={0,0};
113: PetscScalar s;
114: PetscInt i,j,k;
115: PetscDrawLG lg;
116: PEPScale scale;
117: PEPRefine refine;
118: PEPRefineScheme scheme;
122: PEPRegisterAll();
123: PetscObjectOptionsBegin((PetscObject)pep);
124: PetscOptionsFList("-pep_type","Polynomial eigensolver method","PEPSetType",PEPList,(char*)(((PetscObject)pep)->type_name?((PetscObject)pep)->type_name:PEPTOAR),type,sizeof(type),&flg);
125: if (flg) {
126: PEPSetType(pep,type);
127: } else if (!((PetscObject)pep)->type_name) {
128: PEPSetType(pep,PEPTOAR);
129: }
131: PetscOptionsBoolGroupBegin("-pep_general","General polynomial eigenvalue problem","PEPSetProblemType",&flg);
132: if (flg) { PEPSetProblemType(pep,PEP_GENERAL); }
133: PetscOptionsBoolGroup("-pep_hermitian","Hermitian polynomial eigenvalue problem","PEPSetProblemType",&flg);
134: if (flg) { PEPSetProblemType(pep,PEP_HERMITIAN); }
135: PetscOptionsBoolGroup("-pep_hyperbolic","Hyperbolic polynomial eigenvalue problem","PEPSetProblemType",&flg);
136: if (flg) { PEPSetProblemType(pep,PEP_HYPERBOLIC); }
137: PetscOptionsBoolGroupEnd("-pep_gyroscopic","Gyroscopic polynomial eigenvalue problem","PEPSetProblemType",&flg);
138: if (flg) { PEPSetProblemType(pep,PEP_GYROSCOPIC); }
140: scale = pep->scale;
141: PetscOptionsEnum("-pep_scale","Scaling strategy","PEPSetScale",PEPScaleTypes,(PetscEnum)scale,(PetscEnum*)&scale,&flg1);
142: r = pep->sfactor;
143: PetscOptionsReal("-pep_scale_factor","Scale factor","PEPSetScale",pep->sfactor,&r,&flg2);
144: if (!flg2 && r==1.0) r = PETSC_DEFAULT;
145: j = pep->sits;
146: PetscOptionsInt("-pep_scale_its","Number of iterations in diagonal scaling","PEPSetScale",pep->sits,&j,&flg3);
147: t = pep->slambda;
148: PetscOptionsReal("-pep_scale_lambda","Estimate of eigenvalue (modulus) for diagonal scaling","PEPSetScale",pep->slambda,&t,&flg4);
149: if (flg1 || flg2 || flg3 || flg4) { PEPSetScale(pep,scale,r,NULL,NULL,j,t); }
151: PetscOptionsEnum("-pep_extract","Extraction method","PEPSetExtract",PEPExtractTypes,(PetscEnum)pep->extract,(PetscEnum*)&pep->extract,NULL);
153: refine = pep->refine;
154: PetscOptionsEnum("-pep_refine","Iterative refinement method","PEPSetRefine",PEPRefineTypes,(PetscEnum)refine,(PetscEnum*)&refine,&flg1);
155: i = pep->npart;
156: PetscOptionsInt("-pep_refine_partitions","Number of partitions of the communicator for iterative refinement","PEPSetRefine",pep->npart,&i,&flg2);
157: r = pep->rtol;
158: PetscOptionsReal("-pep_refine_tol","Tolerance for iterative refinement","PEPSetRefine",pep->rtol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/1000:pep->rtol,&r,&flg3);
159: j = pep->rits;
160: PetscOptionsInt("-pep_refine_its","Maximum number of iterations for iterative refinement","PEPSetRefine",pep->rits,&j,&flg4);
161: scheme = pep->scheme;
162: PetscOptionsEnum("-pep_refine_scheme","Scheme used for linear systems within iterative refinement","PEPSetRefine",PEPRefineSchemes,(PetscEnum)scheme,(PetscEnum*)&scheme,&flg5);
163: if (flg1 || flg2 || flg3 || flg4 || flg5) { PEPSetRefine(pep,refine,i,r,j,scheme); }
165: i = pep->max_it;
166: PetscOptionsInt("-pep_max_it","Maximum number of iterations","PEPSetTolerances",pep->max_it,&i,&flg1);
167: r = pep->tol;
168: PetscOptionsReal("-pep_tol","Tolerance","PEPSetTolerances",pep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:pep->tol,&r,&flg2);
169: if (flg1 || flg2) { PEPSetTolerances(pep,r,i); }
171: PetscOptionsBoolGroupBegin("-pep_conv_rel","Relative error convergence test","PEPSetConvergenceTest",&flg);
172: if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_REL); }
173: PetscOptionsBoolGroup("-pep_conv_norm","Convergence test relative to the matrix norms","PEPSetConvergenceTest",&flg);
174: if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_NORM); }
175: PetscOptionsBoolGroup("-pep_conv_abs","Absolute error convergence test","PEPSetConvergenceTest",&flg);
176: if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_ABS); }
177: PetscOptionsBoolGroupEnd("-pep_conv_user","User-defined convergence test","PEPSetConvergenceTest",&flg);
178: if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_USER); }
180: PetscOptionsBoolGroupBegin("-pep_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","PEPSetStoppingTest",&flg);
181: if (flg) { PEPSetStoppingTest(pep,PEP_STOP_BASIC); }
182: PetscOptionsBoolGroupEnd("-pep_stop_user","User-defined stopping test","PEPSetStoppingTest",&flg);
183: if (flg) { PEPSetStoppingTest(pep,PEP_STOP_USER); }
185: i = pep->nev;
186: PetscOptionsInt("-pep_nev","Number of eigenvalues to compute","PEPSetDimensions",pep->nev,&i,&flg1);
187: j = pep->ncv;
188: PetscOptionsInt("-pep_ncv","Number of basis vectors","PEPSetDimensions",pep->ncv,&j,&flg2);
189: k = pep->mpd;
190: PetscOptionsInt("-pep_mpd","Maximum dimension of projected problem","PEPSetDimensions",pep->mpd,&k,&flg3);
191: if (flg1 || flg2 || flg3) { PEPSetDimensions(pep,i,j,k); }
193: PetscOptionsEnum("-pep_basis","Polynomial basis","PEPSetBasis",PEPBasisTypes,(PetscEnum)pep->basis,(PetscEnum*)&pep->basis,NULL);
195: PetscOptionsBoolGroupBegin("-pep_largest_magnitude","Compute largest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg);
196: if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_MAGNITUDE); }
197: PetscOptionsBoolGroup("-pep_smallest_magnitude","Compute smallest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg);
198: if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_MAGNITUDE); }
199: PetscOptionsBoolGroup("-pep_largest_real","Compute eigenvalues with largest real parts","PEPSetWhichEigenpairs",&flg);
200: if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_REAL); }
201: PetscOptionsBoolGroup("-pep_smallest_real","Compute eigenvalues with smallest real parts","PEPSetWhichEigenpairs",&flg);
202: if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_REAL); }
203: PetscOptionsBoolGroup("-pep_largest_imaginary","Compute eigenvalues with largest imaginary parts","PEPSetWhichEigenpairs",&flg);
204: if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_IMAGINARY); }
205: PetscOptionsBoolGroup("-pep_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","PEPSetWhichEigenpairs",&flg);
206: if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_IMAGINARY); }
207: PetscOptionsBoolGroup("-pep_target_magnitude","Compute eigenvalues closest to target","PEPSetWhichEigenpairs",&flg);
208: if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE); }
209: PetscOptionsBoolGroup("-pep_target_real","Compute eigenvalues with real parts closest to target","PEPSetWhichEigenpairs",&flg);
210: if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_REAL); }
211: PetscOptionsBoolGroup("-pep_target_imaginary","Compute eigenvalues with imaginary parts closest to target","PEPSetWhichEigenpairs",&flg);
212: if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_IMAGINARY); }
213: PetscOptionsBoolGroupEnd("-pep_all","Compute all eigenvalues in an interval or a region","PEPSetWhichEigenpairs",&flg);
214: if (flg) { PEPSetWhichEigenpairs(pep,PEP_ALL); }
216: PetscOptionsScalar("-pep_target","Value of the target","PEPSetTarget",pep->target,&s,&flg);
217: if (flg) {
218: if (pep->which!=PEP_TARGET_REAL && pep->which!=PEP_TARGET_IMAGINARY) {
219: PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
220: }
221: PEPSetTarget(pep,s);
222: }
224: k = 2;
225: PetscOptionsRealArray("-pep_interval","Computational interval (two real values separated with a comma without spaces)","PEPSetInterval",array,&k,&flg);
226: if (flg) {
227: if (k<2) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_SIZ,"Must pass two values in -pep_interval (comma-separated without spaces)");
228: PEPSetWhichEigenpairs(pep,PEP_ALL);
229: PEPSetInterval(pep,array[0],array[1]);
230: }
232: /* -----------------------------------------------------------------------*/
233: /*
234: Cancels all monitors hardwired into code before call to PEPSetFromOptions()
235: */
236: PetscOptionsBool("-pep_monitor_cancel","Remove any hardwired monitor routines","PEPMonitorCancel",PETSC_FALSE,&flg,&set);
237: if (set && flg) {
238: PEPMonitorCancel(pep);
239: }
240: /*
241: Text monitors
242: */
243: PEPMonitorSetFromOptions(pep,"-pep_monitor","Monitor first unconverged approximate eigenvalue and error estimate","PEPMonitorFirst",PEPMonitorFirst,PETSC_FALSE);
244: PEPConvMonitorSetFromOptions(pep,"-pep_monitor_conv","Monitor approximate eigenvalues and error estimates as they converge","PEPMonitorConverged",PEPMonitorConverged);
245: PEPMonitorSetFromOptions(pep,"-pep_monitor_all","Monitor approximate eigenvalues and error estimates","PEPMonitorAll",PEPMonitorAll,PETSC_TRUE);
246: /*
247: Line graph monitors
248: */
249: PetscOptionsBool("-pep_monitor_lg","Monitor first unconverged approximate error estimate graphically","PEPMonitorSet",PETSC_FALSE,&flg,&set);
250: if (set && flg) {
251: PEPMonitorLGCreate(PetscObjectComm((PetscObject)pep),NULL,"Error estimates",PETSC_DECIDE,PETSC_DECIDE,300,300,&lg);
252: PEPMonitorSet(pep,PEPMonitorLG,lg,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
253: }
254: PetscOptionsBool("-pep_monitor_lg_all","Monitor error estimates graphically","PEPMonitorSet",PETSC_FALSE,&flg,&set);
255: if (set && flg) {
256: PEPMonitorLGCreate(PetscObjectComm((PetscObject)pep),NULL,"Error estimates",PETSC_DECIDE,PETSC_DECIDE,300,300,&lg);
257: PEPMonitorSet(pep,PEPMonitorLGAll,lg,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
258: PEPSetTrackAll(pep,PETSC_TRUE);
259: }
261: /* -----------------------------------------------------------------------*/
262: PetscOptionsName("-pep_view","Print detailed information on solver used","PEPView",NULL);
263: PetscOptionsName("-pep_view_vectors","View computed eigenvectors","PEPVectorsView",NULL);
264: PetscOptionsName("-pep_view_values","View computed eigenvalues","PEPValuesView",NULL);
265: PetscOptionsName("-pep_converged_reason","Print reason for convergence, and number of iterations","PEPConvergedReasonView",NULL);
266: PetscOptionsName("-pep_error_absolute","Print absolute errors of each eigenpair","PEPErrorView",NULL);
267: PetscOptionsName("-pep_error_relative","Print relative errors of each eigenpair","PEPErrorView",NULL);
268: PetscOptionsName("-pep_error_backward","Print backward errors of each eigenpair","PEPErrorView",NULL);
270: if (pep->ops->setfromoptions) {
271: (*pep->ops->setfromoptions)(PetscOptionsObject,pep);
272: }
273: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)pep);
274: PetscOptionsEnd();
276: if (!pep->V) { PEPGetBV(pep,&pep->V); }
277: BVSetFromOptions(pep->V);
278: if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
279: RGSetFromOptions(pep->rg);
280: if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
281: DSSetFromOptions(pep->ds);
282: if (!pep->st) { PEPGetST(pep,&pep->st); }
283: PEPSetDefaultST(pep);
284: STSetFromOptions(pep->st);
285: if (!pep->refineksp) { PEPRefineGetKSP(pep,&pep->refineksp); }
286: KSPSetFromOptions(pep->refineksp);
287: return(0);
288: }
290: /*@C
291: PEPGetTolerances - Gets the tolerance and maximum iteration count used
292: by the PEP convergence tests.
294: Not Collective
296: Input Parameter:
297: . pep - the polynomial eigensolver context
299: Output Parameters:
300: + tol - the convergence tolerance
301: - maxits - maximum number of iterations
303: Notes:
304: The user can specify NULL for any parameter that is not needed.
306: Level: intermediate
308: .seealso: PEPSetTolerances()
309: @*/
310: PetscErrorCode PEPGetTolerances(PEP pep,PetscReal *tol,PetscInt *maxits)311: {
314: if (tol) *tol = pep->tol;
315: if (maxits) *maxits = pep->max_it;
316: return(0);
317: }
319: /*@
320: PEPSetTolerances - Sets the tolerance and maximum iteration count used
321: by the PEP convergence tests.
323: Logically Collective on pep
325: Input Parameters:
326: + pep - the polynomial eigensolver context
327: . tol - the convergence tolerance
328: - maxits - maximum number of iterations to use
330: Options Database Keys:
331: + -pep_tol <tol> - Sets the convergence tolerance
332: - -pep_max_it <maxits> - Sets the maximum number of iterations allowed
334: Notes:
335: Use PETSC_DEFAULT for either argument to assign a reasonably good value.
337: Level: intermediate
339: .seealso: PEPGetTolerances()
340: @*/
341: PetscErrorCode PEPSetTolerances(PEP pep,PetscReal tol,PetscInt maxits)342: {
347: if (tol == PETSC_DEFAULT) {
348: pep->tol = PETSC_DEFAULT;
349: pep->state = PEP_STATE_INITIAL;
350: } else {
351: if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
352: pep->tol = tol;
353: }
354: if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
355: pep->max_it = PETSC_DEFAULT;
356: pep->state = PEP_STATE_INITIAL;
357: } else {
358: if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
359: pep->max_it = maxits;
360: }
361: return(0);
362: }
364: /*@C
365: PEPGetDimensions - Gets the number of eigenvalues to compute
366: and the dimension of the subspace.
368: Not Collective
370: Input Parameter:
371: . pep - the polynomial eigensolver context
373: Output Parameters:
374: + nev - number of eigenvalues to compute
375: . ncv - the maximum dimension of the subspace to be used by the solver
376: - mpd - the maximum dimension allowed for the projected problem
378: Notes:
379: The user can specify NULL for any parameter that is not needed.
381: Level: intermediate
383: .seealso: PEPSetDimensions()
384: @*/
385: PetscErrorCode PEPGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)386: {
389: if (nev) *nev = pep->nev;
390: if (ncv) *ncv = pep->ncv;
391: if (mpd) *mpd = pep->mpd;
392: return(0);
393: }
395: /*@
396: PEPSetDimensions - Sets the number of eigenvalues to compute
397: and the dimension of the subspace.
399: Logically Collective on pep
401: Input Parameters:
402: + pep - the polynomial eigensolver context
403: . nev - number of eigenvalues to compute
404: . ncv - the maximum dimension of the subspace to be used by the solver
405: - mpd - the maximum dimension allowed for the projected problem
407: Options Database Keys:
408: + -pep_nev <nev> - Sets the number of eigenvalues
409: . -pep_ncv <ncv> - Sets the dimension of the subspace
410: - -pep_mpd <mpd> - Sets the maximum projected dimension
412: Notes:
413: Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
414: dependent on the solution method.
416: The parameters ncv and mpd are intimately related, so that the user is advised
417: to set one of them at most. Normal usage is that
418: (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
419: (b) in cases where nev is large, the user sets mpd.
421: The value of ncv should always be between nev and (nev+mpd), typically
422: ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
423: a smaller value should be used.
425: When computing all eigenvalues in an interval, see PEPSetInterval(), these
426: parameters lose relevance, and tuning must be done with PEPSTOARSetDimensions().
428: Level: intermediate
430: .seealso: PEPGetDimensions(), PEPSetInterval(), PEPSTOARSetDimensions()
431: @*/
432: PetscErrorCode PEPSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)433: {
439: if (nev<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
440: pep->nev = nev;
441: if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
442: pep->ncv = PETSC_DEFAULT;
443: } else {
444: if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
445: pep->ncv = ncv;
446: }
447: if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
448: pep->mpd = PETSC_DEFAULT;
449: } else {
450: if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
451: pep->mpd = mpd;
452: }
453: pep->state = PEP_STATE_INITIAL;
454: return(0);
455: }
457: /*@
458: PEPSetWhichEigenpairs - Specifies which portion of the spectrum is
459: to be sought.
461: Logically Collective on pep
463: Input Parameters:
464: + pep - eigensolver context obtained from PEPCreate()
465: - which - the portion of the spectrum to be sought
467: Possible values:
468: The parameter 'which' can have one of these values
470: + PEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
471: . PEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
472: . PEP_LARGEST_REAL - largest real parts
473: . PEP_SMALLEST_REAL - smallest real parts
474: . PEP_LARGEST_IMAGINARY - largest imaginary parts
475: . PEP_SMALLEST_IMAGINARY - smallest imaginary parts
476: . PEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
477: . PEP_TARGET_REAL - eigenvalues with real part closest to target
478: . PEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
479: . PEP_ALL - all eigenvalues contained in a given interval
480: - PEP_WHICH_USER - user defined ordering set with PEPSetEigenvalueComparison()
482: Options Database Keys:
483: + -pep_largest_magnitude - Sets largest eigenvalues in magnitude
484: . -pep_smallest_magnitude - Sets smallest eigenvalues in magnitude
485: . -pep_largest_real - Sets largest real parts
486: . -pep_smallest_real - Sets smallest real parts
487: . -pep_largest_imaginary - Sets largest imaginary parts
488: . -pep_smallest_imaginary - Sets smallest imaginary parts
489: . -pep_target_magnitude - Sets eigenvalues closest to target
490: . -pep_target_real - Sets real parts closest to target
491: . -pep_target_imaginary - Sets imaginary parts closest to target
492: - -pep_all - Sets all eigenvalues in an interval
494: Notes:
495: Not all eigensolvers implemented in PEP account for all the possible values
496: stated above. If SLEPc is compiled for real numbers PEP_LARGEST_IMAGINARY497: and PEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
498: for eigenvalue selection.
500: The target is a scalar value provided with PEPSetTarget().
502: The criterion PEP_TARGET_IMAGINARY is available only in case PETSc and
503: SLEPc have been built with complex scalars.
505: PEP_ALL is intended for use in combination with an interval (see
506: PEPSetInterval()), when all eigenvalues within the interval are requested.
507: In that case, the number of eigenvalues is unknown, so the nev parameter
508: has a different sense, see PEPSetDimensions().
510: Level: intermediate
512: .seealso: PEPGetWhichEigenpairs(), PEPSetTarget(), PEPSetInterval(),
513: PEPSetDimensions(), PEPSetEigenvalueComparison(), PEPWhich514: @*/
515: PetscErrorCode PEPSetWhichEigenpairs(PEP pep,PEPWhich which)516: {
520: switch (which) {
521: case PEP_LARGEST_MAGNITUDE:
522: case PEP_SMALLEST_MAGNITUDE:
523: case PEP_LARGEST_REAL:
524: case PEP_SMALLEST_REAL:
525: case PEP_LARGEST_IMAGINARY:
526: case PEP_SMALLEST_IMAGINARY:
527: case PEP_TARGET_MAGNITUDE:
528: case PEP_TARGET_REAL:
529: #if defined(PETSC_USE_COMPLEX)
530: case PEP_TARGET_IMAGINARY:
531: #endif
532: case PEP_ALL:
533: case PEP_WHICH_USER:
534: if (pep->which != which) {
535: pep->state = PEP_STATE_INITIAL;
536: pep->which = which;
537: }
538: break;
539: #if !defined(PETSC_USE_COMPLEX)
540: case PEP_TARGET_IMAGINARY:
541: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"PEP_TARGET_IMAGINARY can be used only with complex scalars");
542: #endif
543: default:544: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
545: }
546: return(0);
547: }
549: /*@
550: PEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
551: sought.
553: Not Collective
555: Input Parameter:
556: . pep - eigensolver context obtained from PEPCreate()
558: Output Parameter:
559: . which - the portion of the spectrum to be sought
561: Notes:
562: See PEPSetWhichEigenpairs() for possible values of 'which'.
564: Level: intermediate
566: .seealso: PEPSetWhichEigenpairs(), PEPWhich567: @*/
568: PetscErrorCode PEPGetWhichEigenpairs(PEP pep,PEPWhich *which)569: {
573: *which = pep->which;
574: return(0);
575: }
577: /*@C
578: PEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
579: when PEPSetWhichEigenpairs() is set to PEP_WHICH_USER.
581: Logically Collective on pep
583: Input Parameters:
584: + pep - eigensolver context obtained from PEPCreate()
585: . func - a pointer to the comparison function
586: - ctx - a context pointer (the last parameter to the comparison function)
588: Calling Sequence of func:
589: $ func(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
591: + ar - real part of the 1st eigenvalue
592: . ai - imaginary part of the 1st eigenvalue
593: . br - real part of the 2nd eigenvalue
594: . bi - imaginary part of the 2nd eigenvalue
595: . res - result of comparison
596: - ctx - optional context, as set by PEPSetEigenvalueComparison()
598: Note:
599: The returning parameter 'res' can be
600: + negative - if the 1st eigenvalue is preferred to the 2st one
601: . zero - if both eigenvalues are equally preferred
602: - positive - if the 2st eigenvalue is preferred to the 1st one
604: Level: advanced
606: .seealso: PEPSetWhichEigenpairs(), PEPWhich607: @*/
608: PetscErrorCode PEPSetEigenvalueComparison(PEP pep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)609: {
612: pep->sc->comparison = func;
613: pep->sc->comparisonctx = ctx;
614: pep->which = PEP_WHICH_USER;
615: return(0);
616: }
618: /*@
619: PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.
621: Logically Collective on pep
623: Input Parameters:
624: + pep - the polynomial eigensolver context
625: - type - a known type of polynomial eigenvalue problem
627: Options Database Keys:
628: + -pep_general - general problem with no particular structure
629: . -pep_hermitian - problem whose coefficient matrices are Hermitian
630: . -pep_hyperbolic - Hermitian problem that satisfies the definition of hyperbolic
631: - -pep_gyroscopic - problem with Hamiltonian structure
633: Notes:
634: Allowed values for the problem type are: general (PEP_GENERAL), Hermitian
635: (PEP_HERMITIAN), hyperbolic (PEP_HYPERBOLIC), and gyroscopic (PEP_GYROSCOPIC).
637: This function is used to instruct SLEPc to exploit certain structure in
638: the polynomial eigenproblem. By default, no particular structure is assumed.
640: If the problem matrices are Hermitian (symmetric in the real case) or
641: Hermitian/skew-Hermitian then the solver can exploit this fact to perform
642: less operations or provide better stability. Hyperbolic problems are a
643: particular case of Hermitian problems, some solvers may treat them simply as
644: Hermitian.
646: Level: intermediate
648: .seealso: PEPSetOperators(), PEPSetType(), PEPGetProblemType(), PEPProblemType649: @*/
650: PetscErrorCode PEPSetProblemType(PEP pep,PEPProblemType type)651: {
655: if (type!=PEP_GENERAL && type!=PEP_HERMITIAN && type!=PEP_HYPERBOLIC && type!=PEP_GYROSCOPIC) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
656: if (type != pep->problem_type) {
657: pep->problem_type = type;
658: pep->state = PEP_STATE_INITIAL;
659: }
660: return(0);
661: }
663: /*@
664: PEPGetProblemType - Gets the problem type from the PEP object.
666: Not Collective
668: Input Parameter:
669: . pep - the polynomial eigensolver context
671: Output Parameter:
672: . type - the problem type
674: Level: intermediate
676: .seealso: PEPSetProblemType(), PEPProblemType677: @*/
678: PetscErrorCode PEPGetProblemType(PEP pep,PEPProblemType *type)679: {
683: *type = pep->problem_type;
684: return(0);
685: }
687: /*@
688: PEPSetBasis - Specifies the type of polynomial basis used to describe the
689: polynomial eigenvalue problem.
691: Logically Collective on pep
693: Input Parameters:
694: + pep - the polynomial eigensolver context
695: - basis - the type of polynomial basis
697: Options Database Key:
698: . -pep_basis <basis> - Select the basis type
700: Notes:
701: By default, the coefficient matrices passed via PEPSetOperators() are
702: expressed in the monomial basis, i.e.
703: P(lambda) = A_0 + lambda*A_1 + lambda^2*A_2 + ... + lambda^d*A_d.
704: Other polynomial bases may have better numerical behaviour, but the user
705: must then pass the coefficient matrices accordingly.
707: Level: intermediate
709: .seealso: PEPSetOperators(), PEPGetBasis(), PEPBasis710: @*/
711: PetscErrorCode PEPSetBasis(PEP pep,PEPBasis basis)712: {
716: pep->basis = basis;
717: return(0);
718: }
720: /*@
721: PEPGetBasis - Gets the type of polynomial basis from the PEP object.
723: Not Collective
725: Input Parameter:
726: . pep - the polynomial eigensolver context
728: Output Parameter:
729: . basis - the polynomial basis
731: Level: intermediate
733: .seealso: PEPSetBasis(), PEPBasis734: @*/
735: PetscErrorCode PEPGetBasis(PEP pep,PEPBasis *basis)736: {
740: *basis = pep->basis;
741: return(0);
742: }
744: /*@
745: PEPSetTrackAll - Specifies if the solver must compute the residual of all
746: approximate eigenpairs or not.
748: Logically Collective on pep
750: Input Parameters:
751: + pep - the eigensolver context
752: - trackall - whether compute all residuals or not
754: Notes:
755: If the user sets trackall=PETSC_TRUE then the solver explicitly computes
756: the residual for each eigenpair approximation. Computing the residual is
757: usually an expensive operation and solvers commonly compute the associated
758: residual to the first unconverged eigenpair.
759: The options '-pep_monitor_all' and '-pep_monitor_lg_all' automatically
760: activate this option.
762: Level: developer
764: .seealso: PEPGetTrackAll()
765: @*/
766: PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)767: {
771: pep->trackall = trackall;
772: return(0);
773: }
775: /*@
776: PEPGetTrackAll - Returns the flag indicating whether all residual norms must
777: be computed or not.
779: Not Collective
781: Input Parameter:
782: . pep - the eigensolver context
784: Output Parameter:
785: . trackall - the returned flag
787: Level: developer
789: .seealso: PEPSetTrackAll()
790: @*/
791: PetscErrorCode PEPGetTrackAll(PEP pep,PetscBool *trackall)792: {
796: *trackall = pep->trackall;
797: return(0);
798: }
800: /*@C
801: PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
802: used in the convergence test.
804: Logically Collective on pep
806: Input Parameters:
807: + pep - eigensolver context obtained from PEPCreate()
808: . func - a pointer to the convergence test function
809: . ctx - context for private data for the convergence routine (may be null)
810: - destroy - a routine for destroying the context (may be null)
812: Calling Sequence of func:
813: $ func(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
815: + pep - eigensolver context obtained from PEPCreate()
816: . eigr - real part of the eigenvalue
817: . eigi - imaginary part of the eigenvalue
818: . res - residual norm associated to the eigenpair
819: . errest - (output) computed error estimate
820: - ctx - optional context, as set by PEPSetConvergenceTestFunction()
822: Note:
823: If the error estimate returned by the convergence test function is less than
824: the tolerance, then the eigenvalue is accepted as converged.
826: Level: advanced
828: .seealso: PEPSetConvergenceTest(), PEPSetTolerances()
829: @*/
830: PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))831: {
836: if (pep->convergeddestroy) {
837: (*pep->convergeddestroy)(pep->convergedctx);
838: }
839: pep->convergeduser = func;
840: pep->convergeddestroy = destroy;
841: pep->convergedctx = ctx;
842: if (func == PEPConvergedRelative) pep->conv = PEP_CONV_REL;
843: else if (func == PEPConvergedNorm) pep->conv = PEP_CONV_NORM;
844: else if (func == PEPConvergedAbsolute) pep->conv = PEP_CONV_ABS;
845: else {
846: pep->conv = PEP_CONV_USER;
847: pep->converged = pep->convergeduser;
848: }
849: return(0);
850: }
852: /*@
853: PEPSetConvergenceTest - Specifies how to compute the error estimate
854: used in the convergence test.
856: Logically Collective on pep
858: Input Parameters:
859: + pep - eigensolver context obtained from PEPCreate()
860: - conv - the type of convergence test
862: Options Database Keys:
863: + -pep_conv_abs - Sets the absolute convergence test
864: . -pep_conv_rel - Sets the convergence test relative to the eigenvalue
865: . -pep_conv_norm - Sets the convergence test relative to the matrix norms
866: - -pep_conv_user - Selects the user-defined convergence test
868: Note:
869: The parameter 'conv' can have one of these values
870: + PEP_CONV_ABS - absolute error ||r||
871: . PEP_CONV_REL - error relative to the eigenvalue l, ||r||/|l|
872: . PEP_CONV_NORM - error relative matrix norms, ||r||/sum_i(l^i*||A_i||)
873: - PEP_CONV_USER - function set by PEPSetConvergenceTestFunction()
875: Level: intermediate
877: .seealso: PEPGetConvergenceTest(), PEPSetConvergenceTestFunction(), PEPSetStoppingTest(), PEPConv878: @*/
879: PetscErrorCode PEPSetConvergenceTest(PEP pep,PEPConv conv)880: {
884: switch (conv) {
885: case PEP_CONV_ABS: pep->converged = PEPConvergedAbsolute; break;
886: case PEP_CONV_REL: pep->converged = PEPConvergedRelative; break;
887: case PEP_CONV_NORM: pep->converged = PEPConvergedNorm; break;
888: case PEP_CONV_USER:
889: if (!pep->convergeduser) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetConvergenceTestFunction() first");
890: pep->converged = pep->convergeduser;
891: break;
892: default:893: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
894: }
895: pep->conv = conv;
896: return(0);
897: }
899: /*@
900: PEPGetConvergenceTest - Gets the method used to compute the error estimate
901: used in the convergence test.
903: Not Collective
905: Input Parameters:
906: . pep - eigensolver context obtained from PEPCreate()
908: Output Parameters:
909: . conv - the type of convergence test
911: Level: intermediate
913: .seealso: PEPSetConvergenceTest(), PEPConv914: @*/
915: PetscErrorCode PEPGetConvergenceTest(PEP pep,PEPConv *conv)916: {
920: *conv = pep->conv;
921: return(0);
922: }
924: /*@C
925: PEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
926: iteration of the eigensolver.
928: Logically Collective on pep
930: Input Parameters:
931: + pep - eigensolver context obtained from PEPCreate()
932: . func - pointer to the stopping test function
933: . ctx - context for private data for the stopping routine (may be null)
934: - destroy - a routine for destroying the context (may be null)
936: Calling Sequence of func:
937: $ func(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ctx)
939: + pep - eigensolver context obtained from PEPCreate()
940: . its - current number of iterations
941: . max_it - maximum number of iterations
942: . nconv - number of currently converged eigenpairs
943: . nev - number of requested eigenpairs
944: . reason - (output) result of the stopping test
945: - ctx - optional context, as set by PEPSetStoppingTestFunction()
947: Note:
948: Normal usage is to first call the default routine PEPStoppingBasic() and then
949: set reason to PEP_CONVERGED_USER if some user-defined conditions have been
950: met. To let the eigensolver continue iterating, the result must be left as
951: PEP_CONVERGED_ITERATING.
953: Level: advanced
955: .seealso: PEPSetStoppingTest(), PEPStoppingBasic()
956: @*/
957: PetscErrorCode PEPSetStoppingTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))958: {
963: if (pep->stoppingdestroy) {
964: (*pep->stoppingdestroy)(pep->stoppingctx);
965: }
966: pep->stoppinguser = func;
967: pep->stoppingdestroy = destroy;
968: pep->stoppingctx = ctx;
969: if (func == PEPStoppingBasic) pep->stop = PEP_STOP_BASIC;
970: else {
971: pep->stop = PEP_STOP_USER;
972: pep->stopping = pep->stoppinguser;
973: }
974: return(0);
975: }
977: /*@
978: PEPSetStoppingTest - Specifies how to decide the termination of the outer
979: loop of the eigensolver.
981: Logically Collective on pep
983: Input Parameters:
984: + pep - eigensolver context obtained from PEPCreate()
985: - stop - the type of stopping test
987: Options Database Keys:
988: + -pep_stop_basic - Sets the default stopping test
989: - -pep_stop_user - Selects the user-defined stopping test
991: Note:
992: The parameter 'stop' can have one of these values
993: + PEP_STOP_BASIC - default stopping test
994: - PEP_STOP_USER - function set by PEPSetStoppingTestFunction()
996: Level: advanced
998: .seealso: PEPGetStoppingTest(), PEPSetStoppingTestFunction(), PEPSetConvergenceTest(), PEPStop999: @*/
1000: PetscErrorCode PEPSetStoppingTest(PEP pep,PEPStop stop)1001: {
1005: switch (stop) {
1006: case PEP_STOP_BASIC: pep->stopping = PEPStoppingBasic; break;
1007: case PEP_STOP_USER:
1008: if (!pep->stoppinguser) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetStoppingTestFunction() first");
1009: pep->stopping = pep->stoppinguser;
1010: break;
1011: default:1012: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
1013: }
1014: pep->stop = stop;
1015: return(0);
1016: }
1018: /*@
1019: PEPGetStoppingTest - Gets the method used to decide the termination of the outer
1020: loop of the eigensolver.
1022: Not Collective
1024: Input Parameters:
1025: . pep - eigensolver context obtained from PEPCreate()
1027: Output Parameters:
1028: . stop - the type of stopping test
1030: Level: advanced
1032: .seealso: PEPSetStoppingTest(), PEPStop1033: @*/
1034: PetscErrorCode PEPGetStoppingTest(PEP pep,PEPStop *stop)1035: {
1039: *stop = pep->stop;
1040: return(0);
1041: }
1043: /*@
1044: PEPSetScale - Specifies the scaling strategy to be used.
1046: Logically Collective on pep
1048: Input Parameters:
1049: + pep - the eigensolver context
1050: . scale - scaling strategy
1051: . alpha - the scaling factor used in the scalar strategy
1052: . Dl - the left diagonal matrix of the diagonal scaling algorithm
1053: . Dr - the right diagonal matrix of the diagonal scaling algorithm
1054: . its - number of iterations of the diagonal scaling algorithm
1055: - lambda - approximation to wanted eigenvalues (modulus)
1057: Options Database Keys:
1058: + -pep_scale <type> - scaling type, one of <none,scalar,diagonal,both>
1059: . -pep_scale_factor <alpha> - the scaling factor
1060: . -pep_scale_its <its> - number of iterations
1061: - -pep_scale_lambda <lambda> - approximation to eigenvalues
1063: Notes:
1064: There are two non-exclusive scaling strategies: scalar and diagonal.
1066: In the scalar strategy, scaling is applied to the eigenvalue, that is,
1067: mu = lambda/alpha is the new eigenvalue and all matrices are scaled
1068: accordingly. After solving the scaled problem, the original lambda is
1069: recovered. Parameter 'alpha' must be positive. Use PETSC_DEFAULT to let
1070: the solver compute a reasonable scaling factor.
1072: In the diagonal strategy, the solver works implicitly with matrix Dl*A*Dr,
1073: where Dl and Dr are appropriate diagonal matrices. This improves the accuracy
1074: of the computed results in some cases. The user may provide the Dr and Dl
1075: matrices represented as Vec objects storing diagonal elements. If not
1076: provided, these matrices are computed internally. This option requires
1077: that the polynomial coefficient matrices are of MATAIJ type.
1078: The parameter 'its' is the number of iterations performed by the method.
1079: Parameter 'lambda' must be positive. Use PETSC_DEFAULT or set lambda = 1.0 if
1080: no information about eigenvalues is available.
1082: Level: intermediate
1084: .seealso: PEPGetScale()
1085: @*/
1086: PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)1087: {
1093: pep->scale = scale;
1094: if (scale==PEP_SCALE_SCALAR || scale==PEP_SCALE_BOTH) {
1096: if (alpha == PETSC_DEFAULT || alpha == PETSC_DECIDE) {
1097: pep->sfactor = 0.0;
1098: pep->sfactor_set = PETSC_FALSE;
1099: } else {
1100: if (alpha<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
1101: pep->sfactor = alpha;
1102: pep->sfactor_set = PETSC_TRUE;
1103: }
1104: }
1105: if (scale==PEP_SCALE_DIAGONAL || scale==PEP_SCALE_BOTH) {
1106: if (Dl) {
1109: PetscObjectReference((PetscObject)Dl);
1110: VecDestroy(&pep->Dl);
1111: pep->Dl = Dl;
1112: }
1113: if (Dr) {
1116: PetscObjectReference((PetscObject)Dr);
1117: VecDestroy(&pep->Dr);
1118: pep->Dr = Dr;
1119: }
1122: if (its==PETSC_DECIDE || its==PETSC_DEFAULT) pep->sits = 5;
1123: else pep->sits = its;
1124: if (lambda==PETSC_DECIDE || lambda==PETSC_DEFAULT) pep->slambda = 1.0;
1125: else if (lambda<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of lambda. Must be > 0");
1126: else pep->slambda = lambda;
1127: }
1128: pep->state = PEP_STATE_INITIAL;
1129: return(0);
1130: }
1132: /*@C
1133: PEPGetScale - Gets the scaling strategy used by the PEP object, and the
1134: associated parameters.
1136: Not Collectiv, but vectors are shared by all processors that share the PEP1138: Input Parameter:
1139: . pep - the eigensolver context
1141: Output Parameters:
1142: + scale - scaling strategy
1143: . alpha - the scaling factor used in the scalar strategy
1144: . Dl - the left diagonal matrix of the diagonal scaling algorithm
1145: . Dr - the right diagonal matrix of the diagonal scaling algorithm
1146: . its - number of iterations of the diagonal scaling algorithm
1147: - lambda - approximation to wanted eigenvalues (modulus)
1149: Level: intermediate
1151: Note:
1152: The user can specify NULL for any parameter that is not needed.
1154: If Dl or Dr were not set by the user, then the ones computed internally are
1155: returned (or a null pointer if called before PEPSetUp).
1157: .seealso: PEPSetScale(), PEPSetUp()
1158: @*/
1159: PetscErrorCode PEPGetScale(PEP pep,PEPScale *scale,PetscReal *alpha,Vec *Dl,Vec *Dr,PetscInt *its,PetscReal *lambda)1160: {
1163: if (scale) *scale = pep->scale;
1164: if (alpha) *alpha = pep->sfactor;
1165: if (Dl) *Dl = pep->Dl;
1166: if (Dr) *Dr = pep->Dr;
1167: if (its) *its = pep->sits;
1168: if (lambda) *lambda = pep->slambda;
1169: return(0);
1170: }
1172: /*@
1173: PEPSetExtract - Specifies the extraction strategy to be used.
1175: Logically Collective on pep
1177: Input Parameters:
1178: + pep - the eigensolver context
1179: - extract - extraction strategy
1181: Options Database Keys:
1182: . -pep_extract <type> - extraction type, one of <none,norm,residual,structured>
1184: Level: intermediate
1186: .seealso: PEPGetExtract()
1187: @*/
1188: PetscErrorCode PEPSetExtract(PEP pep,PEPExtract extract)1189: {
1193: pep->extract = extract;
1194: return(0);
1195: }
1197: /*@
1198: PEPGetExtract - Gets the extraction strategy used by the PEP object.
1200: Not Collective
1202: Input Parameter:
1203: . pep - the eigensolver context
1205: Output Parameter:
1206: . extract - extraction strategy
1208: Level: intermediate
1210: .seealso: PEPSetExtract(), PEPExtract1211: @*/
1212: PetscErrorCode PEPGetExtract(PEP pep,PEPExtract *extract)1213: {
1217: *extract = pep->extract;
1218: return(0);
1219: }
1221: /*@
1222: PEPSetRefine - Specifies the refinement type (and options) to be used
1223: after the solve.
1225: Logically Collective on pep
1227: Input Parameters:
1228: + pep - the polynomial eigensolver context
1229: . refine - refinement type
1230: . npart - number of partitions of the communicator
1231: . tol - the convergence tolerance
1232: . its - maximum number of refinement iterations
1233: - scheme - which scheme to be used for solving the involved linear systems
1235: Options Database Keys:
1236: + -pep_refine <type> - refinement type, one of <none,simple,multiple>
1237: . -pep_refine_partitions <n> - the number of partitions
1238: . -pep_refine_tol <tol> - the tolerance
1239: . -pep_refine_its <its> - number of iterations
1240: - -pep_refine_scheme - to set the scheme for the linear solves
1242: Notes:
1243: By default, iterative refinement is disabled, since it may be very
1244: costly. There are two possible refinement strategies: simple and multiple.
1245: The simple approach performs iterative refinement on each of the
1246: converged eigenpairs individually, whereas the multiple strategy works
1247: with the invariant pair as a whole, refining all eigenpairs simultaneously.
1248: The latter may be required for the case of multiple eigenvalues.
1250: In some cases, especially when using direct solvers within the
1251: iterative refinement method, it may be helpful for improved scalability
1252: to split the communicator in several partitions. The npart parameter
1253: indicates how many partitions to use (defaults to 1).
1255: The tol and its parameters specify the stopping criterion. In the simple
1256: method, refinement continues until the residual of each eigenpair is
1257: below the tolerance (tol defaults to the PEP tol, but may be set to a
1258: different value). In contrast, the multiple method simply performs its
1259: refinement iterations (just one by default).
1261: The scheme argument is used to change the way in which linear systems are
1262: solved. Possible choices are: explicit, mixed block elimination (MBE),
1263: and Schur complement.
1265: Level: intermediate
1267: .seealso: PEPGetRefine()
1268: @*/
1269: PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PEPRefineScheme scheme)1270: {
1272: PetscMPIInt size;
1281: pep->refine = refine;
1282: if (refine) { /* process parameters only if not REFINE_NONE */
1283: if (npart!=pep->npart) {
1284: PetscSubcommDestroy(&pep->refinesubc);
1285: KSPDestroy(&pep->refineksp);
1286: }
1287: if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
1288: pep->npart = 1;
1289: } else {
1290: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size);
1291: if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
1292: pep->npart = npart;
1293: }
1294: if (tol == PETSC_DEFAULT || tol == PETSC_DECIDE) {
1295: pep->rtol = PETSC_DEFAULT;
1296: } else {
1297: if (tol<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1298: pep->rtol = tol;
1299: }
1300: if (its==PETSC_DECIDE || its==PETSC_DEFAULT) {
1301: pep->rits = PETSC_DEFAULT;
1302: } else {
1303: if (its<0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
1304: pep->rits = its;
1305: }
1306: pep->scheme = scheme;
1307: }
1308: pep->state = PEP_STATE_INITIAL;
1309: return(0);
1310: }
1312: /*@C
1313: PEPGetRefine - Gets the refinement strategy used by the PEP object, and the
1314: associated parameters.
1316: Not Collective
1318: Input Parameter:
1319: . pep - the polynomial eigensolver context
1321: Output Parameters:
1322: + refine - refinement type
1323: . npart - number of partitions of the communicator
1324: . tol - the convergence tolerance
1325: . its - maximum number of refinement iterations
1326: - scheme - the scheme used for solving linear systems
1328: Level: intermediate
1330: Note:
1331: The user can specify NULL for any parameter that is not needed.
1333: .seealso: PEPSetRefine()
1334: @*/
1335: PetscErrorCode PEPGetRefine(PEP pep,PEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,PEPRefineScheme *scheme)1336: {
1339: if (refine) *refine = pep->refine;
1340: if (npart) *npart = pep->npart;
1341: if (tol) *tol = pep->rtol;
1342: if (its) *its = pep->rits;
1343: if (scheme) *scheme = pep->scheme;
1344: return(0);
1345: }
1347: /*@C
1348: PEPSetOptionsPrefix - Sets the prefix used for searching for all
1349: PEP options in the database.
1351: Logically Collective on pep
1353: Input Parameters:
1354: + pep - the polynomial eigensolver context
1355: - prefix - the prefix string to prepend to all PEP option requests
1357: Notes:
1358: A hyphen (-) must NOT be given at the beginning of the prefix name.
1359: The first character of all runtime options is AUTOMATICALLY the
1360: hyphen.
1362: For example, to distinguish between the runtime options for two
1363: different PEP contexts, one could call
1364: .vb
1365: PEPSetOptionsPrefix(pep1,"qeig1_")
1366: PEPSetOptionsPrefix(pep2,"qeig2_")
1367: .ve
1369: Level: advanced
1371: .seealso: PEPAppendOptionsPrefix(), PEPGetOptionsPrefix()
1372: @*/
1373: PetscErrorCode PEPSetOptionsPrefix(PEP pep,const char *prefix)1374: {
1379: if (!pep->st) { PEPGetST(pep,&pep->st); }
1380: STSetOptionsPrefix(pep->st,prefix);
1381: if (!pep->V) { PEPGetBV(pep,&pep->V); }
1382: BVSetOptionsPrefix(pep->V,prefix);
1383: if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1384: DSSetOptionsPrefix(pep->ds,prefix);
1385: if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1386: RGSetOptionsPrefix(pep->rg,prefix);
1387: PetscObjectSetOptionsPrefix((PetscObject)pep,prefix);
1388: return(0);
1389: }
1391: /*@C
1392: PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1393: PEP options in the database.
1395: Logically Collective on pep
1397: Input Parameters:
1398: + pep - the polynomial eigensolver context
1399: - prefix - the prefix string to prepend to all PEP option requests
1401: Notes:
1402: A hyphen (-) must NOT be given at the beginning of the prefix name.
1403: The first character of all runtime options is AUTOMATICALLY the hyphen.
1405: Level: advanced
1407: .seealso: PEPSetOptionsPrefix(), PEPGetOptionsPrefix()
1408: @*/
1409: PetscErrorCode PEPAppendOptionsPrefix(PEP pep,const char *prefix)1410: {
1415: if (!pep->st) { PEPGetST(pep,&pep->st); }
1416: STAppendOptionsPrefix(pep->st,prefix);
1417: if (!pep->V) { PEPGetBV(pep,&pep->V); }
1418: BVAppendOptionsPrefix(pep->V,prefix);
1419: if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1420: DSAppendOptionsPrefix(pep->ds,prefix);
1421: if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1422: RGAppendOptionsPrefix(pep->rg,prefix);
1423: PetscObjectAppendOptionsPrefix((PetscObject)pep,prefix);
1424: return(0);
1425: }
1427: /*@C
1428: PEPGetOptionsPrefix - Gets the prefix used for searching for all
1429: PEP options in the database.
1431: Not Collective
1433: Input Parameters:
1434: . pep - the polynomial eigensolver context
1436: Output Parameters:
1437: . prefix - pointer to the prefix string used is returned
1439: Note:
1440: On the Fortran side, the user should pass in a string 'prefix' of
1441: sufficient length to hold the prefix.
1443: Level: advanced
1445: .seealso: PEPSetOptionsPrefix(), PEPAppendOptionsPrefix()
1446: @*/
1447: PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])1448: {
1454: PetscObjectGetOptionsPrefix((PetscObject)pep,prefix);
1455: return(0);
1456: }