1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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> /*I "slepcpep.h" I*/
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)->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)->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,256,&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? pep->max_it: PETSC_DEFAULT;
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? pep->ncv: PETSC_DEFAULT;
188: PetscOptionsInt("-pep_ncv","Number of basis vectors","PEPSetDimensions",pep->ncv,&j,&flg2);
189: k = pep->mpd? pep->mpd: PETSC_DEFAULT;
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: PetscOptionsBoolGroupEnd("-pep_target_imaginary","Compute eigenvalues with imaginary parts closest to target","PEPSetWhichEigenpairs",&flg);
212: if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_IMAGINARY); }
214: PetscOptionsScalar("-pep_target","Value of the target","PEPSetTarget",pep->target,&s,&flg);
215: if (flg) {
216: if (pep->which!=PEP_TARGET_REAL && pep->which!=PEP_TARGET_IMAGINARY) {
217: PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
218: }
219: PEPSetTarget(pep,s);
220: }
222: k = 2;
223: PetscOptionsRealArray("-pep_interval","Computational interval (two real values separated with a comma without spaces)","PEPSetInterval",array,&k,&flg);
224: if (flg) {
225: if (k<2) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_SIZ,"Must pass two values in -pep_interval (comma-separated without spaces)");
226: PEPSetWhichEigenpairs(pep,PEP_ALL);
227: PEPSetInterval(pep,array[0],array[1]);
228: }
230: /* -----------------------------------------------------------------------*/
231: /*
232: Cancels all monitors hardwired into code before call to PEPSetFromOptions()
233: */
234: PetscOptionsBool("-pep_monitor_cancel","Remove any hardwired monitor routines","PEPMonitorCancel",PETSC_FALSE,&flg,&set);
235: if (set && flg) {
236: PEPMonitorCancel(pep);
237: }
238: /*
239: Text monitors
240: */
241: PEPMonitorSetFromOptions(pep,"-pep_monitor","Monitor first unconverged approximate eigenvalue and error estimate","PEPMonitorFirst",PEPMonitorFirst,PETSC_FALSE);
242: PEPConvMonitorSetFromOptions(pep,"-pep_monitor_conv","Monitor approximate eigenvalues and error estimates as they converge","PEPMonitorConverged",PEPMonitorConverged);
243: PEPMonitorSetFromOptions(pep,"-pep_monitor_all","Monitor approximate eigenvalues and error estimates","PEPMonitorAll",PEPMonitorAll,PETSC_TRUE);
244: /*
245: Line graph monitors
246: */
247: PetscOptionsBool("-pep_monitor_lg","Monitor first unconverged approximate error estimate graphically","PEPMonitorSet",PETSC_FALSE,&flg,&set);
248: if (set && flg) {
249: PEPMonitorLGCreate(PetscObjectComm((PetscObject)pep),NULL,"Error estimates",PETSC_DECIDE,PETSC_DECIDE,300,300,&lg);
250: PEPMonitorSet(pep,PEPMonitorLG,lg,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
251: }
252: PetscOptionsBool("-pep_monitor_lg_all","Monitor error estimates graphically","PEPMonitorSet",PETSC_FALSE,&flg,&set);
253: if (set && flg) {
254: PEPMonitorLGCreate(PetscObjectComm((PetscObject)pep),NULL,"Error estimates",PETSC_DECIDE,PETSC_DECIDE,300,300,&lg);
255: PEPMonitorSet(pep,PEPMonitorLGAll,lg,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
256: PEPSetTrackAll(pep,PETSC_TRUE);
257: }
259: /* -----------------------------------------------------------------------*/
260: PetscOptionsName("-pep_view","Print detailed information on solver used","PEPView",NULL);
261: PetscOptionsName("-pep_view_vectors","View computed eigenvectors","PEPVectorsView",NULL);
262: PetscOptionsName("-pep_view_values","View computed eigenvalues","PEPValuesView",NULL);
263: PetscOptionsName("-pep_converged_reason","Print reason for convergence, and number of iterations","PEPReasonView",NULL);
264: PetscOptionsName("-pep_error_absolute","Print absolute errors of each eigenpair","PEPErrorView",NULL);
265: PetscOptionsName("-pep_error_relative","Print relative errors of each eigenpair","PEPErrorView",NULL);
266: PetscOptionsName("-pep_error_backward","Print backward errors of each eigenpair","PEPErrorView",NULL);
268: if (pep->ops->setfromoptions) {
269: (*pep->ops->setfromoptions)(PetscOptionsObject,pep);
270: }
271: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)pep);
272: PetscOptionsEnd();
274: if (!pep->V) { PEPGetBV(pep,&pep->V); }
275: BVSetFromOptions(pep->V);
276: if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
277: RGSetFromOptions(pep->rg);
278: if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
279: DSSetFromOptions(pep->ds);
280: if (!pep->st) { PEPGetST(pep,&pep->st); }
281: PEPSetDefaultST(pep);
282: STSetFromOptions(pep->st);
283: if (!pep->refineksp) { PEPRefineGetKSP(pep,&pep->refineksp); }
284: KSPSetFromOptions(pep->refineksp);
285: return(0);
286: }
288: /*@C
289: PEPGetTolerances - Gets the tolerance and maximum iteration count used
290: by the PEP convergence tests.
292: Not Collective
294: Input Parameter:
295: . pep - the polynomial eigensolver context
297: Output Parameters:
298: + tol - the convergence tolerance
299: - maxits - maximum number of iterations
301: Notes:
302: The user can specify NULL for any parameter that is not needed.
304: Level: intermediate
306: .seealso: PEPSetTolerances()
307: @*/
308: PetscErrorCode PEPGetTolerances(PEP pep,PetscReal *tol,PetscInt *maxits)309: {
312: if (tol) *tol = pep->tol;
313: if (maxits) *maxits = pep->max_it;
314: return(0);
315: }
317: /*@
318: PEPSetTolerances - Sets the tolerance and maximum iteration count used
319: by the PEP convergence tests.
321: Logically Collective on PEP323: Input Parameters:
324: + pep - the polynomial eigensolver context
325: . tol - the convergence tolerance
326: - maxits - maximum number of iterations to use
328: Options Database Keys:
329: + -pep_tol <tol> - Sets the convergence tolerance
330: - -pep_max_it <maxits> - Sets the maximum number of iterations allowed
332: Notes:
333: Use PETSC_DEFAULT for either argument to assign a reasonably good value.
335: Level: intermediate
337: .seealso: PEPGetTolerances()
338: @*/
339: PetscErrorCode PEPSetTolerances(PEP pep,PetscReal tol,PetscInt maxits)340: {
345: if (tol == PETSC_DEFAULT) {
346: pep->tol = PETSC_DEFAULT;
347: pep->state = PEP_STATE_INITIAL;
348: } else {
349: if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
350: pep->tol = tol;
351: }
352: if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
353: pep->max_it = 0;
354: pep->state = PEP_STATE_INITIAL;
355: } else {
356: if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
357: pep->max_it = maxits;
358: }
359: return(0);
360: }
362: /*@C
363: PEPGetDimensions - Gets the number of eigenvalues to compute
364: and the dimension of the subspace.
366: Not Collective
368: Input Parameter:
369: . pep - the polynomial eigensolver context
371: Output Parameters:
372: + nev - number of eigenvalues to compute
373: . ncv - the maximum dimension of the subspace to be used by the solver
374: - mpd - the maximum dimension allowed for the projected problem
376: Notes:
377: The user can specify NULL for any parameter that is not needed.
379: Level: intermediate
381: .seealso: PEPSetDimensions()
382: @*/
383: PetscErrorCode PEPGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)384: {
387: if (nev) *nev = pep->nev;
388: if (ncv) *ncv = pep->ncv;
389: if (mpd) *mpd = pep->mpd;
390: return(0);
391: }
393: /*@
394: PEPSetDimensions - Sets the number of eigenvalues to compute
395: and the dimension of the subspace.
397: Logically Collective on PEP399: Input Parameters:
400: + pep - the polynomial eigensolver context
401: . nev - number of eigenvalues to compute
402: . ncv - the maximum dimension of the subspace to be used by the solver
403: - mpd - the maximum dimension allowed for the projected problem
405: Options Database Keys:
406: + -pep_nev <nev> - Sets the number of eigenvalues
407: . -pep_ncv <ncv> - Sets the dimension of the subspace
408: - -pep_mpd <mpd> - Sets the maximum projected dimension
410: Notes:
411: Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
412: dependent on the solution method.
414: The parameters ncv and mpd are intimately related, so that the user is advised
415: to set one of them at most. Normal usage is that
416: (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
417: (b) in cases where nev is large, the user sets mpd.
419: The value of ncv should always be between nev and (nev+mpd), typically
420: ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
421: a smaller value should be used.
423: Level: intermediate
425: .seealso: PEPGetDimensions()
426: @*/
427: PetscErrorCode PEPSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)428: {
434: if (nev<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
435: pep->nev = nev;
436: if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
437: pep->ncv = 0;
438: } else {
439: if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
440: pep->ncv = ncv;
441: }
442: if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
443: pep->mpd = 0;
444: } else {
445: if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
446: pep->mpd = mpd;
447: }
448: pep->state = PEP_STATE_INITIAL;
449: return(0);
450: }
452: /*@
453: PEPSetWhichEigenpairs - Specifies which portion of the spectrum is
454: to be sought.
456: Logically Collective on PEP458: Input Parameters:
459: + pep - eigensolver context obtained from PEPCreate()
460: - which - the portion of the spectrum to be sought
462: Possible values:
463: The parameter 'which' can have one of these values
465: + PEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
466: . PEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
467: . PEP_LARGEST_REAL - largest real parts
468: . PEP_SMALLEST_REAL - smallest real parts
469: . PEP_LARGEST_IMAGINARY - largest imaginary parts
470: . PEP_SMALLEST_IMAGINARY - smallest imaginary parts
471: . PEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
472: . PEP_TARGET_REAL - eigenvalues with real part closest to target
473: . PEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
474: - PEP_WHICH_USER - user defined ordering set with PEPSetEigenvalueComparison()
476: Options Database Keys:
477: + -pep_largest_magnitude - Sets largest eigenvalues in magnitude
478: . -pep_smallest_magnitude - Sets smallest eigenvalues in magnitude
479: . -pep_largest_real - Sets largest real parts
480: . -pep_smallest_real - Sets smallest real parts
481: . -pep_largest_imaginary - Sets largest imaginary parts
482: . -pep_smallest_imaginary - Sets smallest imaginary parts
483: . -pep_target_magnitude - Sets eigenvalues closest to target
484: . -pep_target_real - Sets real parts closest to target
485: - -pep_target_imaginary - Sets imaginary parts closest to target
487: Notes:
488: Not all eigensolvers implemented in PEP account for all the possible values
489: stated above. If SLEPc is compiled for real numbers PEP_LARGEST_IMAGINARY490: and PEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
491: for eigenvalue selection.
493: The target is a scalar value provided with PEPSetTarget().
495: The criterion PEP_TARGET_IMAGINARY is available only in case PETSc and
496: SLEPc have been built with complex scalars.
498: Level: intermediate
500: .seealso: PEPGetWhichEigenpairs(), PEPSetTarget(), PEPSetEigenvalueComparison(), PEPWhich501: @*/
502: PetscErrorCode PEPSetWhichEigenpairs(PEP pep,PEPWhich which)503: {
507: switch (which) {
508: case PEP_LARGEST_MAGNITUDE:
509: case PEP_SMALLEST_MAGNITUDE:
510: case PEP_LARGEST_REAL:
511: case PEP_SMALLEST_REAL:
512: case PEP_LARGEST_IMAGINARY:
513: case PEP_SMALLEST_IMAGINARY:
514: case PEP_TARGET_MAGNITUDE:
515: case PEP_TARGET_REAL:
516: #if defined(PETSC_USE_COMPLEX)
517: case PEP_TARGET_IMAGINARY:
518: #endif
519: case PEP_ALL:
520: case PEP_WHICH_USER:
521: if (pep->which != which) {
522: pep->state = PEP_STATE_INITIAL;
523: pep->which = which;
524: }
525: break;
526: default:527: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
528: }
529: return(0);
530: }
532: /*@
533: PEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
534: sought.
536: Not Collective
538: Input Parameter:
539: . pep - eigensolver context obtained from PEPCreate()
541: Output Parameter:
542: . which - the portion of the spectrum to be sought
544: Notes:
545: See PEPSetWhichEigenpairs() for possible values of 'which'.
547: Level: intermediate
549: .seealso: PEPSetWhichEigenpairs(), PEPWhich550: @*/
551: PetscErrorCode PEPGetWhichEigenpairs(PEP pep,PEPWhich *which)552: {
556: *which = pep->which;
557: return(0);
558: }
560: /*@C
561: PEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
562: when PEPSetWhichEigenpairs() is set to PEP_WHICH_USER.
564: Logically Collective on PEP566: Input Parameters:
567: + pep - eigensolver context obtained from PEPCreate()
568: . func - a pointer to the comparison function
569: - ctx - a context pointer (the last parameter to the comparison function)
571: Calling Sequence of func:
572: $ func(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
574: + ar - real part of the 1st eigenvalue
575: . ai - imaginary part of the 1st eigenvalue
576: . br - real part of the 2nd eigenvalue
577: . bi - imaginary part of the 2nd eigenvalue
578: . res - result of comparison
579: - ctx - optional context, as set by PEPSetEigenvalueComparison()
581: Note:
582: The returning parameter 'res' can be
583: + negative - if the 1st eigenvalue is preferred to the 2st one
584: . zero - if both eigenvalues are equally preferred
585: - positive - if the 2st eigenvalue is preferred to the 1st one
587: Level: advanced
589: .seealso: PEPSetWhichEigenpairs(), PEPWhich590: @*/
591: PetscErrorCode PEPSetEigenvalueComparison(PEP pep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)592: {
595: pep->sc->comparison = func;
596: pep->sc->comparisonctx = ctx;
597: pep->which = PEP_WHICH_USER;
598: return(0);
599: }
601: /*@
602: PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.
604: Logically Collective on PEP606: Input Parameters:
607: + pep - the polynomial eigensolver context
608: - type - a known type of polynomial eigenvalue problem
610: Options Database Keys:
611: + -pep_general - general problem with no particular structure
612: . -pep_hermitian - problem whose coefficient matrices are Hermitian
613: . -pep_hyperbolic - Hermitian problem that satisfies the definition of hyperbolic
614: - -pep_gyroscopic - problem with Hamiltonian structure
616: Notes:
617: Allowed values for the problem type are: general (PEP_GENERAL), Hermitian
618: (PEP_HERMITIAN), hyperbolic (PEP_HYPERBOLIC), and gyroscopic (PEP_GYROSCOPIC).
620: This function is used to instruct SLEPc to exploit certain structure in
621: the polynomial eigenproblem. By default, no particular structure is assumed.
623: If the problem matrices are Hermitian (symmetric in the real case) or
624: Hermitian/skew-Hermitian then the solver can exploit this fact to perform
625: less operations or provide better stability. Hyperbolic problems are a
626: particular case of Hermitian problems, some solvers may treat them simply as
627: Hermitian.
629: Level: intermediate
631: .seealso: PEPSetOperators(), PEPSetType(), PEPGetProblemType(), PEPProblemType632: @*/
633: PetscErrorCode PEPSetProblemType(PEP pep,PEPProblemType type)634: {
638: 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");
639: if (type != pep->problem_type) {
640: pep->problem_type = type;
641: pep->state = PEP_STATE_INITIAL;
642: }
643: return(0);
644: }
646: /*@
647: PEPGetProblemType - Gets the problem type from the PEP object.
649: Not Collective
651: Input Parameter:
652: . pep - the polynomial eigensolver context
654: Output Parameter:
655: . type - the problem type
657: Level: intermediate
659: .seealso: PEPSetProblemType(), PEPProblemType660: @*/
661: PetscErrorCode PEPGetProblemType(PEP pep,PEPProblemType *type)662: {
666: *type = pep->problem_type;
667: return(0);
668: }
670: /*@
671: PEPSetBasis - Specifies the type of polynomial basis used to describe the
672: polynomial eigenvalue problem.
674: Logically Collective on PEP676: Input Parameters:
677: + pep - the polynomial eigensolver context
678: - basis - the type of polynomial basis
680: Options Database Key:
681: . -pep_basis <basis> - Select the basis type
683: Notes:
684: By default, the coefficient matrices passed via PEPSetOperators() are
685: expressed in the monomial basis, i.e.
686: P(lambda) = A_0 + lambda*A_1 + lambda^2*A_2 + ... + lambda^d*A_d.
687: Other polynomial bases may have better numerical behaviour, but the user
688: must then pass the coefficient matrices accordingly.
690: Level: intermediate
692: .seealso: PEPSetOperators(), PEPGetBasis(), PEPBasis693: @*/
694: PetscErrorCode PEPSetBasis(PEP pep,PEPBasis basis)695: {
699: pep->basis = basis;
700: return(0);
701: }
703: /*@
704: PEPGetBasis - Gets the type of polynomial basis from the PEP object.
706: Not Collective
708: Input Parameter:
709: . pep - the polynomial eigensolver context
711: Output Parameter:
712: . basis - the polynomial basis
714: Level: intermediate
716: .seealso: PEPSetBasis(), PEPBasis717: @*/
718: PetscErrorCode PEPGetBasis(PEP pep,PEPBasis *basis)719: {
723: *basis = pep->basis;
724: return(0);
725: }
727: /*@
728: PEPSetTrackAll - Specifies if the solver must compute the residual of all
729: approximate eigenpairs or not.
731: Logically Collective on PEP733: Input Parameters:
734: + pep - the eigensolver context
735: - trackall - whether compute all residuals or not
737: Notes:
738: If the user sets trackall=PETSC_TRUE then the solver explicitly computes
739: the residual for each eigenpair approximation. Computing the residual is
740: usually an expensive operation and solvers commonly compute the associated
741: residual to the first unconverged eigenpair.
742: The options '-pep_monitor_all' and '-pep_monitor_lg_all' automatically
743: activate this option.
745: Level: developer
747: .seealso: PEPGetTrackAll()
748: @*/
749: PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)750: {
754: pep->trackall = trackall;
755: return(0);
756: }
758: /*@
759: PEPGetTrackAll - Returns the flag indicating whether all residual norms must
760: be computed or not.
762: Not Collective
764: Input Parameter:
765: . pep - the eigensolver context
767: Output Parameter:
768: . trackall - the returned flag
770: Level: developer
772: .seealso: PEPSetTrackAll()
773: @*/
774: PetscErrorCode PEPGetTrackAll(PEP pep,PetscBool *trackall)775: {
779: *trackall = pep->trackall;
780: return(0);
781: }
783: /*@C
784: PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
785: used in the convergence test.
787: Logically Collective on PEP789: Input Parameters:
790: + pep - eigensolver context obtained from PEPCreate()
791: . func - a pointer to the convergence test function
792: . ctx - context for private data for the convergence routine (may be null)
793: - destroy - a routine for destroying the context (may be null)
795: Calling Sequence of func:
796: $ func(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
798: + pep - eigensolver context obtained from PEPCreate()
799: . eigr - real part of the eigenvalue
800: . eigi - imaginary part of the eigenvalue
801: . res - residual norm associated to the eigenpair
802: . errest - (output) computed error estimate
803: - ctx - optional context, as set by PEPSetConvergenceTestFunction()
805: Note:
806: If the error estimate returned by the convergence test function is less than
807: the tolerance, then the eigenvalue is accepted as converged.
809: Level: advanced
811: .seealso: PEPSetConvergenceTest(), PEPSetTolerances()
812: @*/
813: PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))814: {
819: if (pep->convergeddestroy) {
820: (*pep->convergeddestroy)(pep->convergedctx);
821: }
822: pep->convergeduser = func;
823: pep->convergeddestroy = destroy;
824: pep->convergedctx = ctx;
825: if (func == PEPConvergedRelative) pep->conv = PEP_CONV_REL;
826: else if (func == PEPConvergedNorm) pep->conv = PEP_CONV_NORM;
827: else if (func == PEPConvergedAbsolute) pep->conv = PEP_CONV_ABS;
828: else {
829: pep->conv = PEP_CONV_USER;
830: pep->converged = pep->convergeduser;
831: }
832: return(0);
833: }
835: /*@
836: PEPSetConvergenceTest - Specifies how to compute the error estimate
837: used in the convergence test.
839: Logically Collective on PEP841: Input Parameters:
842: + pep - eigensolver context obtained from PEPCreate()
843: - conv - the type of convergence test
845: Options Database Keys:
846: + -pep_conv_abs - Sets the absolute convergence test
847: . -pep_conv_rel - Sets the convergence test relative to the eigenvalue
848: . -pep_conv_norm - Sets the convergence test relative to the matrix norms
849: - -pep_conv_user - Selects the user-defined convergence test
851: Note:
852: The parameter 'conv' can have one of these values
853: + PEP_CONV_ABS - absolute error ||r||
854: . PEP_CONV_REL - error relative to the eigenvalue l, ||r||/|l|
855: . PEP_CONV_NORM - error relative matrix norms, ||r||/sum_i(l^i*||A_i||)
856: - PEP_CONV_USER - function set by PEPSetConvergenceTestFunction()
858: Level: intermediate
860: .seealso: PEPGetConvergenceTest(), PEPSetConvergenceTestFunction(), PEPSetStoppingTest(), PEPConv861: @*/
862: PetscErrorCode PEPSetConvergenceTest(PEP pep,PEPConv conv)863: {
867: switch (conv) {
868: case PEP_CONV_ABS: pep->converged = PEPConvergedAbsolute; break;
869: case PEP_CONV_REL: pep->converged = PEPConvergedRelative; break;
870: case PEP_CONV_NORM: pep->converged = PEPConvergedNorm; break;
871: case PEP_CONV_USER:
872: if (!pep->convergeduser) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetConvergenceTestFunction() first");
873: pep->converged = pep->convergeduser;
874: break;
875: default:876: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
877: }
878: pep->conv = conv;
879: return(0);
880: }
882: /*@
883: PEPGetConvergenceTest - Gets the method used to compute the error estimate
884: used in the convergence test.
886: Not Collective
888: Input Parameters:
889: . pep - eigensolver context obtained from PEPCreate()
891: Output Parameters:
892: . conv - the type of convergence test
894: Level: intermediate
896: .seealso: PEPSetConvergenceTest(), PEPConv897: @*/
898: PetscErrorCode PEPGetConvergenceTest(PEP pep,PEPConv *conv)899: {
903: *conv = pep->conv;
904: return(0);
905: }
907: /*@C
908: PEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
909: iteration of the eigensolver.
911: Logically Collective on PEP913: Input Parameters:
914: + pep - eigensolver context obtained from PEPCreate()
915: . func - pointer to the stopping test function
916: . ctx - context for private data for the stopping routine (may be null)
917: - destroy - a routine for destroying the context (may be null)
919: Calling Sequence of func:
920: $ func(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ctx)
922: + pep - eigensolver context obtained from PEPCreate()
923: . its - current number of iterations
924: . max_it - maximum number of iterations
925: . nconv - number of currently converged eigenpairs
926: . nev - number of requested eigenpairs
927: . reason - (output) result of the stopping test
928: - ctx - optional context, as set by PEPSetStoppingTestFunction()
930: Note:
931: Normal usage is to first call the default routine PEPStoppingBasic() and then
932: set reason to PEP_CONVERGED_USER if some user-defined conditions have been
933: met. To let the eigensolver continue iterating, the result must be left as
934: PEP_CONVERGED_ITERATING.
936: Level: advanced
938: .seealso: PEPSetStoppingTest(), PEPStoppingBasic()
939: @*/
940: PetscErrorCode PEPSetStoppingTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))941: {
946: if (pep->stoppingdestroy) {
947: (*pep->stoppingdestroy)(pep->stoppingctx);
948: }
949: pep->stoppinguser = func;
950: pep->stoppingdestroy = destroy;
951: pep->stoppingctx = ctx;
952: if (func == PEPStoppingBasic) pep->stop = PEP_STOP_BASIC;
953: else {
954: pep->stop = PEP_STOP_USER;
955: pep->stopping = pep->stoppinguser;
956: }
957: return(0);
958: }
960: /*@
961: PEPSetStoppingTest - Specifies how to decide the termination of the outer
962: loop of the eigensolver.
964: Logically Collective on PEP966: Input Parameters:
967: + pep - eigensolver context obtained from PEPCreate()
968: - stop - the type of stopping test
970: Options Database Keys:
971: + -pep_stop_basic - Sets the default stopping test
972: - -pep_stop_user - Selects the user-defined stopping test
974: Note:
975: The parameter 'stop' can have one of these values
976: + PEP_STOP_BASIC - default stopping test
977: - PEP_STOP_USER - function set by PEPSetStoppingTestFunction()
979: Level: advanced
981: .seealso: PEPGetStoppingTest(), PEPSetStoppingTestFunction(), PEPSetConvergenceTest(), PEPStop982: @*/
983: PetscErrorCode PEPSetStoppingTest(PEP pep,PEPStop stop)984: {
988: switch (stop) {
989: case PEP_STOP_BASIC: pep->stopping = PEPStoppingBasic; break;
990: case PEP_STOP_USER:
991: if (!pep->stoppinguser) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetStoppingTestFunction() first");
992: pep->stopping = pep->stoppinguser;
993: break;
994: default:995: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
996: }
997: pep->stop = stop;
998: return(0);
999: }
1001: /*@
1002: PEPGetStoppingTest - Gets the method used to decide the termination of the outer
1003: loop of the eigensolver.
1005: Not Collective
1007: Input Parameters:
1008: . pep - eigensolver context obtained from PEPCreate()
1010: Output Parameters:
1011: . stop - the type of stopping test
1013: Level: advanced
1015: .seealso: PEPSetStoppingTest(), PEPStop1016: @*/
1017: PetscErrorCode PEPGetStoppingTest(PEP pep,PEPStop *stop)1018: {
1022: *stop = pep->stop;
1023: return(0);
1024: }
1026: /*@C
1027: PEPSetScale - Specifies the scaling strategy to be used.
1029: Logically Collective on PEP1031: Input Parameters:
1032: + pep - the eigensolver context
1033: . scale - scaling strategy
1034: . alpha - the scaling factor used in the scalar strategy
1035: . Dl - the left diagonal matrix of the diagonal scaling algorithm
1036: . Dr - the right diagonal matrix of the diagonal scaling algorithm
1037: . its - number of iterations of the diagonal scaling algorithm
1038: - lambda - approximation to wanted eigenvalues (modulus)
1040: Options Database Keys:
1041: + -pep_scale <type> - scaling type, one of <none,scalar,diagonal,both>
1042: . -pep_scale_factor <alpha> - the scaling factor
1043: . -pep_scale_its <its> - number of iterations
1044: - -pep_scale_lambda <lambda> - approximation to eigenvalues
1046: Notes:
1047: There are two non-exclusive scaling strategies: scalar and diagonal.
1049: In the scalar strategy, scaling is applied to the eigenvalue, that is,
1050: mu = lambda/alpha is the new eigenvalue and all matrices are scaled
1051: accordingly. After solving the scaled problem, the original lambda is
1052: recovered. Parameter 'alpha' must be positive. Use PETSC_DECIDE to let
1053: the solver compute a reasonable scaling factor.
1055: In the diagonal strategy, the solver works implicitly with matrix Dl*A*Dr,
1056: where Dl and Dr are appropriate diagonal matrices. This improves the accuracy
1057: of the computed results in some cases. The user may provide the Dr and Dl
1058: matrices represented as Vec objects storing diagonal elements. If not
1059: provided, these matrices are computed internally. This option requires
1060: that the polynomial coefficient matrices are of MATAIJ type.
1061: The parameter 'its' is the number of iterations performed by the method.
1062: Parameter 'lambda' must be positive. Use PETSC_DECIDE or set lambda = 1.0 if
1063: no information about eigenvalues is available.
1065: Level: intermediate
1067: .seealso: PEPGetScale()
1068: @*/
1069: PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)1070: {
1076: pep->scale = scale;
1077: if (scale==PEP_SCALE_SCALAR || scale==PEP_SCALE_BOTH) {
1079: if (alpha == PETSC_DEFAULT || alpha == PETSC_DECIDE) {
1080: pep->sfactor = 0.0;
1081: pep->sfactor_set = PETSC_FALSE;
1082: } else {
1083: if (alpha<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
1084: pep->sfactor = alpha;
1085: pep->sfactor_set = PETSC_TRUE;
1086: }
1087: }
1088: if (scale==PEP_SCALE_DIAGONAL || scale==PEP_SCALE_BOTH) {
1089: if (Dl) {
1092: PetscObjectReference((PetscObject)Dl);
1093: VecDestroy(&pep->Dl);
1094: pep->Dl = Dl;
1095: }
1096: if (Dr) {
1099: PetscObjectReference((PetscObject)Dr);
1100: VecDestroy(&pep->Dr);
1101: pep->Dr = Dr;
1102: }
1105: if (its==PETSC_DECIDE || its==PETSC_DEFAULT) pep->sits = 5;
1106: else pep->sits = its;
1107: if (lambda==PETSC_DECIDE || lambda==PETSC_DEFAULT) pep->slambda = 1.0;
1108: else if (lambda<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of lambda. Must be > 0");
1109: else pep->slambda = lambda;
1110: }
1111: pep->state = PEP_STATE_INITIAL;
1112: return(0);
1113: }
1115: /*@C
1116: PEPGetScale - Gets the scaling strategy used by the PEP object, and the
1117: associated parameters.
1119: Not Collectiv, but vectors are shared by all processors that share the PEP1121: Input Parameter:
1122: . pep - the eigensolver context
1124: Output Parameters:
1125: + scale - scaling strategy
1126: . alpha - the scaling factor used in the scalar strategy
1127: . Dl - the left diagonal matrix of the diagonal scaling algorithm
1128: . Dr - the right diagonal matrix of the diagonal scaling algorithm
1129: . its - number of iterations of the diagonal scaling algorithm
1130: - lambda - approximation to wanted eigenvalues (modulus)
1132: Level: intermediate
1134: Note:
1135: The user can specify NULL for any parameter that is not needed.
1137: If Dl or Dr were not set by the user, then the ones computed internally are
1138: returned (or a null pointer if called before PEPSetUp).
1140: .seealso: PEPSetScale(), PEPSetUp()
1141: @*/
1142: PetscErrorCode PEPGetScale(PEP pep,PEPScale *scale,PetscReal *alpha,Vec *Dl,Vec *Dr,PetscInt *its,PetscReal *lambda)1143: {
1146: if (scale) *scale = pep->scale;
1147: if (alpha) *alpha = pep->sfactor;
1148: if (Dl) *Dl = pep->Dl;
1149: if (Dr) *Dr = pep->Dr;
1150: if (its) *its = pep->sits;
1151: if (lambda) *lambda = pep->slambda;
1152: return(0);
1153: }
1155: /*@
1156: PEPSetExtract - Specifies the extraction strategy to be used.
1158: Logically Collective on PEP1160: Input Parameters:
1161: + pep - the eigensolver context
1162: - extract - extraction strategy
1164: Options Database Keys:
1165: . -pep_extract <type> - extraction type, one of <none,norm,residual,structured>
1167: Level: intermediate
1169: .seealso: PEPGetExtract()
1170: @*/
1171: PetscErrorCode PEPSetExtract(PEP pep,PEPExtract extract)1172: {
1176: pep->extract = extract;
1177: return(0);
1178: }
1180: /*@
1181: PEPGetExtract - Gets the extraction strategy used by the PEP object.
1183: Not Collective
1185: Input Parameter:
1186: . pep - the eigensolver context
1188: Output Parameter:
1189: . extract - extraction strategy
1191: Level: intermediate
1193: .seealso: PEPSetExtract()
1194: @*/
1195: PetscErrorCode PEPGetExtract(PEP pep,PEPExtract *extract)1196: {
1199: if (extract) *extract = pep->extract;
1200: return(0);
1201: }
1203: /*@
1204: PEPSetRefine - Specifies the refinement type (and options) to be used
1205: after the solve.
1207: Logically Collective on PEP1209: Input Parameters:
1210: + pep - the polynomial eigensolver context
1211: . refine - refinement type
1212: . npart - number of partitions of the communicator
1213: . tol - the convergence tolerance
1214: . its - maximum number of refinement iterations
1215: - scheme - which scheme to be used for solving the involved linear systems
1217: Options Database Keys:
1218: + -pep_refine <type> - refinement type, one of <none,simple,multiple>
1219: . -pep_refine_partitions <n> - the number of partitions
1220: . -pep_refine_tol <tol> - the tolerance
1221: . -pep_refine_its <its> - number of iterations
1222: - -pep_refine_scheme - to set the scheme for the linear solves
1224: Notes:
1225: By default, iterative refinement is disabled, since it may be very
1226: costly. There are two possible refinement strategies: simple and multiple.
1227: The simple approach performs iterative refinement on each of the
1228: converged eigenpairs individually, whereas the multiple strategy works
1229: with the invariant pair as a whole, refining all eigenpairs simultaneously.
1230: The latter may be required for the case of multiple eigenvalues.
1232: In some cases, especially when using direct solvers within the
1233: iterative refinement method, it may be helpful for improved scalability
1234: to split the communicator in several partitions. The npart parameter
1235: indicates how many partitions to use (defaults to 1).
1237: The tol and its parameters specify the stopping criterion. In the simple
1238: method, refinement continues until the residual of each eigenpair is
1239: below the tolerance (tol defaults to the PEP tol, but may be set to a
1240: different value). In contrast, the multiple method simply performs its
1241: refinement iterations (just one by default).
1243: The scheme argument is used to change the way in which linear systems are
1244: solved. Possible choices are: explicit, mixed block elimination (MBE),
1245: and Schur complement.
1247: Level: intermediate
1249: .seealso: PEPGetRefine()
1250: @*/
1251: PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PEPRefineScheme scheme)1252: {
1254: PetscMPIInt size;
1263: pep->refine = refine;
1264: if (refine) { /* process parameters only if not REFINE_NONE */
1265: if (npart!=pep->npart) {
1266: PetscSubcommDestroy(&pep->refinesubc);
1267: KSPDestroy(&pep->refineksp);
1268: }
1269: if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
1270: pep->npart = 1;
1271: } else {
1272: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size);
1273: if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
1274: pep->npart = npart;
1275: }
1276: if (tol == PETSC_DEFAULT || tol == PETSC_DECIDE) {
1277: pep->rtol = PETSC_DEFAULT;
1278: } else {
1279: if (tol<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1280: pep->rtol = tol;
1281: }
1282: if (its==PETSC_DECIDE || its==PETSC_DEFAULT) {
1283: pep->rits = PETSC_DEFAULT;
1284: } else {
1285: if (its<0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
1286: pep->rits = its;
1287: }
1288: pep->scheme = scheme;
1289: }
1290: pep->state = PEP_STATE_INITIAL;
1291: return(0);
1292: }
1294: /*@C
1295: PEPGetRefine - Gets the refinement strategy used by the PEP object, and the
1296: associated parameters.
1298: Not Collective
1300: Input Parameter:
1301: . pep - the polynomial eigensolver context
1303: Output Parameters:
1304: + refine - refinement type
1305: . npart - number of partitions of the communicator
1306: . tol - the convergence tolerance
1307: . its - maximum number of refinement iterations
1308: - scheme - the scheme used for solving linear systems
1310: Level: intermediate
1312: Note:
1313: The user can specify NULL for any parameter that is not needed.
1315: .seealso: PEPSetRefine()
1316: @*/
1317: PetscErrorCode PEPGetRefine(PEP pep,PEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,PEPRefineScheme *scheme)1318: {
1321: if (refine) *refine = pep->refine;
1322: if (npart) *npart = pep->npart;
1323: if (tol) *tol = pep->rtol;
1324: if (its) *its = pep->rits;
1325: if (scheme) *scheme = pep->scheme;
1326: return(0);
1327: }
1329: /*@C
1330: PEPSetOptionsPrefix - Sets the prefix used for searching for all
1331: PEP options in the database.
1333: Logically Collective on PEP1335: Input Parameters:
1336: + pep - the polynomial eigensolver context
1337: - prefix - the prefix string to prepend to all PEP option requests
1339: Notes:
1340: A hyphen (-) must NOT be given at the beginning of the prefix name.
1341: The first character of all runtime options is AUTOMATICALLY the
1342: hyphen.
1344: For example, to distinguish between the runtime options for two
1345: different PEP contexts, one could call
1346: .vb
1347: PEPSetOptionsPrefix(pep1,"qeig1_")
1348: PEPSetOptionsPrefix(pep2,"qeig2_")
1349: .ve
1351: Level: advanced
1353: .seealso: PEPAppendOptionsPrefix(), PEPGetOptionsPrefix()
1354: @*/
1355: PetscErrorCode PEPSetOptionsPrefix(PEP pep,const char *prefix)1356: {
1361: if (!pep->st) { PEPGetST(pep,&pep->st); }
1362: STSetOptionsPrefix(pep->st,prefix);
1363: if (!pep->V) { PEPGetBV(pep,&pep->V); }
1364: BVSetOptionsPrefix(pep->V,prefix);
1365: if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1366: DSSetOptionsPrefix(pep->ds,prefix);
1367: if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1368: RGSetOptionsPrefix(pep->rg,prefix);
1369: PetscObjectSetOptionsPrefix((PetscObject)pep,prefix);
1370: return(0);
1371: }
1373: /*@C
1374: PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1375: PEP options in the database.
1377: Logically Collective on PEP1379: Input Parameters:
1380: + pep - the polynomial eigensolver context
1381: - prefix - the prefix string to prepend to all PEP option requests
1383: Notes:
1384: A hyphen (-) must NOT be given at the beginning of the prefix name.
1385: The first character of all runtime options is AUTOMATICALLY the hyphen.
1387: Level: advanced
1389: .seealso: PEPSetOptionsPrefix(), PEPGetOptionsPrefix()
1390: @*/
1391: PetscErrorCode PEPAppendOptionsPrefix(PEP pep,const char *prefix)1392: {
1397: if (!pep->st) { PEPGetST(pep,&pep->st); }
1398: STAppendOptionsPrefix(pep->st,prefix);
1399: if (!pep->V) { PEPGetBV(pep,&pep->V); }
1400: BVAppendOptionsPrefix(pep->V,prefix);
1401: if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1402: DSAppendOptionsPrefix(pep->ds,prefix);
1403: if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1404: RGAppendOptionsPrefix(pep->rg,prefix);
1405: PetscObjectAppendOptionsPrefix((PetscObject)pep,prefix);
1406: return(0);
1407: }
1409: /*@C
1410: PEPGetOptionsPrefix - Gets the prefix used for searching for all
1411: PEP options in the database.
1413: Not Collective
1415: Input Parameters:
1416: . pep - the polynomial eigensolver context
1418: Output Parameters:
1419: . prefix - pointer to the prefix string used is returned
1421: Note:
1422: On the Fortran side, the user should pass in a string 'prefix' of
1423: sufficient length to hold the prefix.
1425: Level: advanced
1427: .seealso: PEPSetOptionsPrefix(), PEPAppendOptionsPrefix()
1428: @*/
1429: PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])1430: {
1436: PetscObjectGetOptionsPrefix((PetscObject)pep,prefix);
1437: return(0);
1438: }