LCOV - code coverage report
Current view: top level - pep/interface - pepopts.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 457 477 95.8 %
Date: 2024-11-21 00:34:55 Functions: 30 30 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       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             : */
      14             : 
      15             : #include <slepc/private/pepimpl.h>       /*I "slepcpep.h" I*/
      16             : #include <petscdraw.h>
      17             : 
      18             : /*@C
      19             :    PEPMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
      20             :    indicated by the user.
      21             : 
      22             :    Collective
      23             : 
      24             :    Input Parameters:
      25             : +  pep      - the polynomial eigensolver context
      26             : .  opt      - the command line option for this monitor
      27             : .  name     - the monitor type one is seeking
      28             : .  ctx      - an optional user context for the monitor, or NULL
      29             : -  trackall - whether this monitor tracks all eigenvalues or not
      30             : 
      31             :    Level: developer
      32             : 
      33             : .seealso: PEPMonitorSet(), PEPSetTrackAll()
      34             : @*/
      35         477 : PetscErrorCode PEPMonitorSetFromOptions(PEP pep,const char opt[],const char name[],void *ctx,PetscBool trackall)
      36             : {
      37         477 :   PetscErrorCode       (*mfunc)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
      38         477 :   PetscErrorCode       (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
      39         477 :   PetscErrorCode       (*dfunc)(PetscViewerAndFormat**);
      40         477 :   PetscViewerAndFormat *vf;
      41         477 :   PetscViewer          viewer;
      42         477 :   PetscViewerFormat    format;
      43         477 :   PetscViewerType      vtype;
      44         477 :   char                 key[PETSC_MAX_PATH_LEN];
      45         477 :   PetscBool            flg;
      46             : 
      47         477 :   PetscFunctionBegin;
      48         477 :   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->options,((PetscObject)pep)->prefix,opt,&viewer,&format,&flg));
      49         477 :   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
      50             : 
      51           4 :   PetscCall(PetscViewerGetType(viewer,&vtype));
      52           4 :   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
      53           4 :   PetscCall(PetscFunctionListFind(PEPMonitorList,key,&mfunc));
      54           4 :   PetscCheck(mfunc,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Specified viewer and format not supported");
      55           4 :   PetscCall(PetscFunctionListFind(PEPMonitorCreateList,key,&cfunc));
      56           4 :   PetscCall(PetscFunctionListFind(PEPMonitorDestroyList,key,&dfunc));
      57           4 :   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
      58           4 :   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;
      59             : 
      60           4 :   PetscCall((*cfunc)(viewer,format,ctx,&vf));
      61           4 :   PetscCall(PetscViewerDestroy(&viewer));
      62           4 :   PetscCall(PEPMonitorSet(pep,mfunc,vf,(PetscCtxDestroyFn*)dfunc));
      63           4 :   if (trackall) PetscCall(PEPSetTrackAll(pep,PETSC_TRUE));
      64           4 :   PetscFunctionReturn(PETSC_SUCCESS);
      65             : }
      66             : 
      67             : /*@
      68             :    PEPSetFromOptions - Sets PEP options from the options database.
      69             :    This routine must be called before PEPSetUp() if the user is to be
      70             :    allowed to set the solver type.
      71             : 
      72             :    Collective
      73             : 
      74             :    Input Parameters:
      75             : .  pep - the polynomial eigensolver context
      76             : 
      77             :    Notes:
      78             :    To see all options, run your program with the -help option.
      79             : 
      80             :    Level: beginner
      81             : 
      82             : .seealso: PEPSetOptionsPrefix()
      83             : @*/
      84         159 : PetscErrorCode PEPSetFromOptions(PEP pep)
      85             : {
      86         159 :   char            type[256];
      87         159 :   PetscBool       set,flg,flg1,flg2,flg3,flg4,flg5;
      88         159 :   PetscReal       r,t,array[2]={0,0};
      89         159 :   PetscScalar     s;
      90         159 :   PetscInt        i,j,k;
      91         159 :   PEPScale        scale;
      92         159 :   PEPRefine       refine;
      93         159 :   PEPRefineScheme scheme;
      94             : 
      95         159 :   PetscFunctionBegin;
      96         159 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
      97         159 :   PetscCall(PEPRegisterAll());
      98         477 :   PetscObjectOptionsBegin((PetscObject)pep);
      99         169 :     PetscCall(PetscOptionsFList("-pep_type","Polynomial eigensolver method","PEPSetType",PEPList,(char*)(((PetscObject)pep)->type_name?((PetscObject)pep)->type_name:PEPTOAR),type,sizeof(type),&flg));
     100         159 :     if (flg) PetscCall(PEPSetType(pep,type));
     101          67 :     else if (!((PetscObject)pep)->type_name) PetscCall(PEPSetType(pep,PEPTOAR));
     102             : 
     103         159 :     PetscCall(PetscOptionsBoolGroupBegin("-pep_general","General polynomial eigenvalue problem","PEPSetProblemType",&flg));
     104         159 :     if (flg) PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
     105         159 :     PetscCall(PetscOptionsBoolGroup("-pep_hermitian","Hermitian polynomial eigenvalue problem","PEPSetProblemType",&flg));
     106         159 :     if (flg) PetscCall(PEPSetProblemType(pep,PEP_HERMITIAN));
     107         159 :     PetscCall(PetscOptionsBoolGroup("-pep_hyperbolic","Hyperbolic polynomial eigenvalue problem","PEPSetProblemType",&flg));
     108         159 :     if (flg) PetscCall(PEPSetProblemType(pep,PEP_HYPERBOLIC));
     109         159 :     PetscCall(PetscOptionsBoolGroupEnd("-pep_gyroscopic","Gyroscopic polynomial eigenvalue problem","PEPSetProblemType",&flg));
     110         159 :     if (flg) PetscCall(PEPSetProblemType(pep,PEP_GYROSCOPIC));
     111             : 
     112         159 :     scale = pep->scale;
     113         159 :     PetscCall(PetscOptionsEnum("-pep_scale","Scaling strategy","PEPSetScale",PEPScaleTypes,(PetscEnum)scale,(PetscEnum*)&scale,&flg1));
     114         159 :     r = pep->sfactor;
     115         159 :     PetscCall(PetscOptionsReal("-pep_scale_factor","Scale factor","PEPSetScale",pep->sfactor,&r,&flg2));
     116         159 :     if (!flg2 && r==1.0) r = PETSC_DETERMINE;
     117         159 :     j = pep->sits;
     118         159 :     PetscCall(PetscOptionsInt("-pep_scale_its","Number of iterations in diagonal scaling","PEPSetScale",pep->sits,&j,&flg3));
     119         159 :     t = pep->slambda;
     120         159 :     PetscCall(PetscOptionsReal("-pep_scale_lambda","Estimate of eigenvalue (modulus) for diagonal scaling","PEPSetScale",pep->slambda,&t,&flg4));
     121         159 :     if (flg1 || flg2 || flg3 || flg4) PetscCall(PEPSetScale(pep,scale,r,NULL,NULL,j,t));
     122             : 
     123         159 :     PetscCall(PetscOptionsEnum("-pep_extract","Extraction method","PEPSetExtract",PEPExtractTypes,(PetscEnum)pep->extract,(PetscEnum*)&pep->extract,NULL));
     124             : 
     125         159 :     refine = pep->refine;
     126         159 :     PetscCall(PetscOptionsEnum("-pep_refine","Iterative refinement method","PEPSetRefine",PEPRefineTypes,(PetscEnum)refine,(PetscEnum*)&refine,&flg1));
     127         159 :     i = pep->npart;
     128         159 :     PetscCall(PetscOptionsInt("-pep_refine_partitions","Number of partitions of the communicator for iterative refinement","PEPSetRefine",pep->npart,&i,&flg2));
     129         159 :     r = pep->rtol;
     130         160 :     PetscCall(PetscOptionsReal("-pep_refine_tol","Tolerance for iterative refinement","PEPSetRefine",pep->rtol==(PetscReal)PETSC_DETERMINE?SLEPC_DEFAULT_TOL/1000:pep->rtol,&r,&flg3));
     131         159 :     j = pep->rits;
     132         159 :     PetscCall(PetscOptionsInt("-pep_refine_its","Maximum number of iterations for iterative refinement","PEPSetRefine",pep->rits,&j,&flg4));
     133         159 :     scheme = pep->scheme;
     134         159 :     PetscCall(PetscOptionsEnum("-pep_refine_scheme","Scheme used for linear systems within iterative refinement","PEPSetRefine",PEPRefineSchemes,(PetscEnum)scheme,(PetscEnum*)&scheme,&flg5));
     135         159 :     if (flg1 || flg2 || flg3 || flg4 || flg5) PetscCall(PEPSetRefine(pep,refine,i,r,j,scheme));
     136             : 
     137         159 :     i = pep->max_it;
     138         159 :     PetscCall(PetscOptionsInt("-pep_max_it","Maximum number of iterations","PEPSetTolerances",pep->max_it,&i,&flg1));
     139         159 :     r = pep->tol;
     140         256 :     PetscCall(PetscOptionsReal("-pep_tol","Tolerance","PEPSetTolerances",SlepcDefaultTol(pep->tol),&r,&flg2));
     141         159 :     if (flg1 || flg2) PetscCall(PEPSetTolerances(pep,r,i));
     142             : 
     143         159 :     PetscCall(PetscOptionsBoolGroupBegin("-pep_conv_rel","Relative error convergence test","PEPSetConvergenceTest",&flg));
     144         159 :     if (flg) PetscCall(PEPSetConvergenceTest(pep,PEP_CONV_REL));
     145         159 :     PetscCall(PetscOptionsBoolGroup("-pep_conv_norm","Convergence test relative to the matrix norms","PEPSetConvergenceTest",&flg));
     146         159 :     if (flg) PetscCall(PEPSetConvergenceTest(pep,PEP_CONV_NORM));
     147         159 :     PetscCall(PetscOptionsBoolGroup("-pep_conv_abs","Absolute error convergence test","PEPSetConvergenceTest",&flg));
     148         159 :     if (flg) PetscCall(PEPSetConvergenceTest(pep,PEP_CONV_ABS));
     149         159 :     PetscCall(PetscOptionsBoolGroupEnd("-pep_conv_user","User-defined convergence test","PEPSetConvergenceTest",&flg));
     150         159 :     if (flg) PetscCall(PEPSetConvergenceTest(pep,PEP_CONV_USER));
     151             : 
     152         159 :     PetscCall(PetscOptionsBoolGroupBegin("-pep_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","PEPSetStoppingTest",&flg));
     153         159 :     if (flg) PetscCall(PEPSetStoppingTest(pep,PEP_STOP_BASIC));
     154         159 :     PetscCall(PetscOptionsBoolGroupEnd("-pep_stop_user","User-defined stopping test","PEPSetStoppingTest",&flg));
     155         159 :     if (flg) PetscCall(PEPSetStoppingTest(pep,PEP_STOP_USER));
     156             : 
     157         159 :     i = pep->nev;
     158         159 :     PetscCall(PetscOptionsInt("-pep_nev","Number of eigenvalues to compute","PEPSetDimensions",pep->nev,&i,&flg1));
     159         159 :     j = pep->ncv;
     160         159 :     PetscCall(PetscOptionsInt("-pep_ncv","Number of basis vectors","PEPSetDimensions",pep->ncv,&j,&flg2));
     161         159 :     k = pep->mpd;
     162         159 :     PetscCall(PetscOptionsInt("-pep_mpd","Maximum dimension of projected problem","PEPSetDimensions",pep->mpd,&k,&flg3));
     163         159 :     if (flg1 || flg2 || flg3) PetscCall(PEPSetDimensions(pep,i,j,k));
     164             : 
     165         159 :     PetscCall(PetscOptionsEnum("-pep_basis","Polynomial basis","PEPSetBasis",PEPBasisTypes,(PetscEnum)pep->basis,(PetscEnum*)&pep->basis,NULL));
     166             : 
     167         159 :     PetscCall(PetscOptionsBoolGroupBegin("-pep_largest_magnitude","Compute largest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg));
     168         159 :     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_LARGEST_MAGNITUDE));
     169         159 :     PetscCall(PetscOptionsBoolGroup("-pep_smallest_magnitude","Compute smallest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg));
     170         159 :     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_SMALLEST_MAGNITUDE));
     171         159 :     PetscCall(PetscOptionsBoolGroup("-pep_largest_real","Compute eigenvalues with largest real parts","PEPSetWhichEigenpairs",&flg));
     172         159 :     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_LARGEST_REAL));
     173         159 :     PetscCall(PetscOptionsBoolGroup("-pep_smallest_real","Compute eigenvalues with smallest real parts","PEPSetWhichEigenpairs",&flg));
     174         159 :     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_SMALLEST_REAL));
     175         159 :     PetscCall(PetscOptionsBoolGroup("-pep_largest_imaginary","Compute eigenvalues with largest imaginary parts","PEPSetWhichEigenpairs",&flg));
     176         159 :     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_LARGEST_IMAGINARY));
     177         159 :     PetscCall(PetscOptionsBoolGroup("-pep_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","PEPSetWhichEigenpairs",&flg));
     178         159 :     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_SMALLEST_IMAGINARY));
     179         159 :     PetscCall(PetscOptionsBoolGroup("-pep_target_magnitude","Compute eigenvalues closest to target","PEPSetWhichEigenpairs",&flg));
     180         159 :     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE));
     181         159 :     PetscCall(PetscOptionsBoolGroup("-pep_target_real","Compute eigenvalues with real parts closest to target","PEPSetWhichEigenpairs",&flg));
     182         159 :     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_TARGET_REAL));
     183         159 :     PetscCall(PetscOptionsBoolGroup("-pep_target_imaginary","Compute eigenvalues with imaginary parts closest to target","PEPSetWhichEigenpairs",&flg));
     184         159 :     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_TARGET_IMAGINARY));
     185         159 :     PetscCall(PetscOptionsBoolGroupEnd("-pep_all","Compute all eigenvalues in an interval or a region","PEPSetWhichEigenpairs",&flg));
     186         159 :     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_ALL));
     187             : 
     188         159 :     PetscCall(PetscOptionsScalar("-pep_target","Value of the target","PEPSetTarget",pep->target,&s,&flg));
     189         159 :     if (flg) {
     190          54 :       if (pep->which!=PEP_TARGET_REAL && pep->which!=PEP_TARGET_IMAGINARY) PetscCall(PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE));
     191          54 :       PetscCall(PEPSetTarget(pep,s));
     192             :     }
     193             : 
     194         159 :     k = 2;
     195         159 :     PetscCall(PetscOptionsRealArray("-pep_interval","Computational interval (two real values separated with a comma without spaces)","PEPSetInterval",array,&k,&flg));
     196         159 :     if (flg) {
     197           6 :       PetscCheck(k>1,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_SIZ,"Must pass two values in -pep_interval (comma-separated without spaces)");
     198           6 :       PetscCall(PEPSetWhichEigenpairs(pep,PEP_ALL));
     199           6 :       PetscCall(PEPSetInterval(pep,array[0],array[1]));
     200             :     }
     201             : 
     202             :     /* -----------------------------------------------------------------------*/
     203             :     /*
     204             :       Cancels all monitors hardwired into code before call to PEPSetFromOptions()
     205             :     */
     206         159 :     PetscCall(PetscOptionsBool("-pep_monitor_cancel","Remove any hardwired monitor routines","PEPMonitorCancel",PETSC_FALSE,&flg,&set));
     207         159 :     if (set && flg) PetscCall(PEPMonitorCancel(pep));
     208         159 :     PetscCall(PEPMonitorSetFromOptions(pep,"-pep_monitor","first_approximation",NULL,PETSC_FALSE));
     209         159 :     PetscCall(PEPMonitorSetFromOptions(pep,"-pep_monitor_all","all_approximations",NULL,PETSC_TRUE));
     210         159 :     PetscCall(PEPMonitorSetFromOptions(pep,"-pep_monitor_conv","convergence_history",NULL,PETSC_FALSE));
     211             : 
     212             :     /* -----------------------------------------------------------------------*/
     213         159 :     PetscCall(PetscOptionsName("-pep_view","Print detailed information on solver used","PEPView",&set));
     214         159 :     PetscCall(PetscOptionsName("-pep_view_vectors","View computed eigenvectors","PEPVectorsView",&set));
     215         159 :     PetscCall(PetscOptionsName("-pep_view_values","View computed eigenvalues","PEPValuesView",&set));
     216         159 :     PetscCall(PetscOptionsName("-pep_converged_reason","Print reason for convergence, and number of iterations","PEPConvergedReasonView",&set));
     217         159 :     PetscCall(PetscOptionsName("-pep_error_absolute","Print absolute errors of each eigenpair","PEPErrorView",&set));
     218         159 :     PetscCall(PetscOptionsName("-pep_error_relative","Print relative errors of each eigenpair","PEPErrorView",&set));
     219         159 :     PetscCall(PetscOptionsName("-pep_error_backward","Print backward errors of each eigenpair","PEPErrorView",&set));
     220             : 
     221         159 :     PetscTryTypeMethod(pep,setfromoptions,PetscOptionsObject);
     222         159 :     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)pep,PetscOptionsObject));
     223         159 :   PetscOptionsEnd();
     224             : 
     225         159 :   if (!pep->V) PetscCall(PEPGetBV(pep,&pep->V));
     226         159 :   PetscCall(BVSetFromOptions(pep->V));
     227         159 :   if (!pep->rg) PetscCall(PEPGetRG(pep,&pep->rg));
     228         159 :   PetscCall(RGSetFromOptions(pep->rg));
     229         159 :   if (!pep->ds) PetscCall(PEPGetDS(pep,&pep->ds));
     230         159 :   PetscCall(PEPSetDSType(pep));
     231         159 :   PetscCall(DSSetFromOptions(pep->ds));
     232         159 :   if (!pep->st) PetscCall(PEPGetST(pep,&pep->st));
     233         159 :   PetscCall(PEPSetDefaultST(pep));
     234         159 :   PetscCall(STSetFromOptions(pep->st));
     235         159 :   if (!pep->refineksp) PetscCall(PEPRefineGetKSP(pep,&pep->refineksp));
     236         159 :   PetscCall(KSPSetFromOptions(pep->refineksp));
     237         159 :   PetscFunctionReturn(PETSC_SUCCESS);
     238             : }
     239             : 
     240             : /*@
     241             :    PEPGetTolerances - Gets the tolerance and maximum iteration count used
     242             :    by the PEP convergence tests.
     243             : 
     244             :    Not Collective
     245             : 
     246             :    Input Parameter:
     247             : .  pep - the polynomial eigensolver context
     248             : 
     249             :    Output Parameters:
     250             : +  tol - the convergence tolerance
     251             : -  maxits - maximum number of iterations
     252             : 
     253             :    Notes:
     254             :    The user can specify NULL for any parameter that is not needed.
     255             : 
     256             :    Level: intermediate
     257             : 
     258             : .seealso: PEPSetTolerances()
     259             : @*/
     260          27 : PetscErrorCode PEPGetTolerances(PEP pep,PetscReal *tol,PetscInt *maxits)
     261             : {
     262          27 :   PetscFunctionBegin;
     263          27 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     264          27 :   if (tol)    *tol    = pep->tol;
     265          27 :   if (maxits) *maxits = pep->max_it;
     266          27 :   PetscFunctionReturn(PETSC_SUCCESS);
     267             : }
     268             : 
     269             : /*@
     270             :    PEPSetTolerances - Sets the tolerance and maximum iteration count used
     271             :    by the PEP convergence tests.
     272             : 
     273             :    Logically Collective
     274             : 
     275             :    Input Parameters:
     276             : +  pep - the polynomial eigensolver context
     277             : .  tol - the convergence tolerance
     278             : -  maxits - maximum number of iterations to use
     279             : 
     280             :    Options Database Keys:
     281             : +  -pep_tol <tol> - Sets the convergence tolerance
     282             : -  -pep_max_it <maxits> - Sets the maximum number of iterations allowed
     283             : 
     284             :    Notes:
     285             :    Use PETSC_CURRENT to retain the current value of any of the parameters.
     286             :    Use PETSC_DETERMINE for either argument to assign a default value computed
     287             :    internally (may be different in each solver).
     288             :    For maxits use PETSC_UMLIMITED to indicate there is no upper bound on this value.
     289             : 
     290             :    Level: intermediate
     291             : 
     292             : .seealso: PEPGetTolerances()
     293             : @*/
     294         106 : PetscErrorCode PEPSetTolerances(PEP pep,PetscReal tol,PetscInt maxits)
     295             : {
     296         106 :   PetscFunctionBegin;
     297         106 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     298         318 :   PetscValidLogicalCollectiveReal(pep,tol,2);
     299         318 :   PetscValidLogicalCollectiveInt(pep,maxits,3);
     300         106 :   if (tol == (PetscReal)PETSC_DETERMINE) {
     301           0 :     pep->tol   = PETSC_DETERMINE;
     302           0 :     pep->state = PEP_STATE_INITIAL;
     303         106 :   } else if (tol != (PetscReal)PETSC_CURRENT) {
     304         102 :     PetscCheck(tol>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
     305         102 :     pep->tol = tol;
     306             :   }
     307         106 :   if (maxits == PETSC_DETERMINE) {
     308          11 :     pep->max_it = PETSC_DETERMINE;
     309          11 :     pep->state  = PEP_STATE_INITIAL;
     310          95 :   } else if (maxits == PETSC_UNLIMITED) {
     311           0 :     pep->max_it = PETSC_INT_MAX;
     312          95 :   } else if (maxits != PETSC_CURRENT) {
     313          31 :     PetscCheck(maxits>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
     314          31 :     pep->max_it = maxits;
     315             :   }
     316         106 :   PetscFunctionReturn(PETSC_SUCCESS);
     317             : }
     318             : 
     319             : /*@
     320             :    PEPGetDimensions - Gets the number of eigenvalues to compute
     321             :    and the dimension of the subspace.
     322             : 
     323             :    Not Collective
     324             : 
     325             :    Input Parameter:
     326             : .  pep - the polynomial eigensolver context
     327             : 
     328             :    Output Parameters:
     329             : +  nev - number of eigenvalues to compute
     330             : .  ncv - the maximum dimension of the subspace to be used by the solver
     331             : -  mpd - the maximum dimension allowed for the projected problem
     332             : 
     333             :    Notes:
     334             :    The user can specify NULL for any parameter that is not needed.
     335             : 
     336             :    Level: intermediate
     337             : 
     338             : .seealso: PEPSetDimensions()
     339             : @*/
     340          78 : PetscErrorCode PEPGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
     341             : {
     342          78 :   PetscFunctionBegin;
     343          78 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     344          78 :   if (nev) *nev = pep->nev;
     345          78 :   if (ncv) *ncv = pep->ncv;
     346          78 :   if (mpd) *mpd = pep->mpd;
     347          78 :   PetscFunctionReturn(PETSC_SUCCESS);
     348             : }
     349             : 
     350             : /*@
     351             :    PEPSetDimensions - Sets the number of eigenvalues to compute
     352             :    and the dimension of the subspace.
     353             : 
     354             :    Logically Collective
     355             : 
     356             :    Input Parameters:
     357             : +  pep - the polynomial eigensolver context
     358             : .  nev - number of eigenvalues to compute
     359             : .  ncv - the maximum dimension of the subspace to be used by the solver
     360             : -  mpd - the maximum dimension allowed for the projected problem
     361             : 
     362             :    Options Database Keys:
     363             : +  -pep_nev <nev> - Sets the number of eigenvalues
     364             : .  -pep_ncv <ncv> - Sets the dimension of the subspace
     365             : -  -pep_mpd <mpd> - Sets the maximum projected dimension
     366             : 
     367             :    Notes:
     368             :    Use PETSC_DETERMINE for ncv and mpd to assign a reasonably good value, which is
     369             :    dependent on the solution method. For any of the arguments, use PETSC_CURRENT
     370             :    to preserve the current value.
     371             : 
     372             :    The parameters ncv and mpd are intimately related, so that the user is advised
     373             :    to set one of them at most. Normal usage is that
     374             :    (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
     375             :    (b) in cases where nev is large, the user sets mpd.
     376             : 
     377             :    The value of ncv should always be between nev and (nev+mpd), typically
     378             :    ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
     379             :    a smaller value should be used.
     380             : 
     381             :    When computing all eigenvalues in an interval, see PEPSetInterval(), these
     382             :    parameters lose relevance, and tuning must be done with PEPSTOARSetDimensions().
     383             : 
     384             :    Level: intermediate
     385             : 
     386             : .seealso: PEPGetDimensions(), PEPSetInterval(), PEPSTOARSetDimensions()
     387             : @*/
     388         155 : PetscErrorCode PEPSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
     389             : {
     390         155 :   PetscFunctionBegin;
     391         155 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     392         465 :   PetscValidLogicalCollectiveInt(pep,nev,2);
     393         465 :   PetscValidLogicalCollectiveInt(pep,ncv,3);
     394         465 :   PetscValidLogicalCollectiveInt(pep,mpd,4);
     395         155 :   if (nev != PETSC_CURRENT) {
     396         155 :     PetscCheck(nev>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
     397         155 :     pep->nev = nev;
     398             :   }
     399         155 :   if (ncv == PETSC_DETERMINE) {
     400          76 :     pep->ncv = PETSC_DETERMINE;
     401          79 :   } else if (ncv != PETSC_CURRENT) {
     402          79 :     PetscCheck(ncv>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
     403          79 :     pep->ncv = ncv;
     404             :   }
     405         155 :   if (mpd == PETSC_DETERMINE) {
     406         131 :     pep->mpd = PETSC_DETERMINE;
     407          24 :   } else if (mpd != PETSC_CURRENT) {
     408          24 :     PetscCheck(mpd>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
     409          24 :     pep->mpd = mpd;
     410             :   }
     411         155 :   pep->state = PEP_STATE_INITIAL;
     412         155 :   PetscFunctionReturn(PETSC_SUCCESS);
     413             : }
     414             : 
     415             : /*@
     416             :    PEPSetWhichEigenpairs - Specifies which portion of the spectrum is
     417             :    to be sought.
     418             : 
     419             :    Logically Collective
     420             : 
     421             :    Input Parameters:
     422             : +  pep   - eigensolver context obtained from PEPCreate()
     423             : -  which - the portion of the spectrum to be sought
     424             : 
     425             :    Options Database Keys:
     426             : +   -pep_largest_magnitude - Sets largest eigenvalues in magnitude
     427             : .   -pep_smallest_magnitude - Sets smallest eigenvalues in magnitude
     428             : .   -pep_largest_real - Sets largest real parts
     429             : .   -pep_smallest_real - Sets smallest real parts
     430             : .   -pep_largest_imaginary - Sets largest imaginary parts
     431             : .   -pep_smallest_imaginary - Sets smallest imaginary parts
     432             : .   -pep_target_magnitude - Sets eigenvalues closest to target
     433             : .   -pep_target_real - Sets real parts closest to target
     434             : .   -pep_target_imaginary - Sets imaginary parts closest to target
     435             : -   -pep_all - Sets all eigenvalues in an interval or region
     436             : 
     437             :    Notes:
     438             :    The parameter 'which' can have one of these values
     439             : 
     440             : +     PEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
     441             : .     PEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
     442             : .     PEP_LARGEST_REAL - largest real parts
     443             : .     PEP_SMALLEST_REAL - smallest real parts
     444             : .     PEP_LARGEST_IMAGINARY - largest imaginary parts
     445             : .     PEP_SMALLEST_IMAGINARY - smallest imaginary parts
     446             : .     PEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
     447             : .     PEP_TARGET_REAL - eigenvalues with real part closest to target
     448             : .     PEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
     449             : .     PEP_ALL - all eigenvalues contained in a given interval or region
     450             : -     PEP_WHICH_USER - user defined ordering set with PEPSetEigenvalueComparison()
     451             : 
     452             :    Not all eigensolvers implemented in PEP account for all the possible values
     453             :    stated above. If SLEPc is compiled for real numbers PEP_LARGEST_IMAGINARY
     454             :    and PEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
     455             :    for eigenvalue selection.
     456             : 
     457             :    The target is a scalar value provided with PEPSetTarget().
     458             : 
     459             :    The criterion PEP_TARGET_IMAGINARY is available only in case PETSc and
     460             :    SLEPc have been built with complex scalars.
     461             : 
     462             :    PEP_ALL is intended for use in combination with an interval (see
     463             :    PEPSetInterval()), when all eigenvalues within the interval are requested,
     464             :    and also for computing all eigenvalues in a region with the CISS solver.
     465             :    In both cases, the number of eigenvalues is unknown, so the nev parameter
     466             :    has a different sense, see PEPSetDimensions().
     467             : 
     468             :    Level: intermediate
     469             : 
     470             : .seealso: PEPGetWhichEigenpairs(), PEPSetTarget(), PEPSetInterval(),
     471             :           PEPSetDimensions(), PEPSetEigenvalueComparison(), PEPWhich
     472             : @*/
     473          97 : PetscErrorCode PEPSetWhichEigenpairs(PEP pep,PEPWhich which)
     474             : {
     475          97 :   PetscFunctionBegin;
     476          97 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     477         291 :   PetscValidLogicalCollectiveEnum(pep,which,2);
     478          97 :   switch (which) {
     479          97 :     case PEP_LARGEST_MAGNITUDE:
     480             :     case PEP_SMALLEST_MAGNITUDE:
     481             :     case PEP_LARGEST_REAL:
     482             :     case PEP_SMALLEST_REAL:
     483             :     case PEP_LARGEST_IMAGINARY:
     484             :     case PEP_SMALLEST_IMAGINARY:
     485             :     case PEP_TARGET_MAGNITUDE:
     486             :     case PEP_TARGET_REAL:
     487             : #if defined(PETSC_USE_COMPLEX)
     488             :     case PEP_TARGET_IMAGINARY:
     489             : #endif
     490             :     case PEP_ALL:
     491             :     case PEP_WHICH_USER:
     492          97 :       if (pep->which != which) {
     493          93 :         pep->state = PEP_STATE_INITIAL;
     494          93 :         pep->which = which;
     495             :       }
     496          97 :       break;
     497             : #if !defined(PETSC_USE_COMPLEX)
     498           0 :     case PEP_TARGET_IMAGINARY:
     499           0 :       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"PEP_TARGET_IMAGINARY can be used only with complex scalars");
     500             : #endif
     501           0 :     default:
     502           0 :       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
     503             :   }
     504          97 :   PetscFunctionReturn(PETSC_SUCCESS);
     505             : }
     506             : 
     507             : /*@
     508             :     PEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
     509             :     sought.
     510             : 
     511             :     Not Collective
     512             : 
     513             :     Input Parameter:
     514             : .   pep - eigensolver context obtained from PEPCreate()
     515             : 
     516             :     Output Parameter:
     517             : .   which - the portion of the spectrum to be sought
     518             : 
     519             :     Notes:
     520             :     See PEPSetWhichEigenpairs() for possible values of 'which'.
     521             : 
     522             :     Level: intermediate
     523             : 
     524             : .seealso: PEPSetWhichEigenpairs(), PEPWhich
     525             : @*/
     526           1 : PetscErrorCode PEPGetWhichEigenpairs(PEP pep,PEPWhich *which)
     527             : {
     528           1 :   PetscFunctionBegin;
     529           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     530           1 :   PetscAssertPointer(which,2);
     531           1 :   *which = pep->which;
     532           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     533             : }
     534             : 
     535             : /*@C
     536             :    PEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
     537             :    when PEPSetWhichEigenpairs() is set to PEP_WHICH_USER.
     538             : 
     539             :    Logically Collective
     540             : 
     541             :    Input Parameters:
     542             : +  pep  - eigensolver context obtained from PEPCreate()
     543             : .  comp - a pointer to the comparison function
     544             : -  ctx  - a context pointer (the last parameter to the comparison function)
     545             : 
     546             :    Calling sequence of comp:
     547             : $  PetscErrorCode comp(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
     548             : +   ar     - real part of the 1st eigenvalue
     549             : .   ai     - imaginary part of the 1st eigenvalue
     550             : .   br     - real part of the 2nd eigenvalue
     551             : .   bi     - imaginary part of the 2nd eigenvalue
     552             : .   res    - result of comparison
     553             : -   ctx    - optional context, as set by PEPSetEigenvalueComparison()
     554             : 
     555             :    Note:
     556             :    The returning parameter 'res' can be
     557             : +  negative - if the 1st eigenvalue is preferred to the 2st one
     558             : .  zero     - if both eigenvalues are equally preferred
     559             : -  positive - if the 2st eigenvalue is preferred to the 1st one
     560             : 
     561             :    Level: advanced
     562             : 
     563             : .seealso: PEPSetWhichEigenpairs(), PEPWhich
     564             : @*/
     565           1 : PetscErrorCode PEPSetEigenvalueComparison(PEP pep,PetscErrorCode (*comp)(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx),void* ctx)
     566             : {
     567           1 :   PetscFunctionBegin;
     568           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     569           1 :   pep->sc->comparison    = comp;
     570           1 :   pep->sc->comparisonctx = ctx;
     571           1 :   pep->which             = PEP_WHICH_USER;
     572           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     573             : }
     574             : 
     575             : /*@
     576             :    PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.
     577             : 
     578             :    Logically Collective
     579             : 
     580             :    Input Parameters:
     581             : +  pep  - the polynomial eigensolver context
     582             : -  type - a known type of polynomial eigenvalue problem
     583             : 
     584             :    Options Database Keys:
     585             : +  -pep_general - general problem with no particular structure
     586             : .  -pep_hermitian - problem whose coefficient matrices are Hermitian
     587             : .  -pep_hyperbolic - Hermitian problem that satisfies the definition of hyperbolic
     588             : -  -pep_gyroscopic - problem with Hamiltonian structure
     589             : 
     590             :    Notes:
     591             :    Allowed values for the problem type are general (PEP_GENERAL), Hermitian
     592             :    (PEP_HERMITIAN), hyperbolic (PEP_HYPERBOLIC), and gyroscopic (PEP_GYROSCOPIC).
     593             : 
     594             :    This function is used to instruct SLEPc to exploit certain structure in
     595             :    the polynomial eigenproblem. By default, no particular structure is assumed.
     596             : 
     597             :    If the problem matrices are Hermitian (symmetric in the real case) or
     598             :    Hermitian/skew-Hermitian then the solver can exploit this fact to perform
     599             :    less operations or provide better stability. Hyperbolic problems are a
     600             :    particular case of Hermitian problems, some solvers may treat them simply as
     601             :    Hermitian.
     602             : 
     603             :    Level: intermediate
     604             : 
     605             : .seealso: PEPSetOperators(), PEPSetType(), PEPGetProblemType(), PEPProblemType
     606             : @*/
     607         174 : PetscErrorCode PEPSetProblemType(PEP pep,PEPProblemType type)
     608             : {
     609         174 :   PetscFunctionBegin;
     610         174 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     611         522 :   PetscValidLogicalCollectiveEnum(pep,type,2);
     612         174 :   PetscCheck(type==PEP_GENERAL || type==PEP_HERMITIAN || type==PEP_HYPERBOLIC || type==PEP_GYROSCOPIC,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
     613         174 :   if (type != pep->problem_type) {
     614         172 :     pep->problem_type = type;
     615         172 :     pep->state = PEP_STATE_INITIAL;
     616             :   }
     617         174 :   PetscFunctionReturn(PETSC_SUCCESS);
     618             : }
     619             : 
     620             : /*@
     621             :    PEPGetProblemType - Gets the problem type from the PEP object.
     622             : 
     623             :    Not Collective
     624             : 
     625             :    Input Parameter:
     626             : .  pep - the polynomial eigensolver context
     627             : 
     628             :    Output Parameter:
     629             : .  type - the problem type
     630             : 
     631             :    Level: intermediate
     632             : 
     633             : .seealso: PEPSetProblemType(), PEPProblemType
     634             : @*/
     635           3 : PetscErrorCode PEPGetProblemType(PEP pep,PEPProblemType *type)
     636             : {
     637           3 :   PetscFunctionBegin;
     638           3 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     639           3 :   PetscAssertPointer(type,2);
     640           3 :   *type = pep->problem_type;
     641           3 :   PetscFunctionReturn(PETSC_SUCCESS);
     642             : }
     643             : 
     644             : /*@
     645             :    PEPSetBasis - Specifies the type of polynomial basis used to describe the
     646             :    polynomial eigenvalue problem.
     647             : 
     648             :    Logically Collective
     649             : 
     650             :    Input Parameters:
     651             : +  pep   - the polynomial eigensolver context
     652             : -  basis - the type of polynomial basis
     653             : 
     654             :    Options Database Key:
     655             : .  -pep_basis <basis> - Select the basis type
     656             : 
     657             :    Notes:
     658             :    By default, the coefficient matrices passed via PEPSetOperators() are
     659             :    expressed in the monomial basis, i.e.
     660             :    P(lambda) = A_0 + lambda*A_1 + lambda^2*A_2 + ... + lambda^d*A_d.
     661             :    Other polynomial bases may have better numerical behaviour, but the user
     662             :    must then pass the coefficient matrices accordingly.
     663             : 
     664             :    Level: intermediate
     665             : 
     666             : .seealso: PEPSetOperators(), PEPGetBasis(), PEPBasis
     667             : @*/
     668          24 : PetscErrorCode PEPSetBasis(PEP pep,PEPBasis basis)
     669             : {
     670          24 :   PetscFunctionBegin;
     671          24 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     672          72 :   PetscValidLogicalCollectiveEnum(pep,basis,2);
     673          24 :   pep->basis = basis;
     674          24 :   PetscFunctionReturn(PETSC_SUCCESS);
     675             : }
     676             : 
     677             : /*@
     678             :    PEPGetBasis - Gets the type of polynomial basis from the PEP object.
     679             : 
     680             :    Not Collective
     681             : 
     682             :    Input Parameter:
     683             : .  pep - the polynomial eigensolver context
     684             : 
     685             :    Output Parameter:
     686             : .  basis - the polynomial basis
     687             : 
     688             :    Level: intermediate
     689             : 
     690             : .seealso: PEPSetBasis(), PEPBasis
     691             : @*/
     692          14 : PetscErrorCode PEPGetBasis(PEP pep,PEPBasis *basis)
     693             : {
     694          14 :   PetscFunctionBegin;
     695          14 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     696          14 :   PetscAssertPointer(basis,2);
     697          14 :   *basis = pep->basis;
     698          14 :   PetscFunctionReturn(PETSC_SUCCESS);
     699             : }
     700             : 
     701             : /*@
     702             :    PEPSetTrackAll - Specifies if the solver must compute the residual of all
     703             :    approximate eigenpairs or not.
     704             : 
     705             :    Logically Collective
     706             : 
     707             :    Input Parameters:
     708             : +  pep      - the eigensolver context
     709             : -  trackall - whether compute all residuals or not
     710             : 
     711             :    Notes:
     712             :    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
     713             :    the residual for each eigenpair approximation. Computing the residual is
     714             :    usually an expensive operation and solvers commonly compute the associated
     715             :    residual to the first unconverged eigenpair.
     716             : 
     717             :    The option '-pep_monitor_all' automatically activates this option.
     718             : 
     719             :    Level: developer
     720             : 
     721             : .seealso: PEPGetTrackAll()
     722             : @*/
     723          24 : PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)
     724             : {
     725          24 :   PetscFunctionBegin;
     726          24 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     727          72 :   PetscValidLogicalCollectiveBool(pep,trackall,2);
     728          24 :   pep->trackall = trackall;
     729          24 :   PetscFunctionReturn(PETSC_SUCCESS);
     730             : }
     731             : 
     732             : /*@
     733             :    PEPGetTrackAll - Returns the flag indicating whether all residual norms must
     734             :    be computed or not.
     735             : 
     736             :    Not Collective
     737             : 
     738             :    Input Parameter:
     739             : .  pep - the eigensolver context
     740             : 
     741             :    Output Parameter:
     742             : .  trackall - the returned flag
     743             : 
     744             :    Level: developer
     745             : 
     746             : .seealso: PEPSetTrackAll()
     747             : @*/
     748          34 : PetscErrorCode PEPGetTrackAll(PEP pep,PetscBool *trackall)
     749             : {
     750          34 :   PetscFunctionBegin;
     751          34 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     752          34 :   PetscAssertPointer(trackall,2);
     753          34 :   *trackall = pep->trackall;
     754          34 :   PetscFunctionReturn(PETSC_SUCCESS);
     755             : }
     756             : 
     757             : /*@C
     758             :    PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
     759             :    used in the convergence test.
     760             : 
     761             :    Logically Collective
     762             : 
     763             :    Input Parameters:
     764             : +  pep     - eigensolver context obtained from PEPCreate()
     765             : .  conv    - convergence test function, see PEPConvergenceTestFn for the calling sequence
     766             : .  ctx     - context for private data for the convergence routine (may be NULL)
     767             : -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence
     768             : 
     769             :    Note:
     770             :    If the error estimate returned by the convergence test function is less than
     771             :    the tolerance, then the eigenvalue is accepted as converged.
     772             : 
     773             :    Level: advanced
     774             : 
     775             : .seealso: PEPSetConvergenceTest(), PEPSetTolerances()
     776             : @*/
     777           1 : PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PEPConvergenceTestFn *conv,void* ctx,PetscCtxDestroyFn *destroy)
     778             : {
     779           1 :   PetscFunctionBegin;
     780           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     781           1 :   if (pep->convergeddestroy) PetscCall((*pep->convergeddestroy)(&pep->convergedctx));
     782           1 :   pep->convergeduser    = conv;
     783           1 :   pep->convergeddestroy = destroy;
     784           1 :   pep->convergedctx     = ctx;
     785           1 :   if (conv == PEPConvergedRelative) pep->conv = PEP_CONV_REL;
     786           1 :   else if (conv == PEPConvergedNorm) pep->conv = PEP_CONV_NORM;
     787           1 :   else if (conv == PEPConvergedAbsolute) pep->conv = PEP_CONV_ABS;
     788             :   else {
     789           1 :     pep->conv      = PEP_CONV_USER;
     790           1 :     pep->converged = pep->convergeduser;
     791             :   }
     792           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     793             : }
     794             : 
     795             : /*@
     796             :    PEPSetConvergenceTest - Specifies how to compute the error estimate
     797             :    used in the convergence test.
     798             : 
     799             :    Logically Collective
     800             : 
     801             :    Input Parameters:
     802             : +  pep  - eigensolver context obtained from PEPCreate()
     803             : -  conv - the type of convergence test
     804             : 
     805             :    Options Database Keys:
     806             : +  -pep_conv_abs    - Sets the absolute convergence test
     807             : .  -pep_conv_rel    - Sets the convergence test relative to the eigenvalue
     808             : .  -pep_conv_norm   - Sets the convergence test relative to the matrix norms
     809             : -  -pep_conv_user   - Selects the user-defined convergence test
     810             : 
     811             :    Note:
     812             :    The parameter 'conv' can have one of these values
     813             : +     PEP_CONV_ABS    - absolute error ||r||
     814             : .     PEP_CONV_REL    - error relative to the eigenvalue l, ||r||/|l|
     815             : .     PEP_CONV_NORM   - error relative matrix norms, ||r||/sum_i(l^i*||A_i||)
     816             : -     PEP_CONV_USER   - function set by PEPSetConvergenceTestFunction()
     817             : 
     818             :    Level: intermediate
     819             : 
     820             : .seealso: PEPGetConvergenceTest(), PEPSetConvergenceTestFunction(), PEPSetStoppingTest(), PEPConv
     821             : @*/
     822           7 : PetscErrorCode PEPSetConvergenceTest(PEP pep,PEPConv conv)
     823             : {
     824           7 :   PetscFunctionBegin;
     825           7 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     826          21 :   PetscValidLogicalCollectiveEnum(pep,conv,2);
     827           7 :   switch (conv) {
     828           2 :     case PEP_CONV_ABS:  pep->converged = PEPConvergedAbsolute; break;
     829           2 :     case PEP_CONV_REL:  pep->converged = PEPConvergedRelative; break;
     830           3 :     case PEP_CONV_NORM: pep->converged = PEPConvergedNorm; break;
     831           0 :     case PEP_CONV_USER:
     832           0 :       PetscCheck(pep->convergeduser,PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetConvergenceTestFunction() first");
     833           0 :       pep->converged = pep->convergeduser;
     834           0 :       break;
     835           0 :     default:
     836           0 :       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
     837             :   }
     838           7 :   pep->conv = conv;
     839           7 :   PetscFunctionReturn(PETSC_SUCCESS);
     840             : }
     841             : 
     842             : /*@
     843             :    PEPGetConvergenceTest - Gets the method used to compute the error estimate
     844             :    used in the convergence test.
     845             : 
     846             :    Not Collective
     847             : 
     848             :    Input Parameters:
     849             : .  pep   - eigensolver context obtained from PEPCreate()
     850             : 
     851             :    Output Parameters:
     852             : .  conv  - the type of convergence test
     853             : 
     854             :    Level: intermediate
     855             : 
     856             : .seealso: PEPSetConvergenceTest(), PEPConv
     857             : @*/
     858           1 : PetscErrorCode PEPGetConvergenceTest(PEP pep,PEPConv *conv)
     859             : {
     860           1 :   PetscFunctionBegin;
     861           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     862           1 :   PetscAssertPointer(conv,2);
     863           1 :   *conv = pep->conv;
     864           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     865             : }
     866             : 
     867             : /*@C
     868             :    PEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
     869             :    iteration of the eigensolver.
     870             : 
     871             :    Logically Collective
     872             : 
     873             :    Input Parameters:
     874             : +  pep     - eigensolver context obtained from PEPCreate()
     875             : .  stop    - stopping test function, see PEPStoppingTestFn for the calling sequence
     876             : .  ctx     - context for private data for the stopping routine (may be NULL)
     877             : -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence
     878             : 
     879             :    Note:
     880             :    Normal usage is to first call the default routine PEPStoppingBasic() and then
     881             :    set reason to PEP_CONVERGED_USER if some user-defined conditions have been
     882             :    met. To let the eigensolver continue iterating, the result must be left as
     883             :    PEP_CONVERGED_ITERATING.
     884             : 
     885             :    Level: advanced
     886             : 
     887             : .seealso: PEPSetStoppingTest(), PEPStoppingBasic()
     888             : @*/
     889           1 : PetscErrorCode PEPSetStoppingTestFunction(PEP pep,PEPStoppingTestFn *stop,void* ctx,PetscCtxDestroyFn *destroy)
     890             : {
     891           1 :   PetscFunctionBegin;
     892           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     893           1 :   if (pep->stoppingdestroy) PetscCall((*pep->stoppingdestroy)(&pep->stoppingctx));
     894           1 :   pep->stoppinguser    = stop;
     895           1 :   pep->stoppingdestroy = destroy;
     896           1 :   pep->stoppingctx     = ctx;
     897           1 :   if (stop == PEPStoppingBasic) pep->stop = PEP_STOP_BASIC;
     898             :   else {
     899           1 :     pep->stop     = PEP_STOP_USER;
     900           1 :     pep->stopping = pep->stoppinguser;
     901             :   }
     902           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     903             : }
     904             : 
     905             : /*@
     906             :    PEPSetStoppingTest - Specifies how to decide the termination of the outer
     907             :    loop of the eigensolver.
     908             : 
     909             :    Logically Collective
     910             : 
     911             :    Input Parameters:
     912             : +  pep  - eigensolver context obtained from PEPCreate()
     913             : -  stop - the type of stopping test
     914             : 
     915             :    Options Database Keys:
     916             : +  -pep_stop_basic - Sets the default stopping test
     917             : -  -pep_stop_user  - Selects the user-defined stopping test
     918             : 
     919             :    Note:
     920             :    The parameter 'stop' can have one of these values
     921             : +     PEP_STOP_BASIC - default stopping test
     922             : -     PEP_STOP_USER  - function set by PEPSetStoppingTestFunction()
     923             : 
     924             :    Level: advanced
     925             : 
     926             : .seealso: PEPGetStoppingTest(), PEPSetStoppingTestFunction(), PEPSetConvergenceTest(), PEPStop
     927             : @*/
     928           1 : PetscErrorCode PEPSetStoppingTest(PEP pep,PEPStop stop)
     929             : {
     930           1 :   PetscFunctionBegin;
     931           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     932           3 :   PetscValidLogicalCollectiveEnum(pep,stop,2);
     933           1 :   switch (stop) {
     934           1 :     case PEP_STOP_BASIC: pep->stopping = PEPStoppingBasic; break;
     935           0 :     case PEP_STOP_USER:
     936           0 :       PetscCheck(pep->stoppinguser,PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetStoppingTestFunction() first");
     937           0 :       pep->stopping = pep->stoppinguser;
     938           0 :       break;
     939           0 :     default:
     940           0 :       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
     941             :   }
     942           1 :   pep->stop = stop;
     943           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     944             : }
     945             : 
     946             : /*@
     947             :    PEPGetStoppingTest - Gets the method used to decide the termination of the outer
     948             :    loop of the eigensolver.
     949             : 
     950             :    Not Collective
     951             : 
     952             :    Input Parameters:
     953             : .  pep   - eigensolver context obtained from PEPCreate()
     954             : 
     955             :    Output Parameters:
     956             : .  stop  - the type of stopping test
     957             : 
     958             :    Level: advanced
     959             : 
     960             : .seealso: PEPSetStoppingTest(), PEPStop
     961             : @*/
     962           1 : PetscErrorCode PEPGetStoppingTest(PEP pep,PEPStop *stop)
     963             : {
     964           1 :   PetscFunctionBegin;
     965           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     966           1 :   PetscAssertPointer(stop,2);
     967           1 :   *stop = pep->stop;
     968           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     969             : }
     970             : 
     971             : /*@
     972             :    PEPSetScale - Specifies the scaling strategy to be used.
     973             : 
     974             :    Collective
     975             : 
     976             :    Input Parameters:
     977             : +  pep    - the eigensolver context
     978             : .  scale  - scaling strategy
     979             : .  alpha  - the scaling factor used in the scalar strategy
     980             : .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
     981             : .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
     982             : .  its    - number of iterations of the diagonal scaling algorithm
     983             : -  lambda - approximation to wanted eigenvalues (modulus)
     984             : 
     985             :    Options Database Keys:
     986             : +  -pep_scale <type> - scaling type, one of <none,scalar,diagonal,both>
     987             : .  -pep_scale_factor <alpha> - the scaling factor
     988             : .  -pep_scale_its <its> - number of iterations
     989             : -  -pep_scale_lambda <lambda> - approximation to eigenvalues
     990             : 
     991             :    Notes:
     992             :    There are two non-exclusive scaling strategies, scalar and diagonal.
     993             : 
     994             :    In the scalar strategy, scaling is applied to the eigenvalue, that is,
     995             :    mu = lambda/alpha is the new eigenvalue and all matrices are scaled
     996             :    accordingly. After solving the scaled problem, the original lambda is
     997             :    recovered. Parameter 'alpha' must be positive. Use PETSC_DETERMINE to let
     998             :    the solver compute a reasonable scaling factor, and PETSC_CURRENT to
     999             :    retain a previously set value.
    1000             : 
    1001             :    In the diagonal strategy, the solver works implicitly with matrix Dl*A*Dr,
    1002             :    where Dl and Dr are appropriate diagonal matrices. This improves the accuracy
    1003             :    of the computed results in some cases. The user may provide the Dr and Dl
    1004             :    matrices represented as Vec objects storing diagonal elements. If not
    1005             :    provided, these matrices are computed internally. This option requires
    1006             :    that the polynomial coefficient matrices are of MATAIJ type.
    1007             :    The parameter 'its' is the number of iterations performed by the method.
    1008             :    Parameter 'lambda' must be positive. Use PETSC_DETERMINE or set lambda = 1.0
    1009             :    if no information about eigenvalues is available. PETSC_CURRENT can also
    1010             :    be used to leave its and lambda unchanged.
    1011             : 
    1012             :    Level: intermediate
    1013             : 
    1014             : .seealso: PEPGetScale()
    1015             : @*/
    1016          28 : PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)
    1017             : {
    1018          28 :   PetscFunctionBegin;
    1019          28 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
    1020          84 :   PetscValidLogicalCollectiveEnum(pep,scale,2);
    1021          28 :   pep->scale = scale;
    1022          28 :   if (scale==PEP_SCALE_SCALAR || scale==PEP_SCALE_BOTH) {
    1023          36 :     PetscValidLogicalCollectiveReal(pep,alpha,3);
    1024          12 :     if (alpha == (PetscReal)PETSC_DETERMINE) {
    1025          11 :       pep->sfactor = 0.0;
    1026          11 :       pep->sfactor_set = PETSC_FALSE;
    1027           1 :     } else if (alpha != (PetscReal)PETSC_CURRENT) {
    1028           1 :       PetscCheck(alpha>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
    1029           1 :       pep->sfactor = alpha;
    1030           1 :       pep->sfactor_set = PETSC_TRUE;
    1031             :     }
    1032             :   }
    1033          28 :   if (scale==PEP_SCALE_DIAGONAL || scale==PEP_SCALE_BOTH) {
    1034          18 :     if (Dl) {
    1035           1 :       PetscValidHeaderSpecific(Dl,VEC_CLASSID,4);
    1036           1 :       PetscCheckSameComm(pep,1,Dl,4);
    1037           1 :       PetscCall(PetscObjectReference((PetscObject)Dl));
    1038           1 :       PetscCall(VecDestroy(&pep->Dl));
    1039           1 :       pep->Dl = Dl;
    1040             :     }
    1041          18 :     if (Dr) {
    1042           1 :       PetscValidHeaderSpecific(Dr,VEC_CLASSID,5);
    1043           1 :       PetscCheckSameComm(pep,1,Dr,5);
    1044           1 :       PetscCall(PetscObjectReference((PetscObject)Dr));
    1045           1 :       PetscCall(VecDestroy(&pep->Dr));
    1046           1 :       pep->Dr = Dr;
    1047             :     }
    1048          54 :     PetscValidLogicalCollectiveInt(pep,its,6);
    1049          54 :     PetscValidLogicalCollectiveReal(pep,lambda,7);
    1050          18 :     if (its==PETSC_DETERMINE) pep->sits = 5;
    1051          18 :     else if (its!=PETSC_CURRENT) {
    1052          18 :       PetscCheck(its>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be > 0");
    1053          18 :       pep->sits = its;
    1054             :     }
    1055          18 :     if (lambda == (PetscReal)PETSC_DETERMINE) pep->slambda = 1.0;
    1056          18 :     else if (lambda != (PetscReal)PETSC_CURRENT) {
    1057          18 :       PetscCheck(lambda>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of lambda. Must be > 0");
    1058          18 :       pep->slambda = lambda;
    1059             :     }
    1060             :   }
    1061          28 :   pep->state = PEP_STATE_INITIAL;
    1062          28 :   PetscFunctionReturn(PETSC_SUCCESS);
    1063             : }
    1064             : 
    1065             : /*@
    1066             :    PEPGetScale - Gets the scaling strategy used by the PEP object, and the
    1067             :    associated parameters.
    1068             : 
    1069             :    Not Collective
    1070             : 
    1071             :    Input Parameter:
    1072             : .  pep - the eigensolver context
    1073             : 
    1074             :    Output Parameters:
    1075             : +  scale  - scaling strategy
    1076             : .  alpha  - the scaling factor used in the scalar strategy
    1077             : .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
    1078             : .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
    1079             : .  its    - number of iterations of the diagonal scaling algorithm
    1080             : -  lambda - approximation to wanted eigenvalues (modulus)
    1081             : 
    1082             :    Level: intermediate
    1083             : 
    1084             :    Note:
    1085             :    The user can specify NULL for any parameter that is not needed.
    1086             : 
    1087             :    If Dl or Dr were not set by the user, then the ones computed internally are
    1088             :    returned (or a null pointer if called before PEPSetUp).
    1089             : 
    1090             : .seealso: PEPSetScale(), PEPSetUp()
    1091             : @*/
    1092           1 : PetscErrorCode PEPGetScale(PEP pep,PEPScale *scale,PetscReal *alpha,Vec *Dl,Vec *Dr,PetscInt *its,PetscReal *lambda)
    1093             : {
    1094           1 :   PetscFunctionBegin;
    1095           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
    1096           1 :   if (scale)  *scale  = pep->scale;
    1097           1 :   if (alpha)  *alpha  = pep->sfactor;
    1098           1 :   if (Dl)     *Dl     = pep->Dl;
    1099           1 :   if (Dr)     *Dr     = pep->Dr;
    1100           1 :   if (its)    *its    = pep->sits;
    1101           1 :   if (lambda) *lambda = pep->slambda;
    1102           1 :   PetscFunctionReturn(PETSC_SUCCESS);
    1103             : }
    1104             : 
    1105             : /*@
    1106             :    PEPSetExtract - Specifies the extraction strategy to be used.
    1107             : 
    1108             :    Logically Collective
    1109             : 
    1110             :    Input Parameters:
    1111             : +  pep     - the eigensolver context
    1112             : -  extract - extraction strategy
    1113             : 
    1114             :    Options Database Keys:
    1115             : .  -pep_extract <type> - extraction type, one of <none,norm,residual,structured>
    1116             : 
    1117             :    Level: intermediate
    1118             : 
    1119             : .seealso: PEPGetExtract()
    1120             : @*/
    1121           1 : PetscErrorCode PEPSetExtract(PEP pep,PEPExtract extract)
    1122             : {
    1123           1 :   PetscFunctionBegin;
    1124           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
    1125           3 :   PetscValidLogicalCollectiveEnum(pep,extract,2);
    1126           1 :   pep->extract = extract;
    1127           1 :   PetscFunctionReturn(PETSC_SUCCESS);
    1128             : }
    1129             : 
    1130             : /*@
    1131             :    PEPGetExtract - Gets the extraction strategy used by the PEP object.
    1132             : 
    1133             :    Not Collective
    1134             : 
    1135             :    Input Parameter:
    1136             : .  pep - the eigensolver context
    1137             : 
    1138             :    Output Parameter:
    1139             : .  extract - extraction strategy
    1140             : 
    1141             :    Level: intermediate
    1142             : 
    1143             : .seealso: PEPSetExtract(), PEPExtract
    1144             : @*/
    1145           2 : PetscErrorCode PEPGetExtract(PEP pep,PEPExtract *extract)
    1146             : {
    1147           2 :   PetscFunctionBegin;
    1148           2 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
    1149           2 :   PetscAssertPointer(extract,2);
    1150           2 :   *extract = pep->extract;
    1151           2 :   PetscFunctionReturn(PETSC_SUCCESS);
    1152             : }
    1153             : 
    1154             : /*@
    1155             :    PEPSetRefine - Specifies the refinement type (and options) to be used
    1156             :    after the solve.
    1157             : 
    1158             :    Logically Collective
    1159             : 
    1160             :    Input Parameters:
    1161             : +  pep    - the polynomial eigensolver context
    1162             : .  refine - refinement type
    1163             : .  npart  - number of partitions of the communicator
    1164             : .  tol    - the convergence tolerance
    1165             : .  its    - maximum number of refinement iterations
    1166             : -  scheme - which scheme to be used for solving the involved linear systems
    1167             : 
    1168             :    Options Database Keys:
    1169             : +  -pep_refine <type> - refinement type, one of <none,simple,multiple>
    1170             : .  -pep_refine_partitions <n> - the number of partitions
    1171             : .  -pep_refine_tol <tol> - the tolerance
    1172             : .  -pep_refine_its <its> - number of iterations
    1173             : -  -pep_refine_scheme - to set the scheme for the linear solves
    1174             : 
    1175             :    Notes:
    1176             :    By default, iterative refinement is disabled, since it may be very
    1177             :    costly. There are two possible refinement strategies, simple and multiple.
    1178             :    The simple approach performs iterative refinement on each of the
    1179             :    converged eigenpairs individually, whereas the multiple strategy works
    1180             :    with the invariant pair as a whole, refining all eigenpairs simultaneously.
    1181             :    The latter may be required for the case of multiple eigenvalues.
    1182             : 
    1183             :    In some cases, especially when using direct solvers within the
    1184             :    iterative refinement method, it may be helpful for improved scalability
    1185             :    to split the communicator in several partitions. The npart parameter
    1186             :    indicates how many partitions to use (defaults to 1).
    1187             : 
    1188             :    The tol and its parameters specify the stopping criterion. In the simple
    1189             :    method, refinement continues until the residual of each eigenpair is
    1190             :    below the tolerance (tol defaults to the PEP tol, but may be set to a
    1191             :    different value). In contrast, the multiple method simply performs its
    1192             :    refinement iterations (just one by default).
    1193             : 
    1194             :    The scheme argument is used to change the way in which linear systems are
    1195             :    solved. Possible choices are explicit, mixed block elimination (MBE),
    1196             :    and Schur complement.
    1197             : 
    1198             :    Use PETSC_CURRENT to retain the current value of npart, tol or its. Use
    1199             :    PETSC_DETERMINE to assign a default value.
    1200             : 
    1201             :    Level: intermediate
    1202             : 
    1203             : .seealso: PEPGetRefine()
    1204             : @*/
    1205          30 : PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PEPRefineScheme scheme)
    1206             : {
    1207          30 :   PetscMPIInt    size;
    1208             : 
    1209          30 :   PetscFunctionBegin;
    1210          30 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
    1211          90 :   PetscValidLogicalCollectiveEnum(pep,refine,2);
    1212          90 :   PetscValidLogicalCollectiveInt(pep,npart,3);
    1213          90 :   PetscValidLogicalCollectiveReal(pep,tol,4);
    1214          90 :   PetscValidLogicalCollectiveInt(pep,its,5);
    1215          90 :   PetscValidLogicalCollectiveEnum(pep,scheme,6);
    1216          30 :   pep->refine = refine;
    1217          30 :   if (refine) {  /* process parameters only if not REFINE_NONE */
    1218          29 :     if (npart!=pep->npart) {
    1219          12 :       PetscCall(PetscSubcommDestroy(&pep->refinesubc));
    1220          12 :       PetscCall(KSPDestroy(&pep->refineksp));
    1221             :     }
    1222          29 :     if (npart == PETSC_DETERMINE) {
    1223           0 :       pep->npart = 1;
    1224          29 :     } else if (npart != PETSC_CURRENT) {
    1225          29 :       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size));
    1226          29 :       PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
    1227          29 :       pep->npart = npart;
    1228             :     }
    1229          29 :     if (tol == (PetscReal)PETSC_DETERMINE) {
    1230          28 :       pep->rtol = PETSC_DETERMINE;
    1231           1 :     } else if (tol != (PetscReal)PETSC_CURRENT) {
    1232           1 :       PetscCheck(tol>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
    1233           1 :       pep->rtol = tol;
    1234             :     }
    1235          29 :     if (its==PETSC_DETERMINE) {
    1236          28 :       pep->rits = PETSC_DETERMINE;
    1237           1 :     } else if (its != PETSC_CURRENT) {
    1238           1 :       PetscCheck(its>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
    1239           1 :       pep->rits = its;
    1240             :     }
    1241          29 :     pep->scheme = scheme;
    1242             :   }
    1243          30 :   pep->state = PEP_STATE_INITIAL;
    1244          30 :   PetscFunctionReturn(PETSC_SUCCESS);
    1245             : }
    1246             : 
    1247             : /*@
    1248             :    PEPGetRefine - Gets the refinement strategy used by the PEP object, and the
    1249             :    associated parameters.
    1250             : 
    1251             :    Not Collective
    1252             : 
    1253             :    Input Parameter:
    1254             : .  pep - the polynomial eigensolver context
    1255             : 
    1256             :    Output Parameters:
    1257             : +  refine - refinement type
    1258             : .  npart  - number of partitions of the communicator
    1259             : .  tol    - the convergence tolerance
    1260             : .  its    - maximum number of refinement iterations
    1261             : -  scheme - the scheme used for solving linear systems
    1262             : 
    1263             :    Level: intermediate
    1264             : 
    1265             :    Note:
    1266             :    The user can specify NULL for any parameter that is not needed.
    1267             : 
    1268             : .seealso: PEPSetRefine()
    1269             : @*/
    1270           1 : PetscErrorCode PEPGetRefine(PEP pep,PEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,PEPRefineScheme *scheme)
    1271             : {
    1272           1 :   PetscFunctionBegin;
    1273           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
    1274           1 :   if (refine) *refine = pep->refine;
    1275           1 :   if (npart)  *npart  = pep->npart;
    1276           1 :   if (tol)    *tol    = pep->rtol;
    1277           1 :   if (its)    *its    = pep->rits;
    1278           1 :   if (scheme) *scheme = pep->scheme;
    1279           1 :   PetscFunctionReturn(PETSC_SUCCESS);
    1280             : }
    1281             : 
    1282             : /*@
    1283             :    PEPSetOptionsPrefix - Sets the prefix used for searching for all
    1284             :    PEP options in the database.
    1285             : 
    1286             :    Logically Collective
    1287             : 
    1288             :    Input Parameters:
    1289             : +  pep - the polynomial eigensolver context
    1290             : -  prefix - the prefix string to prepend to all PEP option requests
    1291             : 
    1292             :    Notes:
    1293             :    A hyphen (-) must NOT be given at the beginning of the prefix name.
    1294             :    The first character of all runtime options is AUTOMATICALLY the
    1295             :    hyphen.
    1296             : 
    1297             :    For example, to distinguish between the runtime options for two
    1298             :    different PEP contexts, one could call
    1299             : .vb
    1300             :       PEPSetOptionsPrefix(pep1,"qeig1_")
    1301             :       PEPSetOptionsPrefix(pep2,"qeig2_")
    1302             : .ve
    1303             : 
    1304             :    Level: advanced
    1305             : 
    1306             : .seealso: PEPAppendOptionsPrefix(), PEPGetOptionsPrefix()
    1307             : @*/
    1308          24 : PetscErrorCode PEPSetOptionsPrefix(PEP pep,const char *prefix)
    1309             : {
    1310          24 :   PetscFunctionBegin;
    1311          24 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
    1312          24 :   if (!pep->st) PetscCall(PEPGetST(pep,&pep->st));
    1313          24 :   PetscCall(STSetOptionsPrefix(pep->st,prefix));
    1314          24 :   if (!pep->V) PetscCall(PEPGetBV(pep,&pep->V));
    1315          24 :   PetscCall(BVSetOptionsPrefix(pep->V,prefix));
    1316          24 :   if (!pep->ds) PetscCall(PEPGetDS(pep,&pep->ds));
    1317          24 :   PetscCall(DSSetOptionsPrefix(pep->ds,prefix));
    1318          24 :   if (!pep->rg) PetscCall(PEPGetRG(pep,&pep->rg));
    1319          24 :   PetscCall(RGSetOptionsPrefix(pep->rg,prefix));
    1320          24 :   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pep,prefix));
    1321          24 :   PetscFunctionReturn(PETSC_SUCCESS);
    1322             : }
    1323             : 
    1324             : /*@
    1325             :    PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
    1326             :    PEP options in the database.
    1327             : 
    1328             :    Logically Collective
    1329             : 
    1330             :    Input Parameters:
    1331             : +  pep - the polynomial eigensolver context
    1332             : -  prefix - the prefix string to prepend to all PEP option requests
    1333             : 
    1334             :    Notes:
    1335             :    A hyphen (-) must NOT be given at the beginning of the prefix name.
    1336             :    The first character of all runtime options is AUTOMATICALLY the hyphen.
    1337             : 
    1338             :    Level: advanced
    1339             : 
    1340             : .seealso: PEPSetOptionsPrefix(), PEPGetOptionsPrefix()
    1341             : @*/
    1342          24 : PetscErrorCode PEPAppendOptionsPrefix(PEP pep,const char *prefix)
    1343             : {
    1344          24 :   PetscFunctionBegin;
    1345          24 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
    1346          24 :   if (!pep->st) PetscCall(PEPGetST(pep,&pep->st));
    1347          24 :   PetscCall(STAppendOptionsPrefix(pep->st,prefix));
    1348          24 :   if (!pep->V) PetscCall(PEPGetBV(pep,&pep->V));
    1349          24 :   PetscCall(BVAppendOptionsPrefix(pep->V,prefix));
    1350          24 :   if (!pep->ds) PetscCall(PEPGetDS(pep,&pep->ds));
    1351          24 :   PetscCall(DSAppendOptionsPrefix(pep->ds,prefix));
    1352          24 :   if (!pep->rg) PetscCall(PEPGetRG(pep,&pep->rg));
    1353          24 :   PetscCall(RGAppendOptionsPrefix(pep->rg,prefix));
    1354          24 :   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pep,prefix));
    1355          24 :   PetscFunctionReturn(PETSC_SUCCESS);
    1356             : }
    1357             : 
    1358             : /*@
    1359             :    PEPGetOptionsPrefix - Gets the prefix used for searching for all
    1360             :    PEP options in the database.
    1361             : 
    1362             :    Not Collective
    1363             : 
    1364             :    Input Parameters:
    1365             : .  pep - the polynomial eigensolver context
    1366             : 
    1367             :    Output Parameters:
    1368             : .  prefix - pointer to the prefix string used is returned
    1369             : 
    1370             :    Note:
    1371             :    On the Fortran side, the user should pass in a string 'prefix' of
    1372             :    sufficient length to hold the prefix.
    1373             : 
    1374             :    Level: advanced
    1375             : 
    1376             : .seealso: PEPSetOptionsPrefix(), PEPAppendOptionsPrefix()
    1377             : @*/
    1378           1 : PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])
    1379             : {
    1380           1 :   PetscFunctionBegin;
    1381           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
    1382           1 :   PetscAssertPointer(prefix,2);
    1383           1 :   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pep,prefix));
    1384           1 :   PetscFunctionReturn(PETSC_SUCCESS);
    1385             : }

Generated by: LCOV version 1.14