LCOV - code coverage report
Current view: top level - eps/interface - epsopts.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 509 528 96.4 %
Date: 2024-12-18 00:51:33 Functions: 35 35 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             :    EPS routines related to options that can be set via the command-line
      12             :    or procedurally.
      13             : */
      14             : 
      15             : #include <slepc/private/epsimpl.h>   /*I "slepceps.h" I*/
      16             : #include <petscdraw.h>
      17             : 
      18             : /*@C
      19             :    EPSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
      20             :    indicated by the user.
      21             : 
      22             :    Collective
      23             : 
      24             :    Input Parameters:
      25             : +  eps      - the 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: EPSMonitorSet(), EPSSetTrackAll()
      34             : @*/
      35        1713 : PetscErrorCode EPSMonitorSetFromOptions(EPS eps,const char opt[],const char name[],void *ctx,PetscBool trackall)
      36             : {
      37        1713 :   PetscErrorCode       (*mfunc)(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
      38        1713 :   PetscErrorCode       (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
      39        1713 :   PetscErrorCode       (*dfunc)(PetscViewerAndFormat**);
      40        1713 :   PetscViewerAndFormat *vf;
      41        1713 :   PetscViewer          viewer;
      42        1713 :   PetscViewerFormat    format;
      43        1713 :   PetscViewerType      vtype;
      44        1713 :   char                 key[PETSC_MAX_PATH_LEN];
      45        1713 :   PetscBool            flg;
      46             : 
      47        1713 :   PetscFunctionBegin;
      48        1713 :   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)eps),((PetscObject)eps)->options,((PetscObject)eps)->prefix,opt,&viewer,&format,&flg));
      49        1713 :   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
      50             : 
      51           5 :   PetscCall(PetscViewerGetType(viewer,&vtype));
      52           5 :   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
      53           5 :   PetscCall(PetscFunctionListFind(EPSMonitorList,key,&mfunc));
      54           5 :   PetscCheck(mfunc,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Specified viewer and format not supported");
      55           5 :   PetscCall(PetscFunctionListFind(EPSMonitorCreateList,key,&cfunc));
      56           5 :   PetscCall(PetscFunctionListFind(EPSMonitorDestroyList,key,&dfunc));
      57           5 :   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
      58           5 :   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;
      59             : 
      60           5 :   PetscCall((*cfunc)(viewer,format,ctx,&vf));
      61           5 :   PetscCall(PetscViewerDestroy(&viewer));
      62           5 :   PetscCall(EPSMonitorSet(eps,mfunc,vf,(PetscCtxDestroyFn*)dfunc));
      63           5 :   if (trackall) PetscCall(EPSSetTrackAll(eps,PETSC_TRUE));
      64           5 :   PetscFunctionReturn(PETSC_SUCCESS);
      65             : }
      66             : 
      67             : /*@
      68             :    EPSSetFromOptions - Sets EPS options from the options database.
      69             :    This routine must be called before EPSSetUp() if the user is to be
      70             :    allowed to set the solver type.
      71             : 
      72             :    Collective
      73             : 
      74             :    Input Parameters:
      75             : .  eps - the eigensolver context
      76             : 
      77             :    Notes:
      78             :    To see all options, run your program with the -help option.
      79             : 
      80             :    Level: beginner
      81             : 
      82             : .seealso: EPSSetOptionsPrefix()
      83             : @*/
      84         571 : PetscErrorCode EPSSetFromOptions(EPS eps)
      85             : {
      86         571 :   char           type[256];
      87         571 :   PetscBool      set,flg,flg1,flg2,flg3,bval;
      88         571 :   PetscReal      r,array[2]={0,0};
      89         571 :   PetscScalar    s;
      90         571 :   PetscInt       i,j,k;
      91         571 :   EPSBalance     bal;
      92             : 
      93         571 :   PetscFunctionBegin;
      94         571 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
      95         571 :   PetscCall(EPSRegisterAll());
      96        1713 :   PetscObjectOptionsBegin((PetscObject)eps);
      97         627 :     PetscCall(PetscOptionsFList("-eps_type","Eigensolver method","EPSSetType",EPSList,(char*)(((PetscObject)eps)->type_name?((PetscObject)eps)->type_name:EPSKRYLOVSCHUR),type,sizeof(type),&flg));
      98         571 :     if (flg) PetscCall(EPSSetType(eps,type));
      99         319 :     else if (!((PetscObject)eps)->type_name) PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
     100             : 
     101         571 :     PetscCall(PetscOptionsBoolGroupBegin("-eps_hermitian","Hermitian eigenvalue problem","EPSSetProblemType",&flg));
     102         571 :     if (flg) PetscCall(EPSSetProblemType(eps,EPS_HEP));
     103         571 :     PetscCall(PetscOptionsBoolGroup("-eps_gen_hermitian","Generalized Hermitian eigenvalue problem","EPSSetProblemType",&flg));
     104         571 :     if (flg) PetscCall(EPSSetProblemType(eps,EPS_GHEP));
     105         571 :     PetscCall(PetscOptionsBoolGroup("-eps_non_hermitian","Non-Hermitian eigenvalue problem","EPSSetProblemType",&flg));
     106         571 :     if (flg) PetscCall(EPSSetProblemType(eps,EPS_NHEP));
     107         571 :     PetscCall(PetscOptionsBoolGroup("-eps_gen_non_hermitian","Generalized non-Hermitian eigenvalue problem","EPSSetProblemType",&flg));
     108         571 :     if (flg) PetscCall(EPSSetProblemType(eps,EPS_GNHEP));
     109         571 :     PetscCall(PetscOptionsBoolGroup("-eps_pos_gen_non_hermitian","Generalized non-Hermitian eigenvalue problem with positive semi-definite B","EPSSetProblemType",&flg));
     110         571 :     if (flg) PetscCall(EPSSetProblemType(eps,EPS_PGNHEP));
     111         571 :     PetscCall(PetscOptionsBoolGroup("-eps_gen_indefinite","Generalized Hermitian-indefinite eigenvalue problem","EPSSetProblemType",&flg));
     112         571 :     if (flg) PetscCall(EPSSetProblemType(eps,EPS_GHIEP));
     113         571 :     PetscCall(PetscOptionsBoolGroupEnd("-eps_bse","Structured Bethe-Salpeter eigenvalue problem","EPSSetProblemType",&flg));
     114         571 :     if (flg) PetscCall(EPSSetProblemType(eps,EPS_BSE));
     115             : 
     116         571 :     PetscCall(PetscOptionsBoolGroupBegin("-eps_ritz","Rayleigh-Ritz extraction","EPSSetExtraction",&flg));
     117         571 :     if (flg) PetscCall(EPSSetExtraction(eps,EPS_RITZ));
     118         571 :     PetscCall(PetscOptionsBoolGroup("-eps_harmonic","Harmonic Ritz extraction","EPSSetExtraction",&flg));
     119         571 :     if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC));
     120         571 :     PetscCall(PetscOptionsBoolGroup("-eps_harmonic_relative","Relative harmonic Ritz extraction","EPSSetExtraction",&flg));
     121         571 :     if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC_RELATIVE));
     122         571 :     PetscCall(PetscOptionsBoolGroup("-eps_harmonic_right","Right harmonic Ritz extraction","EPSSetExtraction",&flg));
     123         571 :     if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC_RIGHT));
     124         571 :     PetscCall(PetscOptionsBoolGroup("-eps_harmonic_largest","Largest harmonic Ritz extraction","EPSSetExtraction",&flg));
     125         571 :     if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC_LARGEST));
     126         571 :     PetscCall(PetscOptionsBoolGroup("-eps_refined","Refined Ritz extraction","EPSSetExtraction",&flg));
     127         571 :     if (flg) PetscCall(EPSSetExtraction(eps,EPS_REFINED));
     128         571 :     PetscCall(PetscOptionsBoolGroupEnd("-eps_refined_harmonic","Refined harmonic Ritz extraction","EPSSetExtraction",&flg));
     129         571 :     if (flg) PetscCall(EPSSetExtraction(eps,EPS_REFINED_HARMONIC));
     130             : 
     131         571 :     bal = eps->balance;
     132         571 :     PetscCall(PetscOptionsEnum("-eps_balance","Balancing method","EPSSetBalance",EPSBalanceTypes,(PetscEnum)bal,(PetscEnum*)&bal,&flg1));
     133         571 :     j = eps->balance_its;
     134         571 :     PetscCall(PetscOptionsInt("-eps_balance_its","Number of iterations in balancing","EPSSetBalance",eps->balance_its,&j,&flg2));
     135         571 :     r = eps->balance_cutoff;
     136         571 :     PetscCall(PetscOptionsReal("-eps_balance_cutoff","Cutoff value in balancing","EPSSetBalance",eps->balance_cutoff,&r,&flg3));
     137         571 :     if (flg1 || flg2 || flg3) PetscCall(EPSSetBalance(eps,bal,j,r));
     138             : 
     139         571 :     i = eps->max_it;
     140         571 :     PetscCall(PetscOptionsInt("-eps_max_it","Maximum number of iterations","EPSSetTolerances",eps->max_it,&i,&flg1));
     141         571 :     r = eps->tol;
     142         946 :     PetscCall(PetscOptionsReal("-eps_tol","Tolerance","EPSSetTolerances",SlepcDefaultTol(eps->tol),&r,&flg2));
     143         571 :     if (flg1 || flg2) PetscCall(EPSSetTolerances(eps,r,i));
     144             : 
     145         571 :     r = eps->thres;
     146         571 :     PetscCall(PetscOptionsReal("-eps_threshold_absolute","Absolute threshold","EPSSetThreshold",r,&r,&flg));
     147         571 :     if (flg) PetscCall(EPSSetThreshold(eps,r,PETSC_FALSE));
     148         571 :     PetscCall(PetscOptionsReal("-eps_threshold_relative","Relative threshold","EPSSetThreshold",r,&r,&flg));
     149         571 :     if (flg) PetscCall(EPSSetThreshold(eps,r,PETSC_TRUE));
     150             : 
     151         571 :     PetscCall(PetscOptionsBoolGroupBegin("-eps_conv_rel","Relative error convergence test","EPSSetConvergenceTest",&flg));
     152         571 :     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_REL));
     153         571 :     PetscCall(PetscOptionsBoolGroup("-eps_conv_norm","Convergence test relative to the eigenvalue and the matrix norms","EPSSetConvergenceTest",&flg));
     154         571 :     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_NORM));
     155         571 :     PetscCall(PetscOptionsBoolGroup("-eps_conv_abs","Absolute error convergence test","EPSSetConvergenceTest",&flg));
     156         571 :     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_ABS));
     157         571 :     PetscCall(PetscOptionsBoolGroupEnd("-eps_conv_user","User-defined convergence test","EPSSetConvergenceTest",&flg));
     158         571 :     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_USER));
     159             : 
     160         571 :     PetscCall(PetscOptionsBoolGroupBegin("-eps_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","EPSSetStoppingTest",&flg));
     161         571 :     if (flg) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_BASIC));
     162         571 :     PetscCall(PetscOptionsBoolGroup("-eps_stop_threshold","Stop iteration if a converged eigenvalue is below/above the threshold","EPSSetStoppingTest",&flg));
     163         571 :     if (flg) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_THRESHOLD));
     164         571 :     PetscCall(PetscOptionsBoolGroupEnd("-eps_stop_user","User-defined stopping test","EPSSetStoppingTest",&flg));
     165         571 :     if (flg) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_USER));
     166             : 
     167         571 :     i = eps->nev;
     168         571 :     PetscCall(PetscOptionsInt("-eps_nev","Number of eigenvalues to compute","EPSSetDimensions",eps->nev,&i,&flg1));
     169         571 :     if (!flg1) i = PETSC_CURRENT;
     170         571 :     j = eps->ncv;
     171         571 :     PetscCall(PetscOptionsInt("-eps_ncv","Number of basis vectors","EPSSetDimensions",eps->ncv,&j,&flg2));
     172         571 :     k = eps->mpd;
     173         571 :     PetscCall(PetscOptionsInt("-eps_mpd","Maximum dimension of projected problem","EPSSetDimensions",eps->mpd,&k,&flg3));
     174         571 :     if (flg1 || flg2 || flg3) PetscCall(EPSSetDimensions(eps,i,j,k));
     175             : 
     176         571 :     PetscCall(PetscOptionsBoolGroupBegin("-eps_largest_magnitude","Compute largest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg));
     177         571 :     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE));
     178         571 :     PetscCall(PetscOptionsBoolGroup("-eps_smallest_magnitude","Compute smallest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg));
     179         571 :     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE));
     180         571 :     PetscCall(PetscOptionsBoolGroup("-eps_largest_real","Compute eigenvalues with largest real parts","EPSSetWhichEigenpairs",&flg));
     181         571 :     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
     182         571 :     PetscCall(PetscOptionsBoolGroup("-eps_smallest_real","Compute eigenvalues with smallest real parts","EPSSetWhichEigenpairs",&flg));
     183         571 :     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
     184         571 :     PetscCall(PetscOptionsBoolGroup("-eps_largest_imaginary","Compute eigenvalues with largest imaginary parts","EPSSetWhichEigenpairs",&flg));
     185         571 :     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY));
     186         571 :     PetscCall(PetscOptionsBoolGroup("-eps_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","EPSSetWhichEigenpairs",&flg));
     187         571 :     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY));
     188         571 :     PetscCall(PetscOptionsBoolGroup("-eps_target_magnitude","Compute eigenvalues closest to target","EPSSetWhichEigenpairs",&flg));
     189         571 :     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
     190         571 :     PetscCall(PetscOptionsBoolGroup("-eps_target_real","Compute eigenvalues with real parts closest to target","EPSSetWhichEigenpairs",&flg));
     191         571 :     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL));
     192         571 :     PetscCall(PetscOptionsBoolGroup("-eps_target_imaginary","Compute eigenvalues with imaginary parts closest to target","EPSSetWhichEigenpairs",&flg));
     193         571 :     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_IMAGINARY));
     194         571 :     PetscCall(PetscOptionsBoolGroupEnd("-eps_all","Compute all eigenvalues in an interval or a region","EPSSetWhichEigenpairs",&flg));
     195         571 :     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));
     196             : 
     197         571 :     PetscCall(PetscOptionsScalar("-eps_target","Value of the target","EPSSetTarget",eps->target,&s,&flg));
     198         571 :     if (flg) {
     199          37 :       if (eps->which!=EPS_TARGET_REAL && eps->which!=EPS_TARGET_IMAGINARY) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
     200          37 :       PetscCall(EPSSetTarget(eps,s));
     201             :     }
     202             : 
     203         571 :     k = 2;
     204         571 :     PetscCall(PetscOptionsRealArray("-eps_interval","Computational interval (two real values separated with a comma without spaces)","EPSSetInterval",array,&k,&flg));
     205         571 :     if (flg) {
     206          11 :       PetscCheck(k>1,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_SIZ,"Must pass two values in -eps_interval (comma-separated without spaces)");
     207          11 :       PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));
     208          11 :       PetscCall(EPSSetInterval(eps,array[0],array[1]));
     209             :     }
     210             : 
     211         571 :     PetscCall(PetscOptionsBool("-eps_true_residual","Compute true residuals explicitly","EPSSetTrueResidual",eps->trueres,&eps->trueres,NULL));
     212         571 :     PetscCall(PetscOptionsBool("-eps_purify","Postprocess eigenvectors for purification","EPSSetPurify",eps->purify,&bval,&flg));
     213         571 :     if (flg) PetscCall(EPSSetPurify(eps,bval));
     214         571 :     PetscCall(PetscOptionsBool("-eps_two_sided","Use two-sided variant (to compute left eigenvectors)","EPSSetTwoSided",eps->twosided,&bval,&flg));
     215         571 :     if (flg) PetscCall(EPSSetTwoSided(eps,bval));
     216             : 
     217             :     /* -----------------------------------------------------------------------*/
     218             :     /*
     219             :       Cancels all monitors hardwired into code before call to EPSSetFromOptions()
     220             :     */
     221         571 :     PetscCall(PetscOptionsBool("-eps_monitor_cancel","Remove any hardwired monitor routines","EPSMonitorCancel",PETSC_FALSE,&flg,&set));
     222         571 :     if (set && flg) PetscCall(EPSMonitorCancel(eps));
     223         571 :     PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor","first_approximation",NULL,PETSC_FALSE));
     224         571 :     PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor_all","all_approximations",NULL,PETSC_TRUE));
     225         571 :     PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor_conv","convergence_history",NULL,PETSC_FALSE));
     226             : 
     227             :     /* -----------------------------------------------------------------------*/
     228         571 :     PetscCall(PetscOptionsName("-eps_view","Print detailed information on solver used","EPSView",&set));
     229         571 :     PetscCall(PetscOptionsName("-eps_view_vectors","View computed eigenvectors","EPSVectorsView",&set));
     230         571 :     PetscCall(PetscOptionsName("-eps_view_values","View computed eigenvalues","EPSValuesView",&set));
     231         571 :     PetscCall(PetscOptionsName("-eps_converged_reason","Print reason for convergence, and number of iterations","EPSConvergedReasonView",&set));
     232         571 :     PetscCall(PetscOptionsName("-eps_error_absolute","Print absolute errors of each eigenpair","EPSErrorView",&set));
     233         571 :     PetscCall(PetscOptionsName("-eps_error_relative","Print relative errors of each eigenpair","EPSErrorView",&set));
     234         571 :     PetscCall(PetscOptionsName("-eps_error_backward","Print backward errors of each eigenpair","EPSErrorView",&set));
     235             : 
     236         571 :     PetscTryTypeMethod(eps,setfromoptions,PetscOptionsObject);
     237         571 :     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)eps,PetscOptionsObject));
     238         571 :   PetscOptionsEnd();
     239             : 
     240         571 :   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
     241         571 :   PetscCall(BVSetFromOptions(eps->V));
     242         571 :   if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
     243         571 :   PetscCall(RGSetFromOptions(eps->rg));
     244         571 :   if (eps->useds) {
     245         523 :     if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
     246         523 :     PetscCall(EPSSetDSType(eps));
     247         523 :     PetscCall(DSSetFromOptions(eps->ds));
     248             :   }
     249         571 :   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
     250         571 :   PetscCall(EPSSetDefaultST(eps));
     251         571 :   PetscCall(STSetFromOptions(eps->st));
     252         571 :   PetscFunctionReturn(PETSC_SUCCESS);
     253             : }
     254             : 
     255             : /*@
     256             :    EPSGetTolerances - Gets the tolerance and maximum iteration count used
     257             :    by the EPS convergence tests.
     258             : 
     259             :    Not Collective
     260             : 
     261             :    Input Parameter:
     262             : .  eps - the eigensolver context
     263             : 
     264             :    Output Parameters:
     265             : +  tol - the convergence tolerance
     266             : -  maxits - maximum number of iterations
     267             : 
     268             :    Notes:
     269             :    The user can specify NULL for any parameter that is not needed.
     270             : 
     271             :    Level: intermediate
     272             : 
     273             : .seealso: EPSSetTolerances()
     274             : @*/
     275         253 : PetscErrorCode EPSGetTolerances(EPS eps,PetscReal *tol,PetscInt *maxits)
     276             : {
     277         253 :   PetscFunctionBegin;
     278         253 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     279         253 :   if (tol)    *tol    = eps->tol;
     280         253 :   if (maxits) *maxits = eps->max_it;
     281         253 :   PetscFunctionReturn(PETSC_SUCCESS);
     282             : }
     283             : 
     284             : /*@
     285             :    EPSSetTolerances - Sets the tolerance and maximum iteration count used
     286             :    by the EPS convergence tests.
     287             : 
     288             :    Logically Collective
     289             : 
     290             :    Input Parameters:
     291             : +  eps - the eigensolver context
     292             : .  tol - the convergence tolerance
     293             : -  maxits - maximum number of iterations to use
     294             : 
     295             :    Options Database Keys:
     296             : +  -eps_tol <tol> - Sets the convergence tolerance
     297             : -  -eps_max_it <maxits> - Sets the maximum number of iterations allowed
     298             : 
     299             :    Notes:
     300             :    Use PETSC_CURRENT to retain the current value of any of the parameters.
     301             :    Use PETSC_DETERMINE for either argument to assign a default value computed
     302             :    internally (may be different in each solver).
     303             :    For maxits use PETSC_UMLIMITED to indicate there is no upper bound on this value.
     304             : 
     305             :    Level: intermediate
     306             : 
     307             : .seealso: EPSGetTolerances()
     308             : @*/
     309         531 : PetscErrorCode EPSSetTolerances(EPS eps,PetscReal tol,PetscInt maxits)
     310             : {
     311         531 :   PetscFunctionBegin;
     312         531 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     313        1593 :   PetscValidLogicalCollectiveReal(eps,tol,2);
     314        1593 :   PetscValidLogicalCollectiveInt(eps,maxits,3);
     315         531 :   if (tol == (PetscReal)PETSC_DETERMINE) {
     316          19 :     eps->tol   = PETSC_DETERMINE;
     317          19 :     eps->state = EPS_STATE_INITIAL;
     318         512 :   } else if (tol != (PetscReal)PETSC_CURRENT) {
     319         512 :     PetscCheck(tol>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
     320         512 :     eps->tol = tol;
     321             :   }
     322         531 :   if (maxits == PETSC_DETERMINE) {
     323         132 :     eps->max_it = PETSC_DETERMINE;
     324         132 :     eps->state  = EPS_STATE_INITIAL;
     325         399 :   } else if (maxits == PETSC_UNLIMITED) {
     326           0 :     eps->max_it = PETSC_INT_MAX;
     327         399 :   } else if (maxits != PETSC_CURRENT) {
     328         202 :     PetscCheck(maxits>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
     329         202 :     eps->max_it = maxits;
     330             :   }
     331         531 :   PetscFunctionReturn(PETSC_SUCCESS);
     332             : }
     333             : 
     334             : /*@
     335             :    EPSGetDimensions - Gets the number of eigenvalues to compute
     336             :    and the dimension of the subspace.
     337             : 
     338             :    Not Collective
     339             : 
     340             :    Input Parameter:
     341             : .  eps - the eigensolver context
     342             : 
     343             :    Output Parameters:
     344             : +  nev - number of eigenvalues to compute
     345             : .  ncv - the maximum dimension of the subspace to be used by the solver
     346             : -  mpd - the maximum dimension allowed for the projected problem
     347             : 
     348             :    Level: intermediate
     349             : 
     350             : .seealso: EPSSetDimensions()
     351             : @*/
     352         388 : PetscErrorCode EPSGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
     353             : {
     354         388 :   PetscFunctionBegin;
     355         388 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     356         594 :   if (nev) *nev = eps->nev? eps->nev: 1;
     357         388 :   if (ncv) *ncv = eps->ncv;
     358         388 :   if (mpd) *mpd = eps->mpd;
     359         388 :   PetscFunctionReturn(PETSC_SUCCESS);
     360             : }
     361             : 
     362             : /*@
     363             :    EPSSetDimensions - Sets the number of eigenvalues to compute
     364             :    and the dimension of the subspace.
     365             : 
     366             :    Logically Collective
     367             : 
     368             :    Input Parameters:
     369             : +  eps - the eigensolver context
     370             : .  nev - number of eigenvalues to compute
     371             : .  ncv - the maximum dimension of the subspace to be used by the solver
     372             : -  mpd - the maximum dimension allowed for the projected problem
     373             : 
     374             :    Options Database Keys:
     375             : +  -eps_nev <nev> - Sets the number of eigenvalues
     376             : .  -eps_ncv <ncv> - Sets the dimension of the subspace
     377             : -  -eps_mpd <mpd> - Sets the maximum projected dimension
     378             : 
     379             :    Notes:
     380             :    Use PETSC_DETERMINE for ncv and mpd to assign a reasonably good value, which is
     381             :    dependent on the solution method. For any of the arguments, use PETSC_CURRENT
     382             :    to preserve the current value.
     383             : 
     384             :    The parameters ncv and mpd are intimately related, so that the user is advised
     385             :    to set one of them at most. Normal usage is that
     386             :    (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
     387             :    (b) in cases where nev is large, the user sets mpd.
     388             : 
     389             :    The value of ncv should always be between nev and (nev+mpd), typically
     390             :    ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
     391             :    a smaller value should be used.
     392             : 
     393             :    When computing all eigenvalues in an interval, see EPSSetInterval(), these
     394             :    parameters lose relevance, and tuning must be done with
     395             :    EPSKrylovSchurSetDimensions().
     396             : 
     397             :    Level: intermediate
     398             : 
     399             : .seealso: EPSGetDimensions(), EPSSetInterval(), EPSKrylovSchurSetDimensions()
     400             : @*/
     401         538 : PetscErrorCode EPSSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
     402             : {
     403         538 :   PetscFunctionBegin;
     404         538 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     405        1614 :   PetscValidLogicalCollectiveInt(eps,nev,2);
     406        1614 :   PetscValidLogicalCollectiveInt(eps,ncv,3);
     407        1614 :   PetscValidLogicalCollectiveInt(eps,mpd,4);
     408         538 :   if (nev != PETSC_CURRENT) {
     409         518 :     PetscCheck(nev>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
     410         518 :     eps->nev = nev;
     411             :   }
     412         538 :   if (ncv == PETSC_DETERMINE) {
     413         358 :     eps->ncv = PETSC_DETERMINE;
     414         180 :   } else if (ncv != PETSC_CURRENT) {
     415         177 :     PetscCheck(ncv>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
     416         177 :     eps->ncv = ncv;
     417             :   }
     418         538 :   if (mpd == PETSC_DETERMINE) {
     419         484 :     eps->mpd = PETSC_DETERMINE;
     420          54 :   } else if (mpd != PETSC_CURRENT) {
     421          51 :     PetscCheck(mpd>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
     422          51 :     eps->mpd = mpd;
     423             :   }
     424         538 :   eps->state = EPS_STATE_INITIAL;
     425         538 :   PetscFunctionReturn(PETSC_SUCCESS);
     426             : }
     427             : 
     428             : /*@
     429             :    EPSSetWhichEigenpairs - Specifies which portion of the spectrum is
     430             :    to be sought.
     431             : 
     432             :    Logically Collective
     433             : 
     434             :    Input Parameters:
     435             : +  eps   - eigensolver context obtained from EPSCreate()
     436             : -  which - the portion of the spectrum to be sought
     437             : 
     438             :    Options Database Keys:
     439             : +   -eps_largest_magnitude - Sets largest eigenvalues in magnitude
     440             : .   -eps_smallest_magnitude - Sets smallest eigenvalues in magnitude
     441             : .   -eps_largest_real - Sets largest real parts
     442             : .   -eps_smallest_real - Sets smallest real parts
     443             : .   -eps_largest_imaginary - Sets largest imaginary parts
     444             : .   -eps_smallest_imaginary - Sets smallest imaginary parts
     445             : .   -eps_target_magnitude - Sets eigenvalues closest to target
     446             : .   -eps_target_real - Sets real parts closest to target
     447             : .   -eps_target_imaginary - Sets imaginary parts closest to target
     448             : -   -eps_all - Sets all eigenvalues in an interval or region
     449             : 
     450             :    Notes:
     451             :    The parameter 'which' can have one of these values
     452             : 
     453             : +     EPS_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
     454             : .     EPS_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
     455             : .     EPS_LARGEST_REAL - largest real parts
     456             : .     EPS_SMALLEST_REAL - smallest real parts
     457             : .     EPS_LARGEST_IMAGINARY - largest imaginary parts
     458             : .     EPS_SMALLEST_IMAGINARY - smallest imaginary parts
     459             : .     EPS_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
     460             : .     EPS_TARGET_REAL - eigenvalues with real part closest to target
     461             : .     EPS_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
     462             : .     EPS_ALL - all eigenvalues contained in a given interval or region
     463             : -     EPS_WHICH_USER - user defined ordering set with EPSSetEigenvalueComparison()
     464             : 
     465             :    Not all eigensolvers implemented in EPS account for all the possible values
     466             :    stated above. Also, some values make sense only for certain types of
     467             :    problems. If SLEPc is compiled for real numbers EPS_LARGEST_IMAGINARY
     468             :    and EPS_SMALLEST_IMAGINARY use the absolute value of the imaginary part
     469             :    for eigenvalue selection.
     470             : 
     471             :    The target is a scalar value provided with EPSSetTarget().
     472             : 
     473             :    The criterion EPS_TARGET_IMAGINARY is available only in case PETSc and
     474             :    SLEPc have been built with complex scalars.
     475             : 
     476             :    EPS_ALL is intended for use in combination with an interval (see
     477             :    EPSSetInterval()), when all eigenvalues within the interval are requested,
     478             :    or in the context of the CISS solver for computing all eigenvalues in a region.
     479             :    In those cases, the number of eigenvalues is unknown, so the nev parameter
     480             :    has a different sense, see EPSSetDimensions().
     481             : 
     482             :    Level: intermediate
     483             : 
     484             : .seealso: EPSGetWhichEigenpairs(), EPSSetTarget(), EPSSetInterval(),
     485             :           EPSSetDimensions(), EPSSetEigenvalueComparison(), EPSWhich
     486             : @*/
     487         587 : PetscErrorCode EPSSetWhichEigenpairs(EPS eps,EPSWhich which)
     488             : {
     489         587 :   PetscFunctionBegin;
     490         587 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     491        1761 :   PetscValidLogicalCollectiveEnum(eps,which,2);
     492         587 :   switch (which) {
     493         587 :     case EPS_LARGEST_MAGNITUDE:
     494             :     case EPS_SMALLEST_MAGNITUDE:
     495             :     case EPS_LARGEST_REAL:
     496             :     case EPS_SMALLEST_REAL:
     497             :     case EPS_LARGEST_IMAGINARY:
     498             :     case EPS_SMALLEST_IMAGINARY:
     499             :     case EPS_TARGET_MAGNITUDE:
     500             :     case EPS_TARGET_REAL:
     501             : #if defined(PETSC_USE_COMPLEX)
     502             :     case EPS_TARGET_IMAGINARY:
     503             : #endif
     504             :     case EPS_ALL:
     505             :     case EPS_WHICH_USER:
     506         587 :       if (eps->which != which) {
     507         500 :         eps->state = EPS_STATE_INITIAL;
     508         500 :         eps->which = which;
     509             :       }
     510         587 :       break;
     511             : #if !defined(PETSC_USE_COMPLEX)
     512             :     case EPS_TARGET_IMAGINARY:
     513             :       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"EPS_TARGET_IMAGINARY can be used only with complex scalars");
     514             : #endif
     515           0 :     default:
     516           0 :       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
     517             :   }
     518         587 :   PetscFunctionReturn(PETSC_SUCCESS);
     519             : }
     520             : 
     521             : /*@
     522             :    EPSGetWhichEigenpairs - Returns which portion of the spectrum is to be
     523             :    sought.
     524             : 
     525             :    Not Collective
     526             : 
     527             :    Input Parameter:
     528             : .  eps - eigensolver context obtained from EPSCreate()
     529             : 
     530             :    Output Parameter:
     531             : .  which - the portion of the spectrum to be sought
     532             : 
     533             :    Notes:
     534             :    See EPSSetWhichEigenpairs() for possible values of 'which'.
     535             : 
     536             :    Level: intermediate
     537             : 
     538             : .seealso: EPSSetWhichEigenpairs(), EPSWhich
     539             : @*/
     540           1 : PetscErrorCode EPSGetWhichEigenpairs(EPS eps,EPSWhich *which)
     541             : {
     542           1 :   PetscFunctionBegin;
     543           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     544           1 :   PetscAssertPointer(which,2);
     545           1 :   *which = eps->which;
     546           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     547             : }
     548             : 
     549             : /*@
     550             :    EPSSetThreshold - Sets the threshold used in the threshold stopping test.
     551             : 
     552             :    Logically Collective
     553             : 
     554             :    Input Parameters:
     555             : +  eps   - the eigenvalue solver context
     556             : .  thres - the threshold value
     557             : -  rel   - whether the threshold is relative or not
     558             : 
     559             :    Options Database Keys:
     560             : +  -eps_threshold_absolute <thres> - Sets an absolute threshold
     561             : -  -eps_threshold_relative <thres> - Sets a relative threshold
     562             : 
     563             :    Notes:
     564             :    This function internally calls EPSSetStoppingTest() to set a special stopping
     565             :    test based on the threshold, where eigenvalues are computed in sequence until
     566             :    one of the computed eigenvalues is below the threshold (in magnitude). This is
     567             :    the interpretation in case of searching for largest eigenvalues in magnitude,
     568             :    see EPSSetWhichEigenpairs().
     569             : 
     570             :    If the solver is configured to compute smallest magnitude eigenvalues, then the
     571             :    threshold must be interpreted in the opposite direction, i.e., the computation
     572             :    will stop when one of the computed values is above the threshold (in magnitude).
     573             : 
     574             :    The threshold can also be used when computing largest/smallest real eigenvalues
     575             :    (i.e, rightmost or leftmost), in which case the threshold is allowed to be
     576             :    negative. The solver will stop when one of the computed eigenvalues is above
     577             :    or below the threshold (considering the real part of the eigenvalue). This mode
     578             :    is allowed only in problem types whose eigenvalues are always real (e.g., HEP).
     579             : 
     580             :    In the case of largest magnitude eigenvalues, the threshold can be made relative
     581             :    with respect to the dominant eigenvalue. Otherwise, the argument rel should be
     582             :    PETSC_FALSE.
     583             : 
     584             :    An additional use case is with target magnitude selection of eigenvalues (e.g.,
     585             :    with shift-and-invert), but this must be used with caution to avoid unexpected
     586             :    behaviour. With an absolute threshold, the solver will assume that leftmost
     587             :    eigenvalues are being computed (e.g., with target=0 for a problem with real
     588             :    positive eigenvalues). In case of a relative threshold, a value of threshold<1
     589             :    implies that the wanted eigenvalues are the largest ones, and otherwise the
     590             :    solver assumes that smallest eigenvalues are being computed.
     591             : 
     592             :    The test against the threshold is done for converged eigenvalues, which
     593             :    implies that the final number of converged eigenvalues will be at least
     594             :    one more than the actual number of values below/above the threshold.
     595             : 
     596             :    Since the number of computed eigenvalues is not known a priori, the solver
     597             :    will need to reallocate the basis of vectors internally, to have enough room
     598             :    to accommodate all the eigenvectors. Hence, this option must be used with
     599             :    caution to avoid out-of-memory problems. The recommendation is to set the value
     600             :    of ncv to be larger than the estimated number of eigenvalues, to minimize the
     601             :    number of reallocations.
     602             : 
     603             :    If a number of wanted eigenvalues has been set with EPSSetDimensions()
     604             :    it is also taken into account and the solver will stop when one of the two
     605             :    conditions (threshold or number of converged values) is met.
     606             : 
     607             :    Use EPSSetStoppingTest() to return to the usual computation of a fixed number
     608             :    of eigenvalues.
     609             : 
     610             :    Level: advanced
     611             : 
     612             : .seealso: EPSGetThreshold(), EPSSetStoppingTest(), EPSSetDimensions(), EPSSetWhichEigenpairs(), EPSSetProblemType()
     613             : @*/
     614          15 : PetscErrorCode EPSSetThreshold(EPS eps,PetscReal thres,PetscBool rel)
     615             : {
     616          15 :   PetscFunctionBegin;
     617          15 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     618          45 :   PetscValidLogicalCollectiveReal(eps,thres,2);
     619          45 :   PetscValidLogicalCollectiveBool(eps,rel,3);
     620          15 :   if (eps->thres != thres || eps->threlative != rel) {
     621          15 :     eps->thres = thres;
     622          15 :     eps->threlative = rel;
     623          15 :     eps->state = EPS_STATE_INITIAL;
     624          15 :     PetscCall(EPSSetStoppingTest(eps,EPS_STOP_THRESHOLD));
     625             :   }
     626          15 :   PetscFunctionReturn(PETSC_SUCCESS);
     627             : }
     628             : 
     629             : /*@
     630             :    EPSGetThreshold - Gets the threshold used by the threshold stopping test.
     631             : 
     632             :    Not Collective
     633             : 
     634             :    Input Parameter:
     635             : .  eps - the eigenvalue solver context
     636             : 
     637             :    Output Parameters:
     638             : +  thres - the threshold
     639             : -  rel   - whether the threshold is relative or not
     640             : 
     641             :    Level: advanced
     642             : 
     643             : .seealso: EPSSetThreshold()
     644             : @*/
     645           5 : PetscErrorCode EPSGetThreshold(EPS eps,PetscReal *thres,PetscBool *rel)
     646             : {
     647           5 :   PetscFunctionBegin;
     648           5 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     649           5 :   if (thres) *thres = eps->thres;
     650           5 :   if (rel)   *rel   = eps->threlative;
     651           5 :   PetscFunctionReturn(PETSC_SUCCESS);
     652             : }
     653             : 
     654             : /*@C
     655             :    EPSSetEigenvalueComparison - Specifies the eigenvalue comparison function
     656             :    when EPSSetWhichEigenpairs() is set to EPS_WHICH_USER.
     657             : 
     658             :    Logically Collective
     659             : 
     660             :    Input Parameters:
     661             : +  eps  - eigensolver context obtained from EPSCreate()
     662             : .  func - the comparison function, see EPSEigenvalueComparisonFn for the calling sequence
     663             : -  ctx  - a context pointer (the last parameter to the comparison function)
     664             : 
     665             :    Note:
     666             :    The returning parameter 'res' can be
     667             : +  negative - if the 1st eigenvalue is preferred to the 2st one
     668             : .  zero     - if both eigenvalues are equally preferred
     669             : -  positive - if the 2st eigenvalue is preferred to the 1st one
     670             : 
     671             :    Level: advanced
     672             : 
     673             : .seealso: EPSSetWhichEigenpairs(), EPSWhich
     674             : @*/
     675          37 : PetscErrorCode EPSSetEigenvalueComparison(EPS eps,SlepcEigenvalueComparisonFn *func,void* ctx)
     676             : {
     677          37 :   PetscFunctionBegin;
     678          37 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     679          37 :   eps->sc->comparison    = func;
     680          37 :   eps->sc->comparisonctx = ctx;
     681          37 :   eps->which             = EPS_WHICH_USER;
     682          37 :   PetscFunctionReturn(PETSC_SUCCESS);
     683             : }
     684             : 
     685             : /*@C
     686             :    EPSSetArbitrarySelection - Specifies a function intended to look for
     687             :    eigenvalues according to an arbitrary selection criterion. This criterion
     688             :    can be based on a computation involving the current eigenvector approximation.
     689             : 
     690             :    Logically Collective
     691             : 
     692             :    Input Parameters:
     693             : +  eps  - eigensolver context obtained from EPSCreate()
     694             : .  func - the arbitrary selection function, see SlepcArbitrarySelectionFn for a calling sequence
     695             : -  ctx  - a context pointer (the last parameter to the arbitrary selection function)
     696             : 
     697             :    Notes:
     698             :    This provides a mechanism to select eigenpairs by evaluating a user-defined
     699             :    function. When a function has been provided, the default selection based on
     700             :    sorting the eigenvalues is replaced by the sorting of the results of this
     701             :    function (with the same sorting criterion given in EPSSetWhichEigenpairs()).
     702             : 
     703             :    For instance, suppose you want to compute those eigenvectors that maximize
     704             :    a certain computable expression. Then implement the computation using
     705             :    the arguments xr and xi, and return the result in rr. Then set the standard
     706             :    sorting by magnitude so that the eigenpair with largest value of rr is
     707             :    selected.
     708             : 
     709             :    This evaluation function is collective, that is, all processes call it and
     710             :    it can use collective operations; furthermore, the computed result must
     711             :    be the same in all processes.
     712             : 
     713             :    The result of func is expressed as a complex number so that it is possible to
     714             :    use the standard eigenvalue sorting functions, but normally only rr is used.
     715             :    Set ri to zero unless it is meaningful in your application.
     716             : 
     717             :    Level: advanced
     718             : 
     719             : .seealso: EPSSetWhichEigenpairs()
     720             : @*/
     721           8 : PetscErrorCode EPSSetArbitrarySelection(EPS eps,SlepcArbitrarySelectionFn *func,void* ctx)
     722             : {
     723           8 :   PetscFunctionBegin;
     724           8 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     725           8 :   eps->arbitrary    = func;
     726           8 :   eps->arbitraryctx = ctx;
     727           8 :   eps->state        = EPS_STATE_INITIAL;
     728           8 :   PetscFunctionReturn(PETSC_SUCCESS);
     729             : }
     730             : 
     731             : /*@C
     732             :    EPSSetConvergenceTestFunction - Sets a function to compute the error estimate
     733             :    used in the convergence test.
     734             : 
     735             :    Logically Collective
     736             : 
     737             :    Input Parameters:
     738             : +  eps     - eigensolver context obtained from EPSCreate()
     739             : .  func    - convergence test function, see EPSConvergenceTestFn for the calling sequence
     740             : .  ctx     - context for private data for the convergence routine (may be NULL)
     741             : -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence
     742             : 
     743             :    Note:
     744             :    If the error estimate returned by the convergence test function is less than
     745             :    the tolerance, then the eigenvalue is accepted as converged.
     746             : 
     747             :    Level: advanced
     748             : 
     749             : .seealso: EPSSetConvergenceTest(), EPSSetTolerances()
     750             : @*/
     751          33 : PetscErrorCode EPSSetConvergenceTestFunction(EPS eps,EPSConvergenceTestFn *func,void* ctx,PetscCtxDestroyFn *destroy)
     752             : {
     753          33 :   PetscFunctionBegin;
     754          33 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     755          33 :   if (eps->convergeddestroy) PetscCall((*eps->convergeddestroy)(&eps->convergedctx));
     756          33 :   eps->convergeduser    = func;
     757          33 :   eps->convergeddestroy = destroy;
     758          33 :   eps->convergedctx     = ctx;
     759          33 :   if (func == EPSConvergedRelative) eps->conv = EPS_CONV_REL;
     760          33 :   else if (func == EPSConvergedNorm) eps->conv = EPS_CONV_NORM;
     761          33 :   else if (func == EPSConvergedAbsolute) eps->conv = EPS_CONV_ABS;
     762             :   else {
     763          33 :     eps->conv      = EPS_CONV_USER;
     764          33 :     eps->converged = eps->convergeduser;
     765             :   }
     766          33 :   PetscFunctionReturn(PETSC_SUCCESS);
     767             : }
     768             : 
     769             : /*@
     770             :    EPSSetConvergenceTest - Specifies how to compute the error estimate
     771             :    used in the convergence test.
     772             : 
     773             :    Logically Collective
     774             : 
     775             :    Input Parameters:
     776             : +  eps  - eigensolver context obtained from EPSCreate()
     777             : -  conv - the type of convergence test
     778             : 
     779             :    Options Database Keys:
     780             : +  -eps_conv_abs  - Sets the absolute convergence test
     781             : .  -eps_conv_rel  - Sets the convergence test relative to the eigenvalue
     782             : .  -eps_conv_norm - Sets the convergence test relative to the matrix norms
     783             : -  -eps_conv_user - Selects the user-defined convergence test
     784             : 
     785             :    Note:
     786             :    The parameter 'conv' can have one of these values
     787             : +     EPS_CONV_ABS  - absolute error ||r||
     788             : .     EPS_CONV_REL  - error relative to the eigenvalue l, ||r||/|l|
     789             : .     EPS_CONV_NORM - error relative to the matrix norms, ||r||/(||A||+|l|*||B||)
     790             : -     EPS_CONV_USER - function set by EPSSetConvergenceTestFunction()
     791             : 
     792             :    Level: intermediate
     793             : 
     794             : .seealso: EPSGetConvergenceTest(), EPSSetConvergenceTestFunction(), EPSSetStoppingTest(), EPSConv
     795             : @*/
     796         154 : PetscErrorCode EPSSetConvergenceTest(EPS eps,EPSConv conv)
     797             : {
     798         154 :   PetscFunctionBegin;
     799         154 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     800         462 :   PetscValidLogicalCollectiveEnum(eps,conv,2);
     801         154 :   switch (conv) {
     802          19 :     case EPS_CONV_ABS:  eps->converged = EPSConvergedAbsolute; break;
     803          85 :     case EPS_CONV_REL:  eps->converged = EPSConvergedRelative; break;
     804          50 :     case EPS_CONV_NORM: eps->converged = EPSConvergedNorm; break;
     805           0 :     case EPS_CONV_USER:
     806           0 :       PetscCheck(eps->convergeduser,PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetConvergenceTestFunction() first");
     807           0 :       eps->converged = eps->convergeduser;
     808           0 :       break;
     809           0 :     default:
     810           0 :       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
     811             :   }
     812         154 :   eps->conv = conv;
     813         154 :   PetscFunctionReturn(PETSC_SUCCESS);
     814             : }
     815             : 
     816             : /*@
     817             :    EPSGetConvergenceTest - Gets the method used to compute the error estimate
     818             :    used in the convergence test.
     819             : 
     820             :    Not Collective
     821             : 
     822             :    Input Parameters:
     823             : .  eps   - eigensolver context obtained from EPSCreate()
     824             : 
     825             :    Output Parameters:
     826             : .  conv  - the type of convergence test
     827             : 
     828             :    Level: intermediate
     829             : 
     830             : .seealso: EPSSetConvergenceTest(), EPSConv
     831             : @*/
     832           1 : PetscErrorCode EPSGetConvergenceTest(EPS eps,EPSConv *conv)
     833             : {
     834           1 :   PetscFunctionBegin;
     835           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     836           1 :   PetscAssertPointer(conv,2);
     837           1 :   *conv = eps->conv;
     838           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     839             : }
     840             : 
     841             : /*@C
     842             :    EPSSetStoppingTestFunction - Sets a function to decide when to stop the outer
     843             :    iteration of the eigensolver.
     844             : 
     845             :    Logically Collective
     846             : 
     847             :    Input Parameters:
     848             : +  eps     - eigensolver context obtained from EPSCreate()
     849             : .  func    - stopping test function, see EPSStoppingTestFn for the calling sequence
     850             : .  ctx     - context for private data for the stopping routine (may be NULL)
     851             : -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence
     852             : 
     853             :    Note:
     854             :    Normal usage is to first call the default routine EPSStoppingBasic() and then
     855             :    set reason to EPS_CONVERGED_USER if some user-defined conditions have been
     856             :    met. To let the eigensolver continue iterating, the result must be left as
     857             :    EPS_CONVERGED_ITERATING.
     858             : 
     859             :    Level: advanced
     860             : 
     861             : .seealso: EPSSetStoppingTest(), EPSStoppingBasic()
     862             : @*/
     863          17 : PetscErrorCode EPSSetStoppingTestFunction(EPS eps,EPSStoppingTestFn *func,void* ctx,PetscCtxDestroyFn *destroy)
     864             : {
     865          17 :   PetscFunctionBegin;
     866          17 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     867          17 :   if (eps->stoppingdestroy) PetscCall((*eps->stoppingdestroy)(&eps->stoppingctx));
     868          17 :   eps->stoppinguser    = func;
     869          17 :   eps->stoppingdestroy = destroy;
     870          17 :   eps->stoppingctx     = ctx;
     871          17 :   if (func == EPSStoppingBasic) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_BASIC));
     872          17 :   else if (func == EPSStoppingThreshold) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_THRESHOLD));
     873             :   else {
     874           2 :     eps->stop     = EPS_STOP_USER;
     875           2 :     eps->stopping = eps->stoppinguser;
     876             :   }
     877          17 :   PetscFunctionReturn(PETSC_SUCCESS);
     878             : }
     879             : 
     880             : /*@
     881             :    EPSSetStoppingTest - Specifies how to decide the termination of the outer
     882             :    loop of the eigensolver.
     883             : 
     884             :    Logically Collective
     885             : 
     886             :    Input Parameters:
     887             : +  eps  - eigensolver context obtained from EPSCreate()
     888             : -  stop - the type of stopping test
     889             : 
     890             :    Options Database Keys:
     891             : +  -eps_stop_basic     - Sets the default stopping test
     892             : .  -eps_stop_threshold - Sets the threshold stopping test
     893             : -  -eps_stop_user      - Selects the user-defined stopping test
     894             : 
     895             :    Note:
     896             :    The parameter 'stop' can have one of these values
     897             : +     EPS_STOP_BASIC     - default stopping test
     898             : .     EPS_STOP_THRESHOLD - threshold stopping test)
     899             : -     EPS_STOP_USER      - function set by EPSSetStoppingTestFunction()
     900             : 
     901             :    Level: advanced
     902             : 
     903             : .seealso: EPSGetStoppingTest(), EPSSetStoppingTestFunction(), EPSSetConvergenceTest(), EPSStop
     904             : @*/
     905          31 : PetscErrorCode EPSSetStoppingTest(EPS eps,EPSStop stop)
     906             : {
     907          31 :   PetscFunctionBegin;
     908          31 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     909          93 :   PetscValidLogicalCollectiveEnum(eps,stop,2);
     910          31 :   switch (stop) {
     911           1 :     case EPS_STOP_BASIC: eps->stopping = EPSStoppingBasic; break;
     912          30 :     case EPS_STOP_THRESHOLD: eps->stopping = EPSStoppingThreshold; break;
     913           0 :     case EPS_STOP_USER:
     914           0 :       PetscCheck(eps->stoppinguser,PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetStoppingTestFunction() first");
     915           0 :       eps->stopping = eps->stoppinguser;
     916           0 :       break;
     917           0 :     default:
     918           0 :       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
     919             :   }
     920          31 :   eps->stop = stop;
     921          31 :   PetscFunctionReturn(PETSC_SUCCESS);
     922             : }
     923             : 
     924             : /*@
     925             :    EPSGetStoppingTest - Gets the method used to decide the termination of the outer
     926             :    loop of the eigensolver.
     927             : 
     928             :    Not Collective
     929             : 
     930             :    Input Parameters:
     931             : .  eps   - eigensolver context obtained from EPSCreate()
     932             : 
     933             :    Output Parameters:
     934             : .  stop  - the type of stopping test
     935             : 
     936             :    Level: advanced
     937             : 
     938             : .seealso: EPSSetStoppingTest(), EPSStop
     939             : @*/
     940          27 : PetscErrorCode EPSGetStoppingTest(EPS eps,EPSStop *stop)
     941             : {
     942          27 :   PetscFunctionBegin;
     943          27 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     944          27 :   PetscAssertPointer(stop,2);
     945          27 :   *stop = eps->stop;
     946          27 :   PetscFunctionReturn(PETSC_SUCCESS);
     947             : }
     948             : 
     949             : /*@
     950             :    EPSSetProblemType - Specifies the type of the eigenvalue problem.
     951             : 
     952             :    Logically Collective
     953             : 
     954             :    Input Parameters:
     955             : +  eps      - the eigensolver context
     956             : -  type     - a known type of eigenvalue problem
     957             : 
     958             :    Options Database Keys:
     959             : +  -eps_hermitian - Hermitian eigenvalue problem
     960             : .  -eps_gen_hermitian - generalized Hermitian eigenvalue problem
     961             : .  -eps_non_hermitian - non-Hermitian eigenvalue problem
     962             : .  -eps_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
     963             : .  -eps_pos_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
     964             :    with positive semi-definite B
     965             : .  -eps_gen_indefinite - generalized Hermitian-indefinite eigenvalue problem
     966             : -  -eps_bse - structured Bethe-Salpeter eigenvalue problem
     967             : 
     968             :    Notes:
     969             :    This function must be used to instruct SLEPc to exploit symmetry or other
     970             :    kind of structure. If no
     971             :    problem type is specified, by default a non-Hermitian problem is assumed
     972             :    (either standard or generalized). If the user knows that the problem is
     973             :    Hermitian (i.e. A=A^H) or generalized Hermitian (i.e. A=A^H, B=B^H, and
     974             :    B positive definite) then it is recommended to set the problem type so
     975             :    that eigensolver can exploit these properties.
     976             : 
     977             :    If the user does not call this function, the solver will use a reasonable
     978             :    guess.
     979             : 
     980             :    For structured problem types such as EPS_BSE, the matrices passed in via
     981             :    EPSSetOperators() must have been created with the corresponding helper
     982             :    function, i.e., MatCreateBSE().
     983             : 
     984             :    Level: intermediate
     985             : 
     986             : .seealso: EPSSetOperators(), EPSSetType(), EPSGetProblemType(), EPSProblemType
     987             : @*/
     988         667 : PetscErrorCode EPSSetProblemType(EPS eps,EPSProblemType type)
     989             : {
     990         667 :   PetscFunctionBegin;
     991         667 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     992        2001 :   PetscValidLogicalCollectiveEnum(eps,type,2);
     993         667 :   if (type == eps->problem_type) PetscFunctionReturn(PETSC_SUCCESS);
     994         642 :   switch (type) {
     995         253 :     case EPS_HEP:
     996         253 :       eps->isgeneralized = PETSC_FALSE;
     997         253 :       eps->ishermitian = PETSC_TRUE;
     998         253 :       eps->ispositive = PETSC_FALSE;
     999         253 :       eps->isstructured = PETSC_FALSE;
    1000         253 :       break;
    1001         166 :     case EPS_NHEP:
    1002         166 :       eps->isgeneralized = PETSC_FALSE;
    1003         166 :       eps->ishermitian = PETSC_FALSE;
    1004         166 :       eps->ispositive = PETSC_FALSE;
    1005         166 :       eps->isstructured = PETSC_FALSE;
    1006         166 :       break;
    1007         126 :     case EPS_GHEP:
    1008         126 :       eps->isgeneralized = PETSC_TRUE;
    1009         126 :       eps->ishermitian = PETSC_TRUE;
    1010         126 :       eps->ispositive = PETSC_TRUE;
    1011         126 :       eps->isstructured = PETSC_FALSE;
    1012         126 :       break;
    1013          59 :     case EPS_GNHEP:
    1014          59 :       eps->isgeneralized = PETSC_TRUE;
    1015          59 :       eps->ishermitian = PETSC_FALSE;
    1016          59 :       eps->ispositive = PETSC_FALSE;
    1017          59 :       eps->isstructured = PETSC_FALSE;
    1018          59 :       break;
    1019           2 :     case EPS_PGNHEP:
    1020           2 :       eps->isgeneralized = PETSC_TRUE;
    1021           2 :       eps->ishermitian = PETSC_FALSE;
    1022           2 :       eps->ispositive = PETSC_TRUE;
    1023           2 :       eps->isstructured = PETSC_FALSE;
    1024           2 :       break;
    1025          16 :     case EPS_GHIEP:
    1026          16 :       eps->isgeneralized = PETSC_TRUE;
    1027          16 :       eps->ishermitian = PETSC_TRUE;
    1028          16 :       eps->ispositive = PETSC_FALSE;
    1029          16 :       eps->isstructured = PETSC_FALSE;
    1030          16 :       break;
    1031          20 :     case EPS_BSE:
    1032          20 :       eps->isgeneralized = PETSC_FALSE;
    1033          20 :       eps->ishermitian = PETSC_FALSE;
    1034          20 :       eps->ispositive = PETSC_FALSE;
    1035          20 :       eps->isstructured = PETSC_TRUE;
    1036          20 :       break;
    1037           0 :     default:
    1038           0 :       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
    1039             :   }
    1040         642 :   eps->problem_type = type;
    1041         642 :   eps->state = EPS_STATE_INITIAL;
    1042         642 :   PetscFunctionReturn(PETSC_SUCCESS);
    1043             : }
    1044             : 
    1045             : /*@
    1046             :    EPSGetProblemType - Gets the problem type from the EPS object.
    1047             : 
    1048             :    Not Collective
    1049             : 
    1050             :    Input Parameter:
    1051             : .  eps - the eigensolver context
    1052             : 
    1053             :    Output Parameter:
    1054             : .  type - the problem type
    1055             : 
    1056             :    Level: intermediate
    1057             : 
    1058             : .seealso: EPSSetProblemType(), EPSProblemType
    1059             : @*/
    1060          60 : PetscErrorCode EPSGetProblemType(EPS eps,EPSProblemType *type)
    1061             : {
    1062          60 :   PetscFunctionBegin;
    1063          60 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1064          60 :   PetscAssertPointer(type,2);
    1065          60 :   *type = eps->problem_type;
    1066          60 :   PetscFunctionReturn(PETSC_SUCCESS);
    1067             : }
    1068             : 
    1069             : /*@
    1070             :    EPSSetExtraction - Specifies the type of extraction technique to be employed
    1071             :    by the eigensolver.
    1072             : 
    1073             :    Logically Collective
    1074             : 
    1075             :    Input Parameters:
    1076             : +  eps  - the eigensolver context
    1077             : -  extr - a known type of extraction
    1078             : 
    1079             :    Options Database Keys:
    1080             : +  -eps_ritz - Rayleigh-Ritz extraction
    1081             : .  -eps_harmonic - harmonic Ritz extraction
    1082             : .  -eps_harmonic_relative - harmonic Ritz extraction relative to the eigenvalue
    1083             : .  -eps_harmonic_right - harmonic Ritz extraction for rightmost eigenvalues
    1084             : .  -eps_harmonic_largest - harmonic Ritz extraction for largest magnitude
    1085             :    (without target)
    1086             : .  -eps_refined - refined Ritz extraction
    1087             : -  -eps_refined_harmonic - refined harmonic Ritz extraction
    1088             : 
    1089             :    Notes:
    1090             :    Not all eigensolvers support all types of extraction. See the SLEPc
    1091             :    Users Manual for details.
    1092             : 
    1093             :    By default, a standard Rayleigh-Ritz extraction is used. Other extractions
    1094             :    may be useful when computing interior eigenvalues.
    1095             : 
    1096             :    Harmonic-type extractions are used in combination with a 'target'.
    1097             : 
    1098             :    Level: advanced
    1099             : 
    1100             : .seealso: EPSSetTarget(), EPSGetExtraction(), EPSExtraction
    1101             : @*/
    1102          14 : PetscErrorCode EPSSetExtraction(EPS eps,EPSExtraction extr)
    1103             : {
    1104          14 :   PetscFunctionBegin;
    1105          14 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1106          42 :   PetscValidLogicalCollectiveEnum(eps,extr,2);
    1107          14 :   if (eps->extraction != extr) {
    1108          14 :     eps->state      = EPS_STATE_INITIAL;
    1109          14 :     eps->extraction = extr;
    1110             :   }
    1111          14 :   PetscFunctionReturn(PETSC_SUCCESS);
    1112             : }
    1113             : 
    1114             : /*@
    1115             :    EPSGetExtraction - Gets the extraction type used by the EPS object.
    1116             : 
    1117             :    Not Collective
    1118             : 
    1119             :    Input Parameter:
    1120             : .  eps - the eigensolver context
    1121             : 
    1122             :    Output Parameter:
    1123             : .  extr - name of extraction type
    1124             : 
    1125             :    Level: advanced
    1126             : 
    1127             : .seealso: EPSSetExtraction(), EPSExtraction
    1128             : @*/
    1129           2 : PetscErrorCode EPSGetExtraction(EPS eps,EPSExtraction *extr)
    1130             : {
    1131           2 :   PetscFunctionBegin;
    1132           2 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1133           2 :   PetscAssertPointer(extr,2);
    1134           2 :   *extr = eps->extraction;
    1135           2 :   PetscFunctionReturn(PETSC_SUCCESS);
    1136             : }
    1137             : 
    1138             : /*@
    1139             :    EPSSetBalance - Specifies the balancing technique to be employed by the
    1140             :    eigensolver, and some parameters associated to it.
    1141             : 
    1142             :    Logically Collective
    1143             : 
    1144             :    Input Parameters:
    1145             : +  eps    - the eigensolver context
    1146             : .  bal    - the balancing method, one of EPS_BALANCE_NONE, EPS_BALANCE_ONESIDE,
    1147             :             EPS_BALANCE_TWOSIDE, or EPS_BALANCE_USER
    1148             : .  its    - number of iterations of the balancing algorithm
    1149             : -  cutoff - cutoff value
    1150             : 
    1151             :    Options Database Keys:
    1152             : +  -eps_balance <method> - the balancing method, where <method> is one of
    1153             :                            'none', 'oneside', 'twoside', or 'user'
    1154             : .  -eps_balance_its <its> - number of iterations
    1155             : -  -eps_balance_cutoff <cutoff> - cutoff value
    1156             : 
    1157             :    Notes:
    1158             :    When balancing is enabled, the solver works implicitly with matrix DAD^-1,
    1159             :    where D is an appropriate diagonal matrix. This improves the accuracy of
    1160             :    the computed results in some cases. See the SLEPc Users Manual for details.
    1161             : 
    1162             :    Balancing makes sense only for non-Hermitian problems when the required
    1163             :    precision is high (i.e. a small tolerance such as 1e-15).
    1164             : 
    1165             :    By default, balancing is disabled. The two-sided method is much more
    1166             :    effective than the one-sided counterpart, but it requires the system
    1167             :    matrices to have the MatMultTranspose operation defined.
    1168             : 
    1169             :    The parameter 'its' is the number of iterations performed by the method. The
    1170             :    cutoff value is used only in the two-side variant. Use PETSC_DETERMINE to assign
    1171             :    a reasonably good value, or PETSC_CURRENT to leave the value unchanged.
    1172             : 
    1173             :    User-defined balancing is allowed provided that the corresponding matrix
    1174             :    is set via STSetBalanceMatrix.
    1175             : 
    1176             :    Level: intermediate
    1177             : 
    1178             : .seealso: EPSGetBalance(), EPSBalance, STSetBalanceMatrix()
    1179             : @*/
    1180          16 : PetscErrorCode EPSSetBalance(EPS eps,EPSBalance bal,PetscInt its,PetscReal cutoff)
    1181             : {
    1182          16 :   PetscFunctionBegin;
    1183          16 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1184          48 :   PetscValidLogicalCollectiveEnum(eps,bal,2);
    1185          48 :   PetscValidLogicalCollectiveInt(eps,its,3);
    1186          48 :   PetscValidLogicalCollectiveReal(eps,cutoff,4);
    1187          16 :   switch (bal) {
    1188          16 :     case EPS_BALANCE_NONE:
    1189             :     case EPS_BALANCE_ONESIDE:
    1190             :     case EPS_BALANCE_TWOSIDE:
    1191             :     case EPS_BALANCE_USER:
    1192          16 :       if (eps->balance != bal) {
    1193          14 :         eps->state = EPS_STATE_INITIAL;
    1194          14 :         eps->balance = bal;
    1195             :       }
    1196          16 :       break;
    1197           0 :     default:
    1198           0 :       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'bal'");
    1199             :   }
    1200          16 :   if (its==PETSC_DETERMINE) eps->balance_its = 5;
    1201          16 :   else if (its!=PETSC_CURRENT) {
    1202          16 :     PetscCheck(its>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be > 0");
    1203          16 :     eps->balance_its = its;
    1204             :   }
    1205          16 :   if (cutoff==(PetscReal)PETSC_DETERMINE) eps->balance_cutoff = 1e-8;
    1206          16 :   else if (cutoff!=(PetscReal)PETSC_CURRENT) {
    1207          16 :     PetscCheck(cutoff>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of cutoff. Must be > 0");
    1208          16 :     eps->balance_cutoff = cutoff;
    1209             :   }
    1210          16 :   PetscFunctionReturn(PETSC_SUCCESS);
    1211             : }
    1212             : 
    1213             : /*@
    1214             :    EPSGetBalance - Gets the balancing type used by the EPS object, and the
    1215             :    associated parameters.
    1216             : 
    1217             :    Not Collective
    1218             : 
    1219             :    Input Parameter:
    1220             : .  eps - the eigensolver context
    1221             : 
    1222             :    Output Parameters:
    1223             : +  bal    - the balancing method
    1224             : .  its    - number of iterations of the balancing algorithm
    1225             : -  cutoff - cutoff value
    1226             : 
    1227             :    Level: intermediate
    1228             : 
    1229             :    Note:
    1230             :    The user can specify NULL for any parameter that is not needed.
    1231             : 
    1232             : .seealso: EPSSetBalance(), EPSBalance
    1233             : @*/
    1234           1 : PetscErrorCode EPSGetBalance(EPS eps,EPSBalance *bal,PetscInt *its,PetscReal *cutoff)
    1235             : {
    1236           1 :   PetscFunctionBegin;
    1237           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1238           1 :   if (bal)    *bal = eps->balance;
    1239           1 :   if (its)    *its = eps->balance_its;
    1240           1 :   if (cutoff) *cutoff = eps->balance_cutoff;
    1241           1 :   PetscFunctionReturn(PETSC_SUCCESS);
    1242             : }
    1243             : 
    1244             : /*@
    1245             :    EPSSetTwoSided - Sets the solver to use a two-sided variant so that left
    1246             :    eigenvectors are also computed.
    1247             : 
    1248             :    Logically Collective
    1249             : 
    1250             :    Input Parameters:
    1251             : +  eps      - the eigensolver context
    1252             : -  twosided - whether the two-sided variant is to be used or not
    1253             : 
    1254             :    Options Database Keys:
    1255             : .  -eps_two_sided <boolean> - Sets/resets the twosided flag
    1256             : 
    1257             :    Notes:
    1258             :    If the user sets twosided=PETSC_TRUE then the solver uses a variant of
    1259             :    the algorithm that computes both right and left eigenvectors. This is
    1260             :    usually much more costly. This option is not available in all solvers.
    1261             : 
    1262             :    When using two-sided solvers, the problem matrices must have both the
    1263             :    MatMult and MatMultTranspose operations defined.
    1264             : 
    1265             :    Level: advanced
    1266             : 
    1267             : .seealso: EPSGetTwoSided(), EPSGetLeftEigenvector()
    1268             : @*/
    1269          26 : PetscErrorCode EPSSetTwoSided(EPS eps,PetscBool twosided)
    1270             : {
    1271          26 :   PetscFunctionBegin;
    1272          26 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1273          78 :   PetscValidLogicalCollectiveBool(eps,twosided,2);
    1274          26 :   if (twosided!=eps->twosided) {
    1275          17 :     eps->twosided = twosided;
    1276          17 :     eps->state    = EPS_STATE_INITIAL;
    1277             :   }
    1278          26 :   PetscFunctionReturn(PETSC_SUCCESS);
    1279             : }
    1280             : 
    1281             : /*@
    1282             :    EPSGetTwoSided - Returns the flag indicating whether a two-sided variant
    1283             :    of the algorithm is being used or not.
    1284             : 
    1285             :    Not Collective
    1286             : 
    1287             :    Input Parameter:
    1288             : .  eps - the eigensolver context
    1289             : 
    1290             :    Output Parameter:
    1291             : .  twosided - the returned flag
    1292             : 
    1293             :    Level: advanced
    1294             : 
    1295             : .seealso: EPSSetTwoSided()
    1296             : @*/
    1297           5 : PetscErrorCode EPSGetTwoSided(EPS eps,PetscBool *twosided)
    1298             : {
    1299           5 :   PetscFunctionBegin;
    1300           5 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1301           5 :   PetscAssertPointer(twosided,2);
    1302           5 :   *twosided = eps->twosided;
    1303           5 :   PetscFunctionReturn(PETSC_SUCCESS);
    1304             : }
    1305             : 
    1306             : /*@
    1307             :    EPSSetTrueResidual - Specifies if the solver must compute the true residual
    1308             :    explicitly or not.
    1309             : 
    1310             :    Logically Collective
    1311             : 
    1312             :    Input Parameters:
    1313             : +  eps     - the eigensolver context
    1314             : -  trueres - whether true residuals are required or not
    1315             : 
    1316             :    Options Database Keys:
    1317             : .  -eps_true_residual <boolean> - Sets/resets the boolean flag 'trueres'
    1318             : 
    1319             :    Notes:
    1320             :    If the user sets trueres=PETSC_TRUE then the solver explicitly computes
    1321             :    the true residual for each eigenpair approximation, and uses it for
    1322             :    convergence testing. Computing the residual is usually an expensive
    1323             :    operation. Some solvers (e.g., Krylov solvers) can avoid this computation
    1324             :    by using a cheap estimate of the residual norm, but this may sometimes
    1325             :    give inaccurate results (especially if a spectral transform is being
    1326             :    used). On the contrary, preconditioned eigensolvers (e.g., Davidson solvers)
    1327             :    do rely on computing the true residual, so this option is irrelevant for them.
    1328             : 
    1329             :    Level: advanced
    1330             : 
    1331             : .seealso: EPSGetTrueResidual()
    1332             : @*/
    1333           9 : PetscErrorCode EPSSetTrueResidual(EPS eps,PetscBool trueres)
    1334             : {
    1335           9 :   PetscFunctionBegin;
    1336           9 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1337          27 :   PetscValidLogicalCollectiveBool(eps,trueres,2);
    1338           9 :   eps->trueres = trueres;
    1339           9 :   PetscFunctionReturn(PETSC_SUCCESS);
    1340             : }
    1341             : 
    1342             : /*@
    1343             :    EPSGetTrueResidual - Returns the flag indicating whether true
    1344             :    residuals must be computed explicitly or not.
    1345             : 
    1346             :    Not Collective
    1347             : 
    1348             :    Input Parameter:
    1349             : .  eps - the eigensolver context
    1350             : 
    1351             :    Output Parameter:
    1352             : .  trueres - the returned flag
    1353             : 
    1354             :    Level: advanced
    1355             : 
    1356             : .seealso: EPSSetTrueResidual()
    1357             : @*/
    1358           8 : PetscErrorCode EPSGetTrueResidual(EPS eps,PetscBool *trueres)
    1359             : {
    1360           8 :   PetscFunctionBegin;
    1361           8 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1362           8 :   PetscAssertPointer(trueres,2);
    1363           8 :   *trueres = eps->trueres;
    1364           8 :   PetscFunctionReturn(PETSC_SUCCESS);
    1365             : }
    1366             : 
    1367             : /*@
    1368             :    EPSSetTrackAll - Specifies if the solver must compute the residual norm of all
    1369             :    approximate eigenpairs or not.
    1370             : 
    1371             :    Logically Collective
    1372             : 
    1373             :    Input Parameters:
    1374             : +  eps      - the eigensolver context
    1375             : -  trackall - whether to compute all residuals or not
    1376             : 
    1377             :    Notes:
    1378             :    If the user sets trackall=PETSC_TRUE then the solver computes (or estimates)
    1379             :    the residual norm for each eigenpair approximation. Computing the residual is
    1380             :    usually an expensive operation and solvers commonly compute only the residual
    1381             :    associated to the first unconverged eigenpair.
    1382             : 
    1383             :    The option '-eps_monitor_all' automatically activates this option.
    1384             : 
    1385             :    Level: developer
    1386             : 
    1387             : .seealso: EPSGetTrackAll()
    1388             : @*/
    1389         135 : PetscErrorCode EPSSetTrackAll(EPS eps,PetscBool trackall)
    1390             : {
    1391         135 :   PetscFunctionBegin;
    1392         135 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1393         405 :   PetscValidLogicalCollectiveBool(eps,trackall,2);
    1394         135 :   eps->trackall = trackall;
    1395         135 :   PetscFunctionReturn(PETSC_SUCCESS);
    1396             : }
    1397             : 
    1398             : /*@
    1399             :    EPSGetTrackAll - Returns the flag indicating whether all residual norms must
    1400             :    be computed or not.
    1401             : 
    1402             :    Not Collective
    1403             : 
    1404             :    Input Parameter:
    1405             : .  eps - the eigensolver context
    1406             : 
    1407             :    Output Parameter:
    1408             : .  trackall - the returned flag
    1409             : 
    1410             :    Level: developer
    1411             : 
    1412             : .seealso: EPSSetTrackAll()
    1413             : @*/
    1414           1 : PetscErrorCode EPSGetTrackAll(EPS eps,PetscBool *trackall)
    1415             : {
    1416           1 :   PetscFunctionBegin;
    1417           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1418           1 :   PetscAssertPointer(trackall,2);
    1419           1 :   *trackall = eps->trackall;
    1420           1 :   PetscFunctionReturn(PETSC_SUCCESS);
    1421             : }
    1422             : 
    1423             : /*@
    1424             :    EPSSetPurify - Deactivate eigenvector purification (which is activated by default).
    1425             : 
    1426             :    Logically Collective
    1427             : 
    1428             :    Input Parameters:
    1429             : +  eps    - the eigensolver context
    1430             : -  purify - whether purification is required or not
    1431             : 
    1432             :    Options Database Keys:
    1433             : .  -eps_purify <boolean> - Sets/resets the boolean flag 'purify'
    1434             : 
    1435             :    Notes:
    1436             :    By default, eigenvectors of generalized symmetric eigenproblems are purified
    1437             :    in order to purge directions in the nullspace of matrix B. If the user knows
    1438             :    that B is non-singular, then purification can be safely deactivated and some
    1439             :    computational cost is avoided (this is particularly important in interval computations).
    1440             : 
    1441             :    Level: intermediate
    1442             : 
    1443             : .seealso: EPSGetPurify(), EPSSetInterval()
    1444             : @*/
    1445           4 : PetscErrorCode EPSSetPurify(EPS eps,PetscBool purify)
    1446             : {
    1447           4 :   PetscFunctionBegin;
    1448           4 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1449          12 :   PetscValidLogicalCollectiveBool(eps,purify,2);
    1450           4 :   if (purify!=eps->purify) {
    1451           2 :     eps->purify = purify;
    1452           2 :     eps->state  = EPS_STATE_INITIAL;
    1453             :   }
    1454           4 :   PetscFunctionReturn(PETSC_SUCCESS);
    1455             : }
    1456             : 
    1457             : /*@
    1458             :    EPSGetPurify - Returns the flag indicating whether purification is activated
    1459             :    or not.
    1460             : 
    1461             :    Not Collective
    1462             : 
    1463             :    Input Parameter:
    1464             : .  eps - the eigensolver context
    1465             : 
    1466             :    Output Parameter:
    1467             : .  purify - the returned flag
    1468             : 
    1469             :    Level: intermediate
    1470             : 
    1471             : .seealso: EPSSetPurify()
    1472             : @*/
    1473           1 : PetscErrorCode EPSGetPurify(EPS eps,PetscBool *purify)
    1474             : {
    1475           1 :   PetscFunctionBegin;
    1476           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1477           1 :   PetscAssertPointer(purify,2);
    1478           1 :   *purify = eps->purify;
    1479           1 :   PetscFunctionReturn(PETSC_SUCCESS);
    1480             : }
    1481             : 
    1482             : /*@
    1483             :    EPSSetOptionsPrefix - Sets the prefix used for searching for all
    1484             :    EPS options in the database.
    1485             : 
    1486             :    Logically Collective
    1487             : 
    1488             :    Input Parameters:
    1489             : +  eps - the eigensolver context
    1490             : -  prefix - the prefix string to prepend to all EPS option requests
    1491             : 
    1492             :    Notes:
    1493             :    A hyphen (-) must NOT be given at the beginning of the prefix name.
    1494             :    The first character of all runtime options is AUTOMATICALLY the
    1495             :    hyphen.
    1496             : 
    1497             :    For example, to distinguish between the runtime options for two
    1498             :    different EPS contexts, one could call
    1499             : .vb
    1500             :       EPSSetOptionsPrefix(eps1,"eig1_")
    1501             :       EPSSetOptionsPrefix(eps2,"eig2_")
    1502             : .ve
    1503             : 
    1504             :    Level: advanced
    1505             : 
    1506             : .seealso: EPSAppendOptionsPrefix(), EPSGetOptionsPrefix()
    1507             : @*/
    1508         186 : PetscErrorCode EPSSetOptionsPrefix(EPS eps,const char *prefix)
    1509             : {
    1510         186 :   PetscFunctionBegin;
    1511         186 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1512         186 :   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
    1513         186 :   PetscCall(STSetOptionsPrefix(eps->st,prefix));
    1514         186 :   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
    1515         186 :   PetscCall(BVSetOptionsPrefix(eps->V,prefix));
    1516         186 :   if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
    1517         186 :   PetscCall(DSSetOptionsPrefix(eps->ds,prefix));
    1518         186 :   if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
    1519         186 :   PetscCall(RGSetOptionsPrefix(eps->rg,prefix));
    1520         186 :   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)eps,prefix));
    1521         186 :   PetscFunctionReturn(PETSC_SUCCESS);
    1522             : }
    1523             : 
    1524             : /*@
    1525             :    EPSAppendOptionsPrefix - Appends to the prefix used for searching for all
    1526             :    EPS options in the database.
    1527             : 
    1528             :    Logically Collective
    1529             : 
    1530             :    Input Parameters:
    1531             : +  eps - the eigensolver context
    1532             : -  prefix - the prefix string to prepend to all EPS option requests
    1533             : 
    1534             :    Notes:
    1535             :    A hyphen (-) must NOT be given at the beginning of the prefix name.
    1536             :    The first character of all runtime options is AUTOMATICALLY the hyphen.
    1537             : 
    1538             :    Level: advanced
    1539             : 
    1540             : .seealso: EPSSetOptionsPrefix(), EPSGetOptionsPrefix()
    1541             : @*/
    1542         150 : PetscErrorCode EPSAppendOptionsPrefix(EPS eps,const char *prefix)
    1543             : {
    1544         150 :   PetscFunctionBegin;
    1545         150 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1546         150 :   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
    1547         150 :   PetscCall(STAppendOptionsPrefix(eps->st,prefix));
    1548         150 :   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
    1549         150 :   PetscCall(BVAppendOptionsPrefix(eps->V,prefix));
    1550         150 :   if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
    1551         150 :   PetscCall(DSAppendOptionsPrefix(eps->ds,prefix));
    1552         150 :   if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
    1553         150 :   PetscCall(RGAppendOptionsPrefix(eps->rg,prefix));
    1554         150 :   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)eps,prefix));
    1555         150 :   PetscFunctionReturn(PETSC_SUCCESS);
    1556             : }
    1557             : 
    1558             : /*@
    1559             :    EPSGetOptionsPrefix - Gets the prefix used for searching for all
    1560             :    EPS options in the database.
    1561             : 
    1562             :    Not Collective
    1563             : 
    1564             :    Input Parameters:
    1565             : .  eps - the eigensolver context
    1566             : 
    1567             :    Output Parameters:
    1568             : .  prefix - pointer to the prefix string used is returned
    1569             : 
    1570             :    Note:
    1571             :    On the Fortran side, the user should pass in a string 'prefix' of
    1572             :    sufficient length to hold the prefix.
    1573             : 
    1574             :    Level: advanced
    1575             : 
    1576             : .seealso: EPSSetOptionsPrefix(), EPSAppendOptionsPrefix()
    1577             : @*/
    1578          30 : PetscErrorCode EPSGetOptionsPrefix(EPS eps,const char *prefix[])
    1579             : {
    1580          30 :   PetscFunctionBegin;
    1581          30 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1582          30 :   PetscAssertPointer(prefix,2);
    1583          30 :   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)eps,prefix));
    1584          30 :   PetscFunctionReturn(PETSC_SUCCESS);
    1585             : }

Generated by: LCOV version 1.14