LCOV - code coverage report
Current view: top level - eps/impls/davidson/jd - jd.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 214 237 90.3 %
Date: 2024-04-24 00:57:47 Functions: 23 24 95.8 %
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             :    SLEPc eigensolver: "jd"
      12             : 
      13             :    Method: Jacobi-Davidson
      14             : 
      15             :    Algorithm:
      16             : 
      17             :        Jacobi-Davidson with various subspace extraction and
      18             :        restart techniques.
      19             : 
      20             :    References:
      21             : 
      22             :        [1] G.L.G. Sleijpen and H.A. van der Vorst, "A Jacobi-Davidson
      23             :            iteration method for linear eigenvalue problems", SIAM J.
      24             :            Matrix Anal. Appl. 17(2):401-425, 1996.
      25             : 
      26             :        [2] E. Romero and J.E. Roman, "A parallel implementation of
      27             :            Davidson methods for large-scale eigenvalue problems in
      28             :            SLEPc", ACM Trans. Math. Software 40(2), Article 13, 2014.
      29             : */
      30             : 
      31             : #include <slepc/private/epsimpl.h>                /*I "slepceps.h" I*/
      32             : #include <../src/eps/impls/davidson/davidson.h>
      33             : 
      34          20 : static PetscErrorCode EPSSetFromOptions_JD(EPS eps,PetscOptionItems *PetscOptionsObject)
      35             : {
      36          20 :   PetscBool      flg,flg2,op,orth;
      37          20 :   PetscInt       opi,opi0;
      38          20 :   PetscReal      opf;
      39             : 
      40          20 :   PetscFunctionBegin;
      41          20 :   PetscOptionsHeadBegin(PetscOptionsObject,"EPS Jacobi-Davidson (JD) Options");
      42             : 
      43          20 :     PetscCall(EPSJDGetKrylovStart(eps,&op));
      44          20 :     PetscCall(PetscOptionsBool("-eps_jd_krylov_start","Start the search subspace with a Krylov basis","EPSJDSetKrylovStart",op,&op,&flg));
      45          20 :     if (flg) PetscCall(EPSJDSetKrylovStart(eps,op));
      46             : 
      47          20 :     PetscCall(EPSJDGetBOrth(eps,&orth));
      48          20 :     PetscCall(PetscOptionsBool("-eps_jd_borth","Use B-orthogonalization in the search subspace","EPSJDSetBOrth",op,&op,&flg));
      49          20 :     if (flg) PetscCall(EPSJDSetBOrth(eps,op));
      50             : 
      51          20 :     PetscCall(EPSJDGetBlockSize(eps,&opi));
      52          20 :     PetscCall(PetscOptionsInt("-eps_jd_blocksize","Number of vectors to add to the search subspace","EPSJDSetBlockSize",opi,&opi,&flg));
      53          20 :     if (flg) PetscCall(EPSJDSetBlockSize(eps,opi));
      54             : 
      55          20 :     PetscCall(EPSJDGetRestart(eps,&opi,&opi0));
      56          20 :     PetscCall(PetscOptionsInt("-eps_jd_minv","Size of the search subspace after restarting","EPSJDSetRestart",opi,&opi,&flg));
      57          20 :     PetscCall(PetscOptionsInt("-eps_jd_plusk","Number of eigenvectors saved from the previous iteration when restarting","EPSJDSetRestart",opi0,&opi0,&flg2));
      58          20 :     if (flg || flg2) PetscCall(EPSJDSetRestart(eps,opi,opi0));
      59             : 
      60          20 :     PetscCall(EPSJDGetInitialSize(eps,&opi));
      61          20 :     PetscCall(PetscOptionsInt("-eps_jd_initial_size","Initial size of the search subspace","EPSJDSetInitialSize",opi,&opi,&flg));
      62          20 :     if (flg) PetscCall(EPSJDSetInitialSize(eps,opi));
      63             : 
      64          20 :     PetscCall(EPSJDGetFix(eps,&opf));
      65          20 :     PetscCall(PetscOptionsReal("-eps_jd_fix","Tolerance for changing the target in the correction equation","EPSJDSetFix",opf,&opf,&flg));
      66          20 :     if (flg) PetscCall(EPSJDSetFix(eps,opf));
      67             : 
      68          20 :     PetscCall(EPSJDGetConstCorrectionTol(eps,&op));
      69          20 :     PetscCall(PetscOptionsBool("-eps_jd_const_correction_tol","Disable the dynamic stopping criterion when solving the correction equation","EPSJDSetConstCorrectionTol",op,&op,&flg));
      70          20 :     if (flg) PetscCall(EPSJDSetConstCorrectionTol(eps,op));
      71             : 
      72          20 :   PetscOptionsHeadEnd();
      73          20 :   PetscFunctionReturn(PETSC_SUCCESS);
      74             : }
      75             : 
      76          51 : static PetscErrorCode EPSSetDefaultST_JD(EPS eps)
      77             : {
      78          51 :   KSP            ksp;
      79             : 
      80          51 :   PetscFunctionBegin;
      81          51 :   if (!((PetscObject)eps->st)->type_name) {
      82          20 :     PetscCall(STSetType(eps->st,STPRECOND));
      83          20 :     PetscCall(STPrecondSetKSPHasMat(eps->st,PETSC_TRUE));
      84             :   }
      85          51 :   PetscCall(STGetKSP(eps->st,&ksp));
      86          51 :   if (!((PetscObject)ksp)->type_name) {
      87          21 :     PetscCall(KSPSetType(ksp,KSPBCGSL));
      88          21 :     PetscCall(KSPSetTolerances(ksp,1e-4,PETSC_DEFAULT,PETSC_DEFAULT,90));
      89             :   }
      90          51 :   PetscFunctionReturn(PETSC_SUCCESS);
      91             : }
      92             : 
      93          31 : static PetscErrorCode EPSSetUp_JD(EPS eps)
      94             : {
      95          31 :   PetscBool      t;
      96          31 :   KSP            ksp;
      97             : 
      98          31 :   PetscFunctionBegin;
      99             :   /* Setup common for all davidson solvers */
     100          31 :   PetscCall(EPSSetUp_XD(eps));
     101             : 
     102             :   /* Check some constraints */
     103          31 :   PetscCall(STGetKSP(eps->st,&ksp));
     104          31 :   PetscCall(PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&t));
     105          31 :   PetscCheck(!t,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"EPSJD does not work with KSPPREONLY");
     106          31 :   PetscFunctionReturn(PETSC_SUCCESS);
     107             : }
     108             : 
     109           0 : static PetscErrorCode EPSView_JD(EPS eps,PetscViewer viewer)
     110             : {
     111           0 :   PetscBool      isascii,opb;
     112           0 :   PetscReal      opf;
     113           0 :   PetscInt       opi,opi0;
     114             : 
     115           0 :   PetscFunctionBegin;
     116           0 :   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
     117           0 :   if (isascii) {
     118           0 :     PetscCall(EPSXDGetBOrth_XD(eps,&opb));
     119           0 :     if (opb) PetscCall(PetscViewerASCIIPrintf(viewer,"  search subspace is B-orthogonalized\n"));
     120           0 :     else PetscCall(PetscViewerASCIIPrintf(viewer,"  search subspace is orthogonalized\n"));
     121           0 :     PetscCall(EPSXDGetBlockSize_XD(eps,&opi));
     122           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  block size=%" PetscInt_FMT "\n",opi));
     123           0 :     PetscCall(EPSXDGetKrylovStart_XD(eps,&opb));
     124           0 :     if (!opb) PetscCall(PetscViewerASCIIPrintf(viewer,"  type of the initial subspace: non-Krylov\n"));
     125           0 :     else PetscCall(PetscViewerASCIIPrintf(viewer,"  type of the initial subspace: Krylov\n"));
     126           0 :     PetscCall(EPSXDGetRestart_XD(eps,&opi,&opi0));
     127           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  size of the subspace after restarting: %" PetscInt_FMT "\n",opi));
     128           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  number of vectors after restarting from the previous iteration: %" PetscInt_FMT "\n",opi0));
     129             : 
     130           0 :     PetscCall(EPSJDGetFix_JD(eps,&opf));
     131           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  threshold for changing the target in the correction equation (fix): %g\n",(double)opf));
     132             : 
     133           0 :     PetscCall(EPSJDGetConstCorrectionTol_JD(eps,&opb));
     134           0 :     if (!opb) PetscCall(PetscViewerASCIIPrintf(viewer,"  using dynamic tolerance for the correction equation\n"));
     135             :   }
     136           0 :   PetscFunctionReturn(PETSC_SUCCESS);
     137             : }
     138             : 
     139          21 : static PetscErrorCode EPSDestroy_JD(EPS eps)
     140             : {
     141          21 :   PetscFunctionBegin;
     142          21 :   PetscCall(PetscFree(eps->data));
     143          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetKrylovStart_C",NULL));
     144          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetKrylovStart_C",NULL));
     145          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetBlockSize_C",NULL));
     146          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetBlockSize_C",NULL));
     147          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetRestart_C",NULL));
     148          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetRestart_C",NULL));
     149          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetInitialSize_C",NULL));
     150          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetInitialSize_C",NULL));
     151          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetFix_C",NULL));
     152          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetFix_C",NULL));
     153          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetConstCorrectionTol_C",NULL));
     154          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetConstCorrectionTol_C",NULL));
     155          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetBOrth_C",NULL));
     156          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetBOrth_C",NULL));
     157          21 :   PetscFunctionReturn(PETSC_SUCCESS);
     158             : }
     159             : 
     160             : /*@
     161             :    EPSJDSetKrylovStart - Activates or deactivates starting the searching
     162             :    subspace with a Krylov basis.
     163             : 
     164             :    Logically Collective
     165             : 
     166             :    Input Parameters:
     167             : +  eps - the eigenproblem solver context
     168             : -  krylovstart - boolean flag
     169             : 
     170             :    Options Database Key:
     171             : .  -eps_jd_krylov_start - Activates starting the searching subspace with a
     172             :     Krylov basis
     173             : 
     174             :    Level: advanced
     175             : 
     176             : .seealso: EPSJDGetKrylovStart()
     177             : @*/
     178           1 : PetscErrorCode EPSJDSetKrylovStart(EPS eps,PetscBool krylovstart)
     179             : {
     180           1 :   PetscFunctionBegin;
     181           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     182           4 :   PetscValidLogicalCollectiveBool(eps,krylovstart,2);
     183           1 :   PetscTryMethod(eps,"EPSJDSetKrylovStart_C",(EPS,PetscBool),(eps,krylovstart));
     184           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     185             : }
     186             : 
     187             : /*@
     188             :    EPSJDGetKrylovStart - Returns a flag indicating if the searching subspace is started with a
     189             :    Krylov basis.
     190             : 
     191             :    Not Collective
     192             : 
     193             :    Input Parameter:
     194             : .  eps - the eigenproblem solver context
     195             : 
     196             :    Output Parameters:
     197             : .  krylovstart - boolean flag indicating if the searching subspace is started
     198             :    with a Krylov basis
     199             : 
     200             :    Level: advanced
     201             : 
     202             : .seealso: EPSJDSetKrylovStart()
     203             : @*/
     204          20 : PetscErrorCode EPSJDGetKrylovStart(EPS eps,PetscBool *krylovstart)
     205             : {
     206          20 :   PetscFunctionBegin;
     207          20 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     208          20 :   PetscAssertPointer(krylovstart,2);
     209          20 :   PetscUseMethod(eps,"EPSJDGetKrylovStart_C",(EPS,PetscBool*),(eps,krylovstart));
     210          20 :   PetscFunctionReturn(PETSC_SUCCESS);
     211             : }
     212             : 
     213             : /*@
     214             :    EPSJDSetBlockSize - Sets the number of vectors to be added to the searching space
     215             :    in every iteration.
     216             : 
     217             :    Logically Collective
     218             : 
     219             :    Input Parameters:
     220             : +  eps - the eigenproblem solver context
     221             : -  blocksize - number of vectors added to the search space in every iteration
     222             : 
     223             :    Options Database Key:
     224             : .  -eps_jd_blocksize - number of vectors added to the searching space every iteration
     225             : 
     226             :    Level: advanced
     227             : 
     228             : .seealso: EPSJDSetKrylovStart()
     229             : @*/
     230           1 : PetscErrorCode EPSJDSetBlockSize(EPS eps,PetscInt blocksize)
     231             : {
     232           1 :   PetscFunctionBegin;
     233           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     234           4 :   PetscValidLogicalCollectiveInt(eps,blocksize,2);
     235           1 :   PetscTryMethod(eps,"EPSJDSetBlockSize_C",(EPS,PetscInt),(eps,blocksize));
     236           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     237             : }
     238             : 
     239             : /*@
     240             :    EPSJDGetBlockSize - Returns the number of vectors to be added to the searching space
     241             :    in every iteration.
     242             : 
     243             :    Not Collective
     244             : 
     245             :    Input Parameter:
     246             : .  eps - the eigenproblem solver context
     247             : 
     248             :    Output Parameter:
     249             : .  blocksize - number of vectors added to the search space in every iteration
     250             : 
     251             :    Level: advanced
     252             : 
     253             : .seealso: EPSJDSetBlockSize()
     254             : @*/
     255          20 : PetscErrorCode EPSJDGetBlockSize(EPS eps,PetscInt *blocksize)
     256             : {
     257          20 :   PetscFunctionBegin;
     258          20 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     259          20 :   PetscAssertPointer(blocksize,2);
     260          20 :   PetscUseMethod(eps,"EPSJDGetBlockSize_C",(EPS,PetscInt*),(eps,blocksize));
     261          20 :   PetscFunctionReturn(PETSC_SUCCESS);
     262             : }
     263             : 
     264             : /*@
     265             :    EPSJDSetRestart - Sets the number of vectors of the searching space after
     266             :    restarting and the number of vectors saved from the previous iteration.
     267             : 
     268             :    Logically Collective
     269             : 
     270             :    Input Parameters:
     271             : +  eps - the eigenproblem solver context
     272             : .  minv - number of vectors of the searching subspace after restarting
     273             : -  plusk - number of vectors saved from the previous iteration
     274             : 
     275             :    Options Database Keys:
     276             : +  -eps_jd_minv - number of vectors of the searching subspace after restarting
     277             : -  -eps_jd_plusk - number of vectors saved from the previous iteration
     278             : 
     279             :    Level: advanced
     280             : 
     281             : .seealso: EPSJDGetRestart()
     282             : @*/
     283           2 : PetscErrorCode EPSJDSetRestart(EPS eps,PetscInt minv,PetscInt plusk)
     284             : {
     285           2 :   PetscFunctionBegin;
     286           2 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     287           8 :   PetscValidLogicalCollectiveInt(eps,minv,2);
     288           8 :   PetscValidLogicalCollectiveInt(eps,plusk,3);
     289           2 :   PetscTryMethod(eps,"EPSJDSetRestart_C",(EPS,PetscInt,PetscInt),(eps,minv,plusk));
     290           2 :   PetscFunctionReturn(PETSC_SUCCESS);
     291             : }
     292             : 
     293             : /*@
     294             :    EPSJDGetRestart - Gets the number of vectors of the searching space after
     295             :    restarting and the number of vectors saved from the previous iteration.
     296             : 
     297             :    Not Collective
     298             : 
     299             :    Input Parameter:
     300             : .  eps - the eigenproblem solver context
     301             : 
     302             :    Output Parameters:
     303             : +  minv - number of vectors of the searching subspace after restarting
     304             : -  plusk - number of vectors saved from the previous iteration
     305             : 
     306             :    Level: advanced
     307             : 
     308             : .seealso: EPSJDSetRestart()
     309             : @*/
     310          20 : PetscErrorCode EPSJDGetRestart(EPS eps,PetscInt *minv,PetscInt *plusk)
     311             : {
     312          20 :   PetscFunctionBegin;
     313          20 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     314          20 :   PetscUseMethod(eps,"EPSJDGetRestart_C",(EPS,PetscInt*,PetscInt*),(eps,minv,plusk));
     315          20 :   PetscFunctionReturn(PETSC_SUCCESS);
     316             : }
     317             : 
     318             : /*@
     319             :    EPSJDSetInitialSize - Sets the initial size of the searching space.
     320             : 
     321             :    Logically Collective
     322             : 
     323             :    Input Parameters:
     324             : +  eps - the eigenproblem solver context
     325             : -  initialsize - number of vectors of the initial searching subspace
     326             : 
     327             :    Options Database Key:
     328             : .  -eps_jd_initial_size - number of vectors of the initial searching subspace
     329             : 
     330             :    Notes:
     331             :    If EPSJDGetKrylovStart() is PETSC_FALSE and the user provides vectors with
     332             :    EPSSetInitialSpace(), up to initialsize vectors will be used; and if the
     333             :    provided vectors are not enough, the solver completes the subspace with
     334             :    random vectors. In the case of EPSJDGetKrylovStart() being PETSC_TRUE, the solver
     335             :    gets the first vector provided by the user or, if not available, a random vector,
     336             :    and expands the Krylov basis up to initialsize vectors.
     337             : 
     338             :    Level: advanced
     339             : 
     340             : .seealso: EPSJDGetInitialSize(), EPSJDGetKrylovStart()
     341             : @*/
     342           1 : PetscErrorCode EPSJDSetInitialSize(EPS eps,PetscInt initialsize)
     343             : {
     344           1 :   PetscFunctionBegin;
     345           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     346           4 :   PetscValidLogicalCollectiveInt(eps,initialsize,2);
     347           1 :   PetscTryMethod(eps,"EPSJDSetInitialSize_C",(EPS,PetscInt),(eps,initialsize));
     348           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     349             : }
     350             : 
     351             : /*@
     352             :    EPSJDGetInitialSize - Returns the initial size of the searching space.
     353             : 
     354             :    Not Collective
     355             : 
     356             :    Input Parameter:
     357             : .  eps - the eigenproblem solver context
     358             : 
     359             :    Output Parameter:
     360             : .  initialsize - number of vectors of the initial searching subspace
     361             : 
     362             :    Notes:
     363             :    If EPSJDGetKrylovStart() is PETSC_FALSE and the user provides vectors with
     364             :    EPSSetInitialSpace(), up to initialsize vectors will be used; and if the
     365             :    provided vectors are not enough, the solver completes the subspace with
     366             :    random vectors. In the case of EPSJDGetKrylovStart() being PETSC_TRUE, the solver
     367             :    gets the first vector provided by the user or, if not available, a random vector,
     368             :    and expands the Krylov basis up to initialsize vectors.
     369             : 
     370             :    Level: advanced
     371             : 
     372             : .seealso: EPSJDSetInitialSize(), EPSJDGetKrylovStart()
     373             : @*/
     374          20 : PetscErrorCode EPSJDGetInitialSize(EPS eps,PetscInt *initialsize)
     375             : {
     376          20 :   PetscFunctionBegin;
     377          20 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     378          20 :   PetscAssertPointer(initialsize,2);
     379          20 :   PetscUseMethod(eps,"EPSJDGetInitialSize_C",(EPS,PetscInt*),(eps,initialsize));
     380          20 :   PetscFunctionReturn(PETSC_SUCCESS);
     381             : }
     382             : 
     383           1 : static PetscErrorCode EPSJDSetFix_JD(EPS eps,PetscReal fix)
     384             : {
     385           1 :   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
     386             : 
     387           1 :   PetscFunctionBegin;
     388           1 :   if (fix == (PetscReal)PETSC_DEFAULT || fix == (PetscReal)PETSC_DECIDE) fix = 0.01;
     389           1 :   PetscCheck(fix>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value, must be >0");
     390           1 :   data->fix = fix;
     391           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     392             : }
     393             : 
     394             : /*@
     395             :    EPSJDSetFix - Sets the threshold for changing the target in the correction
     396             :    equation.
     397             : 
     398             :    Logically Collective
     399             : 
     400             :    Input Parameters:
     401             : +  eps - the eigenproblem solver context
     402             : -  fix - threshold for changing the target
     403             : 
     404             :    Options Database Key:
     405             : .  -eps_jd_fix - the fix value
     406             : 
     407             :    Note:
     408             :    The target in the correction equation is fixed at the first iterations.
     409             :    When the norm of the residual vector is lower than the fix value,
     410             :    the target is set to the corresponding eigenvalue.
     411             : 
     412             :    Level: advanced
     413             : 
     414             : .seealso: EPSJDGetFix()
     415             : @*/
     416           1 : PetscErrorCode EPSJDSetFix(EPS eps,PetscReal fix)
     417             : {
     418           1 :   PetscFunctionBegin;
     419           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     420           4 :   PetscValidLogicalCollectiveReal(eps,fix,2);
     421           1 :   PetscTryMethod(eps,"EPSJDSetFix_C",(EPS,PetscReal),(eps,fix));
     422           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     423             : }
     424             : 
     425          20 : PetscErrorCode EPSJDGetFix_JD(EPS eps,PetscReal *fix)
     426             : {
     427          20 :   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
     428             : 
     429          20 :   PetscFunctionBegin;
     430          20 :   *fix = data->fix;
     431          20 :   PetscFunctionReturn(PETSC_SUCCESS);
     432             : }
     433             : 
     434             : /*@
     435             :    EPSJDGetFix - Returns the threshold for changing the target in the correction
     436             :    equation.
     437             : 
     438             :    Not Collective
     439             : 
     440             :    Input Parameter:
     441             : .  eps - the eigenproblem solver context
     442             : 
     443             :    Output Parameter:
     444             : .  fix - threshold for changing the target
     445             : 
     446             :    Note:
     447             :    The target in the correction equation is fixed at the first iterations.
     448             :    When the norm of the residual vector is lower than the fix value,
     449             :    the target is set to the corresponding eigenvalue.
     450             : 
     451             :    Level: advanced
     452             : 
     453             : .seealso: EPSJDSetFix()
     454             : @*/
     455          20 : PetscErrorCode EPSJDGetFix(EPS eps,PetscReal *fix)
     456             : {
     457          20 :   PetscFunctionBegin;
     458          20 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     459          20 :   PetscAssertPointer(fix,2);
     460          20 :   PetscUseMethod(eps,"EPSJDGetFix_C",(EPS,PetscReal*),(eps,fix));
     461          20 :   PetscFunctionReturn(PETSC_SUCCESS);
     462             : }
     463             : 
     464           1 : static PetscErrorCode EPSJDSetConstCorrectionTol_JD(EPS eps,PetscBool constant)
     465             : {
     466           1 :   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
     467             : 
     468           1 :   PetscFunctionBegin;
     469           1 :   data->dynamic = PetscNot(constant);
     470           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     471             : }
     472             : 
     473             : /*@
     474             :    EPSJDSetConstCorrectionTol - If true, deactivates the dynamic stopping criterion
     475             :    (also called Newton) that sets the KSP relative tolerance
     476             :    to 0.5**i, where i is the number of EPS iterations from the last converged value.
     477             : 
     478             :    Logically Collective
     479             : 
     480             :    Input Parameters:
     481             : +  eps - the eigenproblem solver context
     482             : -  constant - if false, the KSP relative tolerance is set to 0.5**i.
     483             : 
     484             :    Options Database Key:
     485             : .  -eps_jd_const_correction_tol - Deactivates the dynamic stopping criterion.
     486             : 
     487             :    Level: advanced
     488             : 
     489             : .seealso: EPSJDGetConstCorrectionTol()
     490             : @*/
     491           1 : PetscErrorCode EPSJDSetConstCorrectionTol(EPS eps,PetscBool constant)
     492             : {
     493           1 :   PetscFunctionBegin;
     494           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     495           4 :   PetscValidLogicalCollectiveBool(eps,constant,2);
     496           1 :   PetscTryMethod(eps,"EPSJDSetConstCorrectionTol_C",(EPS,PetscBool),(eps,constant));
     497           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     498             : }
     499             : 
     500          20 : PetscErrorCode EPSJDGetConstCorrectionTol_JD(EPS eps,PetscBool *constant)
     501             : {
     502          20 :   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
     503             : 
     504          20 :   PetscFunctionBegin;
     505          20 :   *constant = PetscNot(data->dynamic);
     506          20 :   PetscFunctionReturn(PETSC_SUCCESS);
     507             : }
     508             : 
     509             : /*@
     510             :    EPSJDGetConstCorrectionTol - Returns a flag indicating if the dynamic stopping is being used for
     511             :    solving the correction equation.
     512             : 
     513             :    Not Collective
     514             : 
     515             :    Input Parameter:
     516             : .  eps - the eigenproblem solver context
     517             : 
     518             :    Output Parameter:
     519             : .  constant - boolean flag indicating if the dynamic stopping criterion is not being used.
     520             : 
     521             :    Notes:
     522             :    If the flag is false the KSP relative tolerance is set to 0.5**i, where i is the number
     523             :    of EPS iterations from the last converged value.
     524             : 
     525             :    Level: advanced
     526             : 
     527             : .seealso: EPSJDSetConstCorrectionTol()
     528             : @*/
     529          20 : PetscErrorCode EPSJDGetConstCorrectionTol(EPS eps,PetscBool *constant)
     530             : {
     531          20 :   PetscFunctionBegin;
     532          20 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     533          20 :   PetscAssertPointer(constant,2);
     534          20 :   PetscUseMethod(eps,"EPSJDGetConstCorrectionTol_C",(EPS,PetscBool*),(eps,constant));
     535          20 :   PetscFunctionReturn(PETSC_SUCCESS);
     536             : }
     537             : 
     538             : /*@
     539             :    EPSJDSetBOrth - Selects the orthogonalization that will be used in the search
     540             :    subspace in case of generalized Hermitian problems.
     541             : 
     542             :    Logically Collective
     543             : 
     544             :    Input Parameters:
     545             : +  eps   - the eigenproblem solver context
     546             : -  borth - whether to B-orthogonalize the search subspace
     547             : 
     548             :    Options Database Key:
     549             : .  -eps_jd_borth - Set the orthogonalization used in the search subspace
     550             : 
     551             :    Level: advanced
     552             : 
     553             : .seealso: EPSJDGetBOrth()
     554             : @*/
     555           1 : PetscErrorCode EPSJDSetBOrth(EPS eps,PetscBool borth)
     556             : {
     557           1 :   PetscFunctionBegin;
     558           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     559           4 :   PetscValidLogicalCollectiveBool(eps,borth,2);
     560           1 :   PetscTryMethod(eps,"EPSJDSetBOrth_C",(EPS,PetscBool),(eps,borth));
     561           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     562             : }
     563             : 
     564             : /*@
     565             :    EPSJDGetBOrth - Returns the orthogonalization used in the search
     566             :    subspace in case of generalized Hermitian problems.
     567             : 
     568             :    Not Collective
     569             : 
     570             :    Input Parameter:
     571             : .  eps - the eigenproblem solver context
     572             : 
     573             :    Output Parameters:
     574             : .  borth - whether to B-orthogonalize the search subspace
     575             : 
     576             :    Level: advanced
     577             : 
     578             : .seealso: EPSJDSetBOrth()
     579             : @*/
     580          20 : PetscErrorCode EPSJDGetBOrth(EPS eps,PetscBool *borth)
     581             : {
     582          20 :   PetscFunctionBegin;
     583          20 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     584          20 :   PetscAssertPointer(borth,2);
     585          20 :   PetscUseMethod(eps,"EPSJDGetBOrth_C",(EPS,PetscBool*),(eps,borth));
     586          20 :   PetscFunctionReturn(PETSC_SUCCESS);
     587             : }
     588             : 
     589          21 : SLEPC_EXTERN PetscErrorCode EPSCreate_JD(EPS eps)
     590             : {
     591          21 :   EPS_DAVIDSON   *data;
     592             : 
     593          21 :   PetscFunctionBegin;
     594          21 :   PetscCall(PetscNew(&data));
     595          21 :   eps->data = (void*)data;
     596             : 
     597          21 :   data->blocksize   = 1;
     598          21 :   data->initialsize = 0;
     599          21 :   data->minv        = 0;
     600          21 :   data->plusk       = PETSC_DEFAULT;
     601          21 :   data->ipB         = PETSC_TRUE;
     602          21 :   data->fix         = 0.01;
     603          21 :   data->krylovstart = PETSC_FALSE;
     604          21 :   data->dynamic     = PETSC_FALSE;
     605             : 
     606          21 :   eps->useds = PETSC_TRUE;
     607          21 :   eps->categ = EPS_CATEGORY_PRECOND;
     608             : 
     609          21 :   eps->ops->solve          = EPSSolve_XD;
     610          21 :   eps->ops->setup          = EPSSetUp_JD;
     611          21 :   eps->ops->setupsort      = EPSSetUpSort_Default;
     612          21 :   eps->ops->setfromoptions = EPSSetFromOptions_JD;
     613          21 :   eps->ops->destroy        = EPSDestroy_JD;
     614          21 :   eps->ops->reset          = EPSReset_XD;
     615          21 :   eps->ops->view           = EPSView_JD;
     616          21 :   eps->ops->backtransform  = EPSBackTransform_Default;
     617          21 :   eps->ops->computevectors = EPSComputeVectors_XD;
     618          21 :   eps->ops->setdefaultst   = EPSSetDefaultST_JD;
     619             : 
     620          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetKrylovStart_C",EPSXDSetKrylovStart_XD));
     621          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetKrylovStart_C",EPSXDGetKrylovStart_XD));
     622          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetBlockSize_C",EPSXDSetBlockSize_XD));
     623          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetBlockSize_C",EPSXDGetBlockSize_XD));
     624          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetRestart_C",EPSXDSetRestart_XD));
     625          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetRestart_C",EPSXDGetRestart_XD));
     626          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetInitialSize_C",EPSXDSetInitialSize_XD));
     627          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetInitialSize_C",EPSXDGetInitialSize_XD));
     628          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetFix_C",EPSJDSetFix_JD));
     629          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetFix_C",EPSJDGetFix_JD));
     630          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetConstCorrectionTol_C",EPSJDSetConstCorrectionTol_JD));
     631          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetConstCorrectionTol_C",EPSJDGetConstCorrectionTol_JD));
     632          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDSetBOrth_C",EPSXDSetBOrth_XD));
     633          21 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSJDGetBOrth_C",EPSXDGetBOrth_XD));
     634          21 :   PetscFunctionReturn(PETSC_SUCCESS);
     635             : }

Generated by: LCOV version 1.14