xref: /petsc/src/tao/bound/impls/bncg/bncg.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1ac9112b8SAlp Dener #include <petsctaolinesearch.h>
2414d97d3SAlp Dener #include <../src/tao/bound/impls/bncg/bncg.h> /*I "petsctao.h" I*/
350b47da0SAdam Denchfield #include <petscksp.h>
4ac9112b8SAlp Dener 
5c8bcdf1eSAdam Denchfield #define CG_GradientDescent    0
6c8bcdf1eSAdam Denchfield #define CG_HestenesStiefel    1
7c8bcdf1eSAdam Denchfield #define CG_FletcherReeves     2
850b47da0SAdam Denchfield #define CG_PolakRibierePolyak 3
9c8bcdf1eSAdam Denchfield #define CG_PolakRibierePlus   4
10c8bcdf1eSAdam Denchfield #define CG_DaiYuan            5
11c8bcdf1eSAdam Denchfield #define CG_HagerZhang         6
12c8bcdf1eSAdam Denchfield #define CG_DaiKou             7
13c8bcdf1eSAdam Denchfield #define CG_KouDai             8
14c8bcdf1eSAdam Denchfield #define CG_SSML_BFGS          9
15c8bcdf1eSAdam Denchfield #define CG_SSML_DFP           10
16c8bcdf1eSAdam Denchfield #define CG_SSML_BROYDEN       11
17484c7b14SAdam Denchfield #define CG_PCGradientDescent  12
18484c7b14SAdam Denchfield #define CGTypes               13
19ac9112b8SAlp Dener 
20484c7b14SAdam Denchfield static const char *CG_Table[64] = {"gd", "hs", "fr", "pr", "prp", "dy", "hz", "dk", "kd", "ssml_bfgs", "ssml_dfp", "ssml_brdn", "pcgd"};
21ac9112b8SAlp Dener 
2261be54a6SAlp Dener #define CG_AS_NONE      0
2361be54a6SAlp Dener #define CG_AS_BERTSEKAS 1
2461be54a6SAlp Dener #define CG_AS_SIZE      2
25ac9112b8SAlp Dener 
2661be54a6SAlp Dener static const char *CG_AS_TYPE[64] = {"none", "bertsekas"};
27ac9112b8SAlp Dener 
289371c9d4SSatish Balay PetscErrorCode TaoBNCGEstimateActiveSet(Tao tao, PetscInt asType) {
2961be54a6SAlp Dener   TAO_BNCG *cg = (TAO_BNCG *)tao->data;
3061be54a6SAlp Dener 
3161be54a6SAlp Dener   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->inactive_old));
3361be54a6SAlp Dener   if (cg->inactive_idx) {
349566063dSJacob Faibussowitsch     PetscCall(ISDuplicate(cg->inactive_idx, &cg->inactive_old));
359566063dSJacob Faibussowitsch     PetscCall(ISCopy(cg->inactive_idx, cg->inactive_old));
3661be54a6SAlp Dener   }
3761be54a6SAlp Dener   switch (asType) {
3861be54a6SAlp Dener   case CG_AS_NONE:
399566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cg->inactive_idx));
409566063dSJacob Faibussowitsch     PetscCall(VecWhichInactive(tao->XL, tao->solution, cg->unprojected_gradient, tao->XU, PETSC_TRUE, &cg->inactive_idx));
419566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cg->active_idx));
429566063dSJacob Faibussowitsch     PetscCall(ISComplementVec(cg->inactive_idx, tao->solution, &cg->active_idx));
4361be54a6SAlp Dener     break;
4461be54a6SAlp Dener 
4561be54a6SAlp Dener   case CG_AS_BERTSEKAS:
4661be54a6SAlp Dener     /* Use gradient descent to estimate the active set */
479566063dSJacob Faibussowitsch     PetscCall(VecCopy(cg->unprojected_gradient, cg->W));
489566063dSJacob Faibussowitsch     PetscCall(VecScale(cg->W, -1.0));
499371c9d4SSatish Balay     PetscCall(TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, cg->unprojected_gradient, cg->W, cg->work, cg->as_step, &cg->as_tol, &cg->active_lower, &cg->active_upper, &cg->active_fixed, &cg->active_idx, &cg->inactive_idx));
50c4b75bccSAlp Dener     break;
519371c9d4SSatish Balay   default: break;
5261be54a6SAlp Dener   }
5361be54a6SAlp Dener   PetscFunctionReturn(0);
5461be54a6SAlp Dener }
5561be54a6SAlp Dener 
569371c9d4SSatish Balay PetscErrorCode TaoBNCGBoundStep(Tao tao, PetscInt asType, Vec step) {
5761be54a6SAlp Dener   TAO_BNCG *cg = (TAO_BNCG *)tao->data;
5861be54a6SAlp Dener 
5961be54a6SAlp Dener   PetscFunctionBegin;
60a1318120SAlp Dener   switch (asType) {
619371c9d4SSatish Balay   case CG_AS_NONE: PetscCall(VecISSet(step, cg->active_idx, 0.0)); break;
6261be54a6SAlp Dener 
639371c9d4SSatish Balay   case CG_AS_BERTSEKAS: PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, cg->active_lower, cg->active_upper, cg->active_fixed, 1.0, step)); break;
6461be54a6SAlp Dener 
659371c9d4SSatish Balay   default: break;
6661be54a6SAlp Dener   }
6761be54a6SAlp Dener   PetscFunctionReturn(0);
6861be54a6SAlp Dener }
6961be54a6SAlp Dener 
709371c9d4SSatish Balay static PetscErrorCode TaoSolve_BNCG(Tao tao) {
71ac9112b8SAlp Dener   TAO_BNCG *cg   = (TAO_BNCG *)tao->data;
72484c7b14SAdam Denchfield   PetscReal step = 1.0, gnorm, gnorm2, resnorm;
73c4b75bccSAlp Dener   PetscInt  nDiff;
74ac9112b8SAlp Dener 
75ac9112b8SAlp Dener   PetscFunctionBegin;
76ac9112b8SAlp Dener   /*   Project the current point onto the feasible set */
779566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
789566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
79ac9112b8SAlp Dener 
80c8bcdf1eSAdam Denchfield   /* Project the initial point onto the feasible region */
819566063dSJacob Faibussowitsch   PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
82484c7b14SAdam Denchfield 
8348a46eb9SPierre Jolivet   if (nDiff > 0 || !tao->recycle) PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient));
849566063dSJacob Faibussowitsch   PetscCall(VecNorm(cg->unprojected_gradient, NORM_2, &gnorm));
853c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(cg->f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
86ac9112b8SAlp Dener 
8761be54a6SAlp Dener   /* Estimate the active set and compute the projected gradient */
889566063dSJacob Faibussowitsch   PetscCall(TaoBNCGEstimateActiveSet(tao, cg->as_type));
8961be54a6SAlp Dener 
90ac9112b8SAlp Dener   /* Project the gradient and calculate the norm */
919566063dSJacob Faibussowitsch   PetscCall(VecCopy(cg->unprojected_gradient, tao->gradient));
929566063dSJacob Faibussowitsch   PetscCall(VecISSet(tao->gradient, cg->active_idx, 0.0));
939566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
94ac9112b8SAlp Dener   gnorm2 = gnorm * gnorm;
95ac9112b8SAlp Dener 
96c8bcdf1eSAdam Denchfield   /* Initialize counters */
97e031d6f5SAlp Dener   tao->niter   = 0;
9850b47da0SAdam Denchfield   cg->ls_fails = cg->descent_error = 0;
99c8bcdf1eSAdam Denchfield   cg->resets                       = -1;
100484c7b14SAdam Denchfield   cg->skipped_updates = cg->pure_gd_steps = 0;
101c8bcdf1eSAdam Denchfield   cg->iter_quad                           = 0;
102c8bcdf1eSAdam Denchfield 
103c8bcdf1eSAdam Denchfield   /* Convergence test at the starting point. */
104ac9112b8SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
105484c7b14SAdam Denchfield 
1069566063dSJacob Faibussowitsch   PetscCall(VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W));
1079566063dSJacob Faibussowitsch   PetscCall(VecNorm(cg->W, NORM_2, &resnorm));
1083c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
1099566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its));
1109566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step));
111dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
112ac9112b8SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
113484c7b14SAdam Denchfield   /* Calculate initial direction. */
114414d97d3SAlp Dener   if (!tao->recycle) {
115484c7b14SAdam Denchfield     /* We are not recycling a solution/history from a past TaoSolve */
1169566063dSJacob Faibussowitsch     PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
117ac9112b8SAlp Dener   }
118c8bcdf1eSAdam Denchfield   /* Initial gradient descent step. Scaling by 1.0 also does a decent job for some problems. */
119c8bcdf1eSAdam Denchfield   while (1) {
120e1e80dc8SAlp Dener     /* Call general purpose update function */
121e1e80dc8SAlp Dener     if (tao->ops->update) {
122dbbe0bcdSBarry Smith       PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
1237494f0b1SStefano Zampini       PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient));
124e1e80dc8SAlp Dener     }
1259566063dSJacob Faibussowitsch     PetscCall(TaoBNCGConductIteration(tao, gnorm));
126c8bcdf1eSAdam Denchfield     if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
127ac9112b8SAlp Dener   }
128ac9112b8SAlp Dener }
129ac9112b8SAlp Dener 
1309371c9d4SSatish Balay static PetscErrorCode TaoSetUp_BNCG(Tao tao) {
131ac9112b8SAlp Dener   TAO_BNCG *cg = (TAO_BNCG *)tao->data;
132ac9112b8SAlp Dener 
133ac9112b8SAlp Dener   PetscFunctionBegin;
13448a46eb9SPierre Jolivet   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
13548a46eb9SPierre Jolivet   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
13648a46eb9SPierre Jolivet   if (!cg->W) PetscCall(VecDuplicate(tao->solution, &cg->W));
13748a46eb9SPierre Jolivet   if (!cg->work) PetscCall(VecDuplicate(tao->solution, &cg->work));
13848a46eb9SPierre Jolivet   if (!cg->sk) PetscCall(VecDuplicate(tao->solution, &cg->sk));
13948a46eb9SPierre Jolivet   if (!cg->yk) PetscCall(VecDuplicate(tao->gradient, &cg->yk));
14048a46eb9SPierre Jolivet   if (!cg->X_old) PetscCall(VecDuplicate(tao->solution, &cg->X_old));
14148a46eb9SPierre Jolivet   if (!cg->G_old) PetscCall(VecDuplicate(tao->gradient, &cg->G_old));
142c8bcdf1eSAdam Denchfield   if (cg->diag_scaling) {
1439566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &cg->d_work));
1449566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &cg->y_work));
1459566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &cg->g_work));
146c4b75bccSAlp Dener   }
14748a46eb9SPierre Jolivet   if (!cg->unprojected_gradient) PetscCall(VecDuplicate(tao->gradient, &cg->unprojected_gradient));
14848a46eb9SPierre Jolivet   if (!cg->unprojected_gradient_old) PetscCall(VecDuplicate(tao->gradient, &cg->unprojected_gradient_old));
1499566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(cg->B, cg->sk, cg->yk));
1501baa6e33SBarry Smith   if (cg->pc) PetscCall(MatLMVMSetJ0(cg->B, cg->pc));
151ac9112b8SAlp Dener   PetscFunctionReturn(0);
152ac9112b8SAlp Dener }
153ac9112b8SAlp Dener 
1549371c9d4SSatish Balay static PetscErrorCode TaoDestroy_BNCG(Tao tao) {
155ac9112b8SAlp Dener   TAO_BNCG *cg = (TAO_BNCG *)tao->data;
156ac9112b8SAlp Dener 
157ac9112b8SAlp Dener   PetscFunctionBegin;
158ac9112b8SAlp Dener   if (tao->setupcalled) {
1599566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->W));
1609566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->work));
1619566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->X_old));
1629566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->G_old));
1639566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->unprojected_gradient));
1649566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->unprojected_gradient_old));
1659566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->g_work));
1669566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->d_work));
1679566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->y_work));
1689566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->sk));
1699566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->yk));
170ac9112b8SAlp Dener   }
1719566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->active_lower));
1729566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->active_upper));
1739566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->active_fixed));
1749566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->active_idx));
1759566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->inactive_idx));
1769566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->inactive_old));
1779566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->new_inactives));
1789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&cg->B));
17948a46eb9SPierre Jolivet   if (cg->pc) PetscCall(MatDestroy(&cg->pc));
1809566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
181ac9112b8SAlp Dener   PetscFunctionReturn(0);
182ac9112b8SAlp Dener }
183ac9112b8SAlp Dener 
1849371c9d4SSatish Balay static PetscErrorCode TaoSetFromOptions_BNCG(Tao tao, PetscOptionItems *PetscOptionsObject) {
185ac9112b8SAlp Dener   TAO_BNCG *cg = (TAO_BNCG *)tao->data;
186ac9112b8SAlp Dener 
187ac9112b8SAlp Dener   PetscFunctionBegin;
188d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Nonlinear Conjugate Gradient method for unconstrained optimization");
1899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_bncg_type", "cg formula", "", CG_Table, CGTypes, CG_Table[cg->cg_type], &cg->cg_type, NULL));
1908ebe3e4eSStefano Zampini   if (cg->cg_type != CG_SSML_BFGS) cg->alpha = -1.0; /* Setting defaults for non-BFGS methods. User can change it below. */
191484c7b14SAdam Denchfield   if (CG_GradientDescent == cg->cg_type) {
192484c7b14SAdam Denchfield     cg->cg_type          = CG_PCGradientDescent;
193484c7b14SAdam Denchfield     /* Set scaling equal to none or, at best, scalar scaling. */
194484c7b14SAdam Denchfield     cg->unscaled_restart = PETSC_TRUE;
195484c7b14SAdam Denchfield     cg->diag_scaling     = PETSC_FALSE;
196484c7b14SAdam Denchfield   }
1979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_bncg_as_type", "active set estimation method", "", CG_AS_TYPE, CG_AS_SIZE, CG_AS_TYPE[cg->cg_type], &cg->cg_type, NULL));
1989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_hz_eta", "(developer) cutoff tolerance for HZ", "", cg->hz_eta, &cg->hz_eta, NULL));
1999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_eps", "(developer) cutoff value for restarts", "", cg->epsilon, &cg->epsilon, NULL));
2009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_dk_eta", "(developer) cutoff tolerance for DK", "", cg->dk_eta, &cg->dk_eta, NULL));
2019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_xi", "(developer) Parameter in the KD method", "", cg->xi, &cg->xi, NULL));
2029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_theta", "(developer) update parameter for the Broyden method", "", cg->theta, &cg->theta, NULL));
2039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_hz_theta", "(developer) parameter for the HZ (2006) method", "", cg->hz_theta, &cg->hz_theta, NULL));
2049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_alpha", "(developer) parameter for the scalar scaling", "", cg->alpha, &cg->alpha, NULL));
2059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_bfgs_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->bfgs_scale, &cg->bfgs_scale, NULL));
2069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_dfp_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->dfp_scale, &cg->dfp_scale, NULL));
2079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_bncg_diag_scaling", "Enable diagonal Broyden-like preconditioning", "", cg->diag_scaling, &cg->diag_scaling, NULL));
2089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_bncg_dynamic_restart", "(developer) use dynamic restarts as in HZ, DK, KD", "", cg->use_dynamic_restart, &cg->use_dynamic_restart, NULL));
2099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_bncg_unscaled_restart", "(developer) use unscaled gradient restarts", "", cg->unscaled_restart, &cg->unscaled_restart, NULL));
2109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_zeta", "(developer) Free parameter for the Kou-Dai method", "", cg->zeta, &cg->zeta, NULL));
2119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_bncg_min_quad", "(developer) Number of iterations with approximate quadratic behavior needed for restart", "", cg->min_quad, &cg->min_quad, NULL));
2129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_bncg_min_restart_num", "(developer) Number of iterations between restarts (times dimension)", "", cg->min_restart_num, &cg->min_restart_num, NULL));
2139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_bncg_spaced_restart", "(developer) Enable regular steepest descent restarting every fixed number of iterations", "", cg->spaced_restart, &cg->spaced_restart, NULL));
2149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_bncg_no_scaling", "Disable all scaling except in restarts", "", cg->no_scaling, &cg->no_scaling, NULL));
215484c7b14SAdam Denchfield   if (cg->no_scaling) {
216484c7b14SAdam Denchfield     cg->diag_scaling = PETSC_FALSE;
217484c7b14SAdam Denchfield     cg->alpha        = -1.0;
218484c7b14SAdam Denchfield   }
219b474139fSKarl Rupp   if (cg->alpha == -1.0 && cg->cg_type == CG_KouDai && !cg->diag_scaling) { /* Some more default options that appear to be good. */
220484c7b14SAdam Denchfield     cg->neg_xi = PETSC_TRUE;
221484c7b14SAdam Denchfield   }
2229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_bncg_neg_xi", "(developer) Use negative xi when it might be a smaller descent direction than necessary", "", cg->neg_xi, &cg->neg_xi, NULL));
2239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", cg->as_tol, &cg->as_tol, NULL));
2249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_as_step", "(developer) step length used when estimating actively bounded variables", "", cg->as_step, &cg->as_step, NULL));
2259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_delta_min", "(developer) minimum scaling factor used for scaled gradient restarts", "", cg->delta_min, &cg->delta_min, NULL));
2269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_delta_max", "(developer) maximum scaling factor used for scaled gradient restarts", "", cg->delta_max, &cg->delta_max, NULL));
22750b47da0SAdam Denchfield 
228d0609cedSBarry Smith   PetscOptionsHeadEnd();
2299566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(cg->B, ((PetscObject)tao)->prefix));
2309566063dSJacob Faibussowitsch   PetscCall(MatAppendOptionsPrefix(cg->B, "tao_bncg_"));
2319566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(cg->B));
232ac9112b8SAlp Dener   PetscFunctionReturn(0);
233ac9112b8SAlp Dener }
234ac9112b8SAlp Dener 
2359371c9d4SSatish Balay static PetscErrorCode TaoView_BNCG(Tao tao, PetscViewer viewer) {
236ac9112b8SAlp Dener   PetscBool isascii;
237ac9112b8SAlp Dener   TAO_BNCG *cg = (TAO_BNCG *)tao->data;
238ac9112b8SAlp Dener 
239ac9112b8SAlp Dener   PetscFunctionBegin;
2409566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
241ac9112b8SAlp Dener   if (isascii) {
2429566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
2439566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cg->cg_type]));
24463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Skipped Stepdirection Updates: %" PetscInt_FMT "\n", cg->skipped_updates));
24563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", cg->resets));
24663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Pure gradient steps: %" PetscInt_FMT "\n", cg->pure_gd_steps));
24763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Not a descent direction: %" PetscInt_FMT "\n", cg->descent_error));
24863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Line search fails: %" PetscInt_FMT "\n", cg->ls_fails));
249484c7b14SAdam Denchfield     if (cg->diag_scaling) {
2509566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
251484c7b14SAdam Denchfield       if (isascii) {
2529566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
2539566063dSJacob Faibussowitsch         PetscCall(MatView(cg->B, viewer));
2549566063dSJacob Faibussowitsch         PetscCall(PetscViewerPopFormat(viewer));
255484c7b14SAdam Denchfield       }
256484c7b14SAdam Denchfield     }
2579566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
258ac9112b8SAlp Dener   }
259ac9112b8SAlp Dener   PetscFunctionReturn(0);
260ac9112b8SAlp Dener }
261ac9112b8SAlp Dener 
2629371c9d4SSatish Balay PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal yty, PetscReal yts, PetscReal sts, PetscReal *scale, PetscReal alpha) {
263c8bcdf1eSAdam Denchfield   PetscReal a, b, c, sig1, sig2;
264c8bcdf1eSAdam Denchfield 
265c8bcdf1eSAdam Denchfield   PetscFunctionBegin;
266c8bcdf1eSAdam Denchfield   *scale = 0.0;
2678ebe3e4eSStefano Zampini   if (1.0 == alpha) *scale = yts / yty;
2688ebe3e4eSStefano Zampini   else if (0.0 == alpha) *scale = sts / yts;
26950b47da0SAdam Denchfield   else if (-1.0 == alpha) *scale = 1.0;
270c8bcdf1eSAdam Denchfield   else {
271c8bcdf1eSAdam Denchfield     a = yty;
272c8bcdf1eSAdam Denchfield     b = yts;
273c8bcdf1eSAdam Denchfield     c = sts;
274c8bcdf1eSAdam Denchfield     a *= alpha;
275c8bcdf1eSAdam Denchfield     b *= -(2.0 * alpha - 1.0);
276c8bcdf1eSAdam Denchfield     c *= alpha - 1.0;
277c8bcdf1eSAdam Denchfield     sig1 = (-b + PetscSqrtReal(b * b - 4.0 * a * c)) / (2.0 * a);
278c8bcdf1eSAdam Denchfield     sig2 = (-b - PetscSqrtReal(b * b - 4.0 * a * c)) / (2.0 * a);
279c8bcdf1eSAdam Denchfield     /* accept the positive root as the scalar */
2808ebe3e4eSStefano Zampini     if (sig1 > 0.0) *scale = sig1;
2818ebe3e4eSStefano Zampini     else if (sig2 > 0.0) *scale = sig2;
2828ebe3e4eSStefano Zampini     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Cannot find positive scalar");
283c8bcdf1eSAdam Denchfield   }
284c8bcdf1eSAdam Denchfield   PetscFunctionReturn(0);
285c8bcdf1eSAdam Denchfield }
286c8bcdf1eSAdam Denchfield 
287ac9112b8SAlp Dener /*MC
288ac9112b8SAlp Dener   TAOBNCG - Bound-constrained Nonlinear Conjugate Gradient method.
289ac9112b8SAlp Dener 
290ac9112b8SAlp Dener   Options Database Keys:
29150b47da0SAdam Denchfield + -tao_bncg_recycle - enable recycling the latest calculated gradient vector in subsequent TaoSolve() calls (currently disabled)
292c4b75bccSAlp Dener . -tao_bncg_eta <r> - restart tolerance
29361be54a6SAlp Dener . -tao_bncg_type <taocg_type> - cg formula
294c4b75bccSAlp Dener . -tao_bncg_as_type <none,bertsekas> - active set estimation method
295c4b75bccSAlp Dener . -tao_bncg_as_tol <r> - tolerance used in Bertsekas active-set estimation
296c4b75bccSAlp Dener . -tao_bncg_as_step <r> - trial step length used in Bertsekas active-set estimation
29750b47da0SAdam Denchfield . -tao_bncg_eps <r> - cutoff used for determining whether or not we restart based on steplength each iteration, as well as determining whether or not we continue using the last stepdirection. Defaults to machine precision.
29850b47da0SAdam Denchfield . -tao_bncg_theta <r> - convex combination parameter for the Broyden method
29950b47da0SAdam Denchfield . -tao_bncg_hz_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods
30050b47da0SAdam Denchfield . -tao_bncg_dk_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods
30150b47da0SAdam Denchfield . -tao_bncg_xi <r> - Multiplicative constant of the gamma term in the KD method
30250b47da0SAdam Denchfield . -tao_bncg_hz_theta <r> - Multiplicative constant of the theta term for the HZ method
30350b47da0SAdam Denchfield . -tao_bncg_bfgs_scale <r> - Scaling parameter of the bfgs contribution to the scalar Broyden method
30450b47da0SAdam Denchfield . -tao_bncg_dfp_scale <r> - Scaling parameter of the dfp contribution to the scalar Broyden method
30550b47da0SAdam Denchfield . -tao_bncg_diag_scaling <b> - Whether or not to use diagonal initialization/preconditioning for the CG methods. Default True.
30650b47da0SAdam Denchfield . -tao_bncg_dynamic_restart <b> - use dynamic restart strategy in the HZ, DK, KD methods
30750b47da0SAdam Denchfield . -tao_bncg_unscaled_restart <b> - whether or not to scale the gradient when doing gradient descent restarts
30850b47da0SAdam Denchfield . -tao_bncg_zeta <r> - Scaling parameter in the KD method
309484c7b14SAdam Denchfield . -tao_bncg_delta_min <r> - Minimum bound for rescaling during restarted gradient descent steps
310484c7b14SAdam Denchfield . -tao_bncg_delta_max <r> - Maximum bound for rescaling during restarted gradient descent steps
31150b47da0SAdam Denchfield . -tao_bncg_min_quad <i> - Number of quadratic-like steps in a row necessary to do a dynamic restart
31250b47da0SAdam Denchfield . -tao_bncg_min_restart_num <i> - This number, x, makes sure there is a gradient descent step every x*n iterations, where n is the dimension of the problem
31350b47da0SAdam Denchfield . -tao_bncg_spaced_restart <b> - whether or not to do gradient descent steps every x*n iterations
314484c7b14SAdam Denchfield . -tao_bncg_no_scaling <b> - If true, eliminates all scaling, including defaults.
3153850be85SAlp Dener - -tao_bncg_neg_xi <b> - Whether or not to use negative xi in the KD method under certain conditions
316ac9112b8SAlp Dener 
317ac9112b8SAlp Dener   Notes:
318ac9112b8SAlp Dener     CG formulas are:
3193850be85SAlp Dener + "gd" - Gradient Descent
3203850be85SAlp Dener . "fr" - Fletcher-Reeves
3213850be85SAlp Dener . "pr" - Polak-Ribiere-Polyak
3223850be85SAlp Dener . "prp" - Polak-Ribiere-Plus
3233850be85SAlp Dener . "hs" - Hestenes-Steifel
3243850be85SAlp Dener . "dy" - Dai-Yuan
3253850be85SAlp Dener . "ssml_bfgs" - Self-Scaling Memoryless BFGS
3263850be85SAlp Dener . "ssml_dfp"  - Self-Scaling Memoryless DFP
3273850be85SAlp Dener . "ssml_brdn" - Self-Scaling Memoryless Broyden
3283850be85SAlp Dener . "hz" - Hager-Zhang (CG_DESCENT 5.3)
3293850be85SAlp Dener . "dk" - Dai-Kou (2013)
3303850be85SAlp Dener - "kd" - Kou-Dai (2015)
3319abc5736SPatrick Sanan 
332ac9112b8SAlp Dener   Level: beginner
333ac9112b8SAlp Dener M*/
334ac9112b8SAlp Dener 
3359371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_BNCG(Tao tao) {
336ac9112b8SAlp Dener   TAO_BNCG   *cg;
337ac9112b8SAlp Dener   const char *morethuente_type = TAOLINESEARCHMT;
338ac9112b8SAlp Dener 
339ac9112b8SAlp Dener   PetscFunctionBegin;
340ac9112b8SAlp Dener   tao->ops->setup          = TaoSetUp_BNCG;
341ac9112b8SAlp Dener   tao->ops->solve          = TaoSolve_BNCG;
342ac9112b8SAlp Dener   tao->ops->view           = TaoView_BNCG;
343ac9112b8SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNCG;
344ac9112b8SAlp Dener   tao->ops->destroy        = TaoDestroy_BNCG;
345ac9112b8SAlp Dener 
346ac9112b8SAlp Dener   /* Override default settings (unless already changed) */
347ac9112b8SAlp Dener   if (!tao->max_it_changed) tao->max_it = 2000;
348ac9112b8SAlp Dener   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
349ac9112b8SAlp Dener 
350ac9112b8SAlp Dener   /*  Note: nondefault values should be used for nonlinear conjugate gradient  */
351ac9112b8SAlp Dener   /*  method.  In particular, gtol should be less that 0.5; the value used in  */
352ac9112b8SAlp Dener   /*  Nocedal and Wright is 0.10.  We use the default values for the  */
353ac9112b8SAlp Dener   /*  linesearch because it seems to work better. */
3549566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
3559566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
3569566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
3579566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
358ac9112b8SAlp Dener 
359*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&cg));
360ac9112b8SAlp Dener   tao->data = (void *)cg;
3619566063dSJacob Faibussowitsch   PetscCall(KSPInitializePackage());
3629566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &cg->B));
3639566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)cg->B, (PetscObject)tao, 1));
3649566063dSJacob Faibussowitsch   PetscCall(MatSetType(cg->B, MATLMVMDIAGBROYDEN));
36550b47da0SAdam Denchfield 
366484c7b14SAdam Denchfield   cg->pc = NULL;
367484c7b14SAdam Denchfield 
36850b47da0SAdam Denchfield   cg->dk_eta           = 0.5;
36950b47da0SAdam Denchfield   cg->hz_eta           = 0.4;
370c8bcdf1eSAdam Denchfield   cg->dynamic_restart  = PETSC_FALSE;
371c8bcdf1eSAdam Denchfield   cg->unscaled_restart = PETSC_FALSE;
372484c7b14SAdam Denchfield   cg->no_scaling       = PETSC_FALSE;
373484c7b14SAdam Denchfield   cg->delta_min        = 1e-7;
374484c7b14SAdam Denchfield   cg->delta_max        = 100;
375c8bcdf1eSAdam Denchfield   cg->theta            = 1.0;
376c8bcdf1eSAdam Denchfield   cg->hz_theta         = 1.0;
377c8bcdf1eSAdam Denchfield   cg->dfp_scale        = 1.0;
378c8bcdf1eSAdam Denchfield   cg->bfgs_scale       = 1.0;
37950b47da0SAdam Denchfield   cg->zeta             = 0.1;
38050b47da0SAdam Denchfield   cg->min_quad         = 6;
381c8bcdf1eSAdam Denchfield   cg->min_restart_num  = 6; /* As in CG_DESCENT and KD2015*/
382c8bcdf1eSAdam Denchfield   cg->xi               = 1.0;
38350b47da0SAdam Denchfield   cg->neg_xi           = PETSC_TRUE;
384c8bcdf1eSAdam Denchfield   cg->spaced_restart   = PETSC_FALSE;
385c8bcdf1eSAdam Denchfield   cg->tol_quad         = 1e-8;
38661be54a6SAlp Dener   cg->as_step          = 0.001;
38761be54a6SAlp Dener   cg->as_tol           = 0.001;
38850b47da0SAdam Denchfield   cg->eps_23           = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0); /* Just a little tighter*/
38961be54a6SAlp Dener   cg->as_type          = CG_AS_BERTSEKAS;
390c8bcdf1eSAdam Denchfield   cg->cg_type          = CG_SSML_BFGS;
391c8bcdf1eSAdam Denchfield   cg->alpha            = 1.0;
392c8bcdf1eSAdam Denchfield   cg->diag_scaling     = PETSC_TRUE;
393c8bcdf1eSAdam Denchfield   PetscFunctionReturn(0);
394c8bcdf1eSAdam Denchfield }
395c8bcdf1eSAdam Denchfield 
3969371c9d4SSatish Balay PetscErrorCode TaoBNCGResetUpdate(Tao tao, PetscReal gnormsq) {
397c8bcdf1eSAdam Denchfield   TAO_BNCG *cg = (TAO_BNCG *)tao->data;
398c8bcdf1eSAdam Denchfield   PetscReal scaling;
399c8bcdf1eSAdam Denchfield 
400c8bcdf1eSAdam Denchfield   PetscFunctionBegin;
401c8bcdf1eSAdam Denchfield   ++cg->resets;
402c8bcdf1eSAdam Denchfield   scaling = 2.0 * PetscMax(1.0, PetscAbsScalar(cg->f)) / PetscMax(gnormsq, cg->eps_23);
403484c7b14SAdam Denchfield   scaling = PetscMin(cg->delta_max, PetscMax(cg->delta_min, scaling));
404484c7b14SAdam Denchfield   if (cg->unscaled_restart) {
405484c7b14SAdam Denchfield     scaling = 1.0;
406484c7b14SAdam Denchfield     ++cg->pure_gd_steps;
407484c7b14SAdam Denchfield   }
4089566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(tao->stepdirection, -scaling, 0.0, tao->gradient));
409c8bcdf1eSAdam Denchfield   /* Also want to reset our diagonal scaling with each restart */
4101baa6e33SBarry Smith   if (cg->diag_scaling) PetscCall(MatLMVMReset(cg->B, PETSC_FALSE));
411c8bcdf1eSAdam Denchfield   PetscFunctionReturn(0);
412c8bcdf1eSAdam Denchfield }
413c8bcdf1eSAdam Denchfield 
4149371c9d4SSatish Balay PetscErrorCode TaoBNCGCheckDynamicRestart(Tao tao, PetscReal stepsize, PetscReal gd, PetscReal gd_old, PetscBool *dynrestart, PetscReal fold) {
415c8bcdf1eSAdam Denchfield   TAO_BNCG *cg = (TAO_BNCG *)tao->data;
416c8bcdf1eSAdam Denchfield   PetscReal quadinterp;
417c8bcdf1eSAdam Denchfield 
418c8bcdf1eSAdam Denchfield   PetscFunctionBegin;
41950b47da0SAdam Denchfield   if (cg->f < cg->min_quad / 10) {
42050b47da0SAdam Denchfield     *dynrestart = PETSC_FALSE;
42150b47da0SAdam Denchfield     PetscFunctionReturn(0);
42250b47da0SAdam Denchfield   } /* just skip this since this strategy doesn't work well for functions near zero */
423484c7b14SAdam Denchfield   quadinterp = 2.0 * (cg->f - fold) / (stepsize * (gd + gd_old));
42450b47da0SAdam Denchfield   if (PetscAbs(quadinterp - 1.0) < cg->tol_quad) ++cg->iter_quad;
425c8bcdf1eSAdam Denchfield   else {
426c8bcdf1eSAdam Denchfield     cg->iter_quad = 0;
427c8bcdf1eSAdam Denchfield     *dynrestart   = PETSC_FALSE;
428c8bcdf1eSAdam Denchfield   }
429c8bcdf1eSAdam Denchfield   if (cg->iter_quad >= cg->min_quad) {
430c8bcdf1eSAdam Denchfield     cg->iter_quad = 0;
431c8bcdf1eSAdam Denchfield     *dynrestart   = PETSC_TRUE;
432c8bcdf1eSAdam Denchfield   }
433c8bcdf1eSAdam Denchfield   PetscFunctionReturn(0);
434c8bcdf1eSAdam Denchfield }
435c8bcdf1eSAdam Denchfield 
4369371c9d4SSatish Balay PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao tao, PetscReal gnorm2, PetscReal step, PetscReal fold, PetscReal gnorm2_old, PetscReal dnorm, PetscBool pcgd_fallback) {
437c8bcdf1eSAdam Denchfield   TAO_BNCG *cg    = (TAO_BNCG *)tao->data;
43850b47da0SAdam Denchfield   PetscReal gamma = 1.0, tau_k, beta;
439484c7b14SAdam Denchfield   PetscReal tmp = 1.0, ynorm, ynorm2 = 1.0, snorm = 1.0, dk_yk = 1.0, gd;
44050b47da0SAdam Denchfield   PetscReal gkp1_yk, gd_old, tau_bfgs, tau_dfp, gkp1D_yk, gtDg;
441c8bcdf1eSAdam Denchfield   PetscInt  dim;
442484c7b14SAdam Denchfield   PetscBool cg_restart = PETSC_FALSE;
443c8bcdf1eSAdam Denchfield   PetscFunctionBegin;
444c8bcdf1eSAdam Denchfield 
44550b47da0SAdam Denchfield   /* Local curvature check to see if we need to restart */
446414d97d3SAlp Dener   if (tao->niter >= 1 || tao->recycle) {
4479566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient));
4489566063dSJacob Faibussowitsch     PetscCall(VecNorm(cg->yk, NORM_2, &ynorm));
449c8bcdf1eSAdam Denchfield     ynorm2 = ynorm * ynorm;
4509566063dSJacob Faibussowitsch     PetscCall(VecDot(cg->yk, tao->stepdirection, &dk_yk));
451484c7b14SAdam Denchfield     if (step * dnorm < PETSC_MACHINE_EPSILON || step * dk_yk < PETSC_MACHINE_EPSILON) {
452e2570530SAlp Dener       cg_restart = PETSC_TRUE;
453484c7b14SAdam Denchfield       ++cg->skipped_updates;
454484c7b14SAdam Denchfield     }
45550b47da0SAdam Denchfield     if (cg->spaced_restart) {
4569566063dSJacob Faibussowitsch       PetscCall(VecGetSize(tao->gradient, &dim));
457e2570530SAlp Dener       if (tao->niter % (dim * cg->min_restart_num)) cg_restart = PETSC_TRUE;
45850b47da0SAdam Denchfield     }
45950b47da0SAdam Denchfield   }
46050b47da0SAdam Denchfield   /* If the user wants regular restarts, do it every 6n iterations, where n=dimension */
46150b47da0SAdam Denchfield   if (cg->spaced_restart) {
4629566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->gradient, &dim));
463e2570530SAlp Dener     if (0 == tao->niter % (6 * dim)) cg_restart = PETSC_TRUE;
46450b47da0SAdam Denchfield   }
46550b47da0SAdam Denchfield   /* Compute the diagonal scaling vector if applicable */
4661baa6e33SBarry Smith   if (cg->diag_scaling) PetscCall(MatLMVMUpdate(cg->B, tao->solution, tao->gradient));
46750b47da0SAdam Denchfield 
468484c7b14SAdam Denchfield   /* A note on diagonal scaling (to be added to paper):
469484c7b14SAdam Denchfield    For the FR, PR, PRP, and DY methods, the diagonally scaled versions
470484c7b14SAdam Denchfield    must be derived as a preconditioned CG method rather than as
471484c7b14SAdam Denchfield    a Hessian initialization like in the Broyden methods. */
47250b47da0SAdam Denchfield 
473484c7b14SAdam Denchfield   /* In that case, one writes the objective function as
474484c7b14SAdam Denchfield    f(x) \equiv f(Ay). Gradient evaluations yield g(x_k) = A g(Ay_k) = A g(x_k).
475484c7b14SAdam Denchfield    Furthermore, the direction d_k \equiv (x_k - x_{k-1})/step according to
476484c7b14SAdam Denchfield    HZ (2006) becomes A^{-1} d_k, such that d_k^T g_k remains the
477484c7b14SAdam Denchfield    same under preconditioning. Note that A is diagonal, such that A^T = A. */
47850b47da0SAdam Denchfield 
479484c7b14SAdam Denchfield   /* This yields questions like what the dot product d_k^T y_k
480484c7b14SAdam Denchfield    should look like. HZ mistakenly treats that as the same under
481484c7b14SAdam Denchfield    preconditioning, but that is not necessarily true. */
48250b47da0SAdam Denchfield 
483484c7b14SAdam Denchfield   /* Observe y_k \equiv g_k - g_{k-1}, and under the P.C. transformation,
484484c7b14SAdam Denchfield    we get d_k^T y_k = (d_k^T A_k^{-T} A_k g_k - d_k^T A_k^{-T} A_{k-1} g_{k-1}),
485484c7b14SAdam Denchfield    yielding d_k^T y_k = d_k^T g_k - d_k^T (A_k^{-T} A_{k-1} g_{k-1}), which is
486484c7b14SAdam Denchfield    NOT the same if our preconditioning matrix is updated between iterations.
487484c7b14SAdam Denchfield    This same issue is found when considering dot products of the form g_{k+1}^T y_k. */
48850b47da0SAdam Denchfield 
48950b47da0SAdam Denchfield   /* Compute CG step direction */
49050b47da0SAdam Denchfield   if (cg_restart) {
4919566063dSJacob Faibussowitsch     PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
492484c7b14SAdam Denchfield   } else if (pcgd_fallback) {
493484c7b14SAdam Denchfield     /* Just like preconditioned CG */
4949566063dSJacob Faibussowitsch     PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
4959566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work));
49650b47da0SAdam Denchfield   } else if (ynorm2 > PETSC_MACHINE_EPSILON) {
49750b47da0SAdam Denchfield     switch (cg->cg_type) {
498484c7b14SAdam Denchfield     case CG_PCGradientDescent:
49950b47da0SAdam Denchfield       if (!cg->diag_scaling) {
500484c7b14SAdam Denchfield         if (!cg->no_scaling) {
50150b47da0SAdam Denchfield           cg->sts = step * step * dnorm * dnorm;
5029566063dSJacob Faibussowitsch           PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, cg->sts, &tau_k, cg->alpha));
503484c7b14SAdam Denchfield         } else {
504484c7b14SAdam Denchfield           tau_k = 1.0;
505484c7b14SAdam Denchfield           ++cg->pure_gd_steps;
506484c7b14SAdam Denchfield         }
5079566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, 0.0, tao->gradient));
50850b47da0SAdam Denchfield       } else {
5099566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5109566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work));
51150b47da0SAdam Denchfield       }
51250b47da0SAdam Denchfield       break;
513484c7b14SAdam Denchfield 
51450b47da0SAdam Denchfield     case CG_HestenesStiefel:
51550b47da0SAdam Denchfield       /* Classic Hestenes-Stiefel method, modified with scalar and diagonal preconditioning. */
51650b47da0SAdam Denchfield       if (!cg->diag_scaling) {
51750b47da0SAdam Denchfield         cg->sts = step * step * dnorm * dnorm;
5189566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
5199566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, cg->sts, &tau_k, cg->alpha));
52050b47da0SAdam Denchfield         beta = tau_k * gkp1_yk / dk_yk;
5219566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
52250b47da0SAdam Denchfield       } else {
5239566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5249566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, cg->g_work, &gkp1_yk));
52550b47da0SAdam Denchfield         beta = gkp1_yk / dk_yk;
5269566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
52750b47da0SAdam Denchfield       }
528c8bcdf1eSAdam Denchfield       break;
529484c7b14SAdam Denchfield 
530c8bcdf1eSAdam Denchfield     case CG_FletcherReeves:
5319566063dSJacob Faibussowitsch       PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old));
5329566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient));
5339566063dSJacob Faibussowitsch       PetscCall(VecNorm(cg->yk, NORM_2, &ynorm));
53450b47da0SAdam Denchfield       ynorm2 = ynorm * ynorm;
5359566063dSJacob Faibussowitsch       PetscCall(VecDot(cg->yk, tao->stepdirection, &dk_yk));
53650b47da0SAdam Denchfield       if (!cg->diag_scaling) {
5379566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, step * step * dnorm * dnorm, &tau_k, cg->alpha));
53850b47da0SAdam Denchfield         beta = tau_k * gnorm2 / gnorm2_old;
5399566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
54050b47da0SAdam Denchfield       } else {
5419566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old)); /* Before it's updated */
5429566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5439566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->gradient, cg->g_work, &tmp));
54450b47da0SAdam Denchfield         beta = tmp / gnorm2_old;
5459566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
54650b47da0SAdam Denchfield       }
547c8bcdf1eSAdam Denchfield       break;
548484c7b14SAdam Denchfield 
54950b47da0SAdam Denchfield     case CG_PolakRibierePolyak:
55050b47da0SAdam Denchfield       snorm = step * dnorm;
55150b47da0SAdam Denchfield       if (!cg->diag_scaling) {
5529566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old));
5539566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
5549566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, snorm * snorm, &tau_k, cg->alpha));
55550b47da0SAdam Denchfield         beta = tau_k * gkp1_yk / gnorm2_old;
5569566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
55750b47da0SAdam Denchfield       } else {
5589566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old));
5599566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5609566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
56150b47da0SAdam Denchfield         beta = gkp1_yk / gnorm2_old;
5629566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
56350b47da0SAdam Denchfield       }
564c8bcdf1eSAdam Denchfield       break;
565484c7b14SAdam Denchfield 
566c8bcdf1eSAdam Denchfield     case CG_PolakRibierePlus:
5679566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient));
5689566063dSJacob Faibussowitsch       PetscCall(VecNorm(cg->yk, NORM_2, &ynorm));
56950b47da0SAdam Denchfield       ynorm2 = ynorm * ynorm;
57050b47da0SAdam Denchfield       if (!cg->diag_scaling) {
5719566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old));
5729566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
5739566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, snorm * snorm, &tau_k, cg->alpha));
57450b47da0SAdam Denchfield         beta = tau_k * gkp1_yk / gnorm2_old;
57550b47da0SAdam Denchfield         beta = PetscMax(beta, 0.0);
5769566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
57750b47da0SAdam Denchfield       } else {
5789566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old)); /* Old gtDg */
5799566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5809566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
58150b47da0SAdam Denchfield         beta = gkp1_yk / gnorm2_old;
58250b47da0SAdam Denchfield         beta = PetscMax(beta, 0.0);
5839566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
58450b47da0SAdam Denchfield       }
585c8bcdf1eSAdam Denchfield       break;
586484c7b14SAdam Denchfield 
587484c7b14SAdam Denchfield     case CG_DaiYuan:
588484c7b14SAdam Denchfield       /* Dai, Yu-Hong, and Yaxiang Yuan. "A nonlinear conjugate gradient method with a strong global convergence property."
589484c7b14SAdam Denchfield          SIAM Journal on optimization 10, no. 1 (1999): 177-182. */
59050b47da0SAdam Denchfield       if (!cg->diag_scaling) {
5919566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->stepdirection, tao->gradient, &gd));
5929566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
5939566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, cg->yts, &tau_k, cg->alpha));
59450b47da0SAdam Denchfield         beta = tau_k * gnorm2 / (gd - gd_old);
5959566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
59650b47da0SAdam Denchfield       } else {
5979566063dSJacob Faibussowitsch         PetscCall(MatMult(cg->B, tao->stepdirection, cg->d_work));
5989566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5999566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, tao->gradient, &gtDg));
6009566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->stepdirection, cg->G_old, &gd_old));
6019566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->d_work, cg->g_work, &dk_yk));
60250b47da0SAdam Denchfield         dk_yk = dk_yk - gd_old;
60350b47da0SAdam Denchfield         beta  = gtDg / dk_yk;
6049566063dSJacob Faibussowitsch         PetscCall(VecScale(cg->d_work, beta));
6059566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(tao->stepdirection, -1.0, cg->g_work, cg->d_work));
60650b47da0SAdam Denchfield       }
607c8bcdf1eSAdam Denchfield       break;
608484c7b14SAdam Denchfield 
609c8bcdf1eSAdam Denchfield     case CG_HagerZhang:
610484c7b14SAdam Denchfield       /* Hager, William W., and Hongchao Zhang. "Algorithm 851: CG_DESCENT, a conjugate gradient method with guaranteed descent."
611484c7b14SAdam Denchfield          ACM Transactions on Mathematical Software (TOMS) 32, no. 1 (2006): 113-137. */
6129566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
6139566063dSJacob Faibussowitsch       PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
6149566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
61550b47da0SAdam Denchfield       snorm   = dnorm * step;
61650b47da0SAdam Denchfield       cg->yts = step * dk_yk;
61748a46eb9SPierre Jolivet       if (cg->use_dynamic_restart) PetscCall(TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold));
61850b47da0SAdam Denchfield       if (cg->dynamic_restart) {
6199566063dSJacob Faibussowitsch         PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
620c8bcdf1eSAdam Denchfield       } else {
621c8bcdf1eSAdam Denchfield         if (!cg->diag_scaling) {
6229566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
6239566063dSJacob Faibussowitsch           PetscCall(TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm * snorm, &tau_k, cg->alpha));
624c8bcdf1eSAdam Denchfield           /* Supplying cg->alpha = -1.0 will give the CG_DESCENT 5.3 special case of tau_k = 1.0 */
625c8bcdf1eSAdam Denchfield           tmp  = gd / dk_yk;
626c8bcdf1eSAdam Denchfield           beta = tau_k * (gkp1_yk / dk_yk - ynorm2 * gd / (dk_yk * dk_yk));
627c8bcdf1eSAdam Denchfield           /* Bound beta as in CG_DESCENT 5.3, as implemented, with the third comparison from DK 2013 */
62850b47da0SAdam Denchfield           beta = PetscMax(PetscMax(beta, cg->hz_eta * tau_k * gd_old / (dnorm * dnorm)), cg->dk_eta * tau_k * gd / (dnorm * dnorm));
629c8bcdf1eSAdam Denchfield           /* d <- -t*g + beta*t*d */
6309566063dSJacob Faibussowitsch           PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
631c8bcdf1eSAdam Denchfield         } else {
632c8bcdf1eSAdam Denchfield           /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
633c8bcdf1eSAdam Denchfield           cg->yty = ynorm2;
634c8bcdf1eSAdam Denchfield           cg->sts = snorm * snorm;
63550b47da0SAdam Denchfield           /* Apply the diagonal scaling to all my vectors */
6369566063dSJacob Faibussowitsch           PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
6379566063dSJacob Faibussowitsch           PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
6389566063dSJacob Faibussowitsch           PetscCall(MatSolve(cg->B, tao->stepdirection, cg->d_work));
639c8bcdf1eSAdam Denchfield           /* Construct the constant ytDgkp1 */
6409566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, cg->g_work, &gkp1_yk));
641c8bcdf1eSAdam Denchfield           /* Construct the constant for scaling Dkyk in the update */
642c8bcdf1eSAdam Denchfield           tmp = gd / dk_yk;
6439566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, cg->y_work, &tau_k));
644c8bcdf1eSAdam Denchfield           tau_k = -tau_k * gd / (dk_yk * dk_yk);
645c8bcdf1eSAdam Denchfield           /* beta is the constant which adds the dk contribution */
646484c7b14SAdam Denchfield           beta  = gkp1_yk / dk_yk + cg->hz_theta * tau_k; /* HZ; (1.15) from DK 2013 */
647c8bcdf1eSAdam Denchfield           /* From HZ2013, modified to account for diagonal scaling*/
6489566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->G_old, cg->d_work, &gd_old));
6499566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->stepdirection, cg->g_work, &gd));
65050b47da0SAdam Denchfield           beta = PetscMax(PetscMax(beta, cg->hz_eta * gd_old / (dnorm * dnorm)), cg->dk_eta * gd / (dnorm * dnorm));
651c8bcdf1eSAdam Denchfield           /* Do the update */
6529566063dSJacob Faibussowitsch           PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
65350b47da0SAdam Denchfield         }
65450b47da0SAdam Denchfield       }
655c8bcdf1eSAdam Denchfield       break;
656484c7b14SAdam Denchfield 
657c8bcdf1eSAdam Denchfield     case CG_DaiKou:
658484c7b14SAdam Denchfield       /* Dai, Yu-Hong, and Cai-Xia Kou. "A nonlinear conjugate gradient algorithm with an optimal property and an improved Wolfe line search."
659484c7b14SAdam Denchfield          SIAM Journal on Optimization 23, no. 1 (2013): 296-320. */
6609566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
6619566063dSJacob Faibussowitsch       PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
6629566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
66350b47da0SAdam Denchfield       snorm   = step * dnorm;
66450b47da0SAdam Denchfield       cg->yts = dk_yk * step;
665c8bcdf1eSAdam Denchfield       if (!cg->diag_scaling) {
6669566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
6679566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm * snorm, &tau_k, cg->alpha));
668c8bcdf1eSAdam Denchfield         /* Use cg->alpha = -1.0 to get tau_k = 1.0 as in CG_DESCENT 5.3 */
669c8bcdf1eSAdam Denchfield         tmp  = gd / dk_yk;
67050b47da0SAdam Denchfield         beta = tau_k * (gkp1_yk / dk_yk - ynorm2 * gd / (dk_yk * dk_yk) + gd / (dnorm * dnorm)) - step * gd / dk_yk;
67150b47da0SAdam Denchfield         beta = PetscMax(PetscMax(beta, cg->hz_eta * tau_k * gd_old / (dnorm * dnorm)), cg->dk_eta * tau_k * gd / (dnorm * dnorm));
672c8bcdf1eSAdam Denchfield         /* d <- -t*g + beta*t*d */
6739566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, 0.0, beta, tao->gradient, cg->yk));
674c8bcdf1eSAdam Denchfield       } else {
675c8bcdf1eSAdam Denchfield         /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
676c8bcdf1eSAdam Denchfield         cg->yty = ynorm2;
677c8bcdf1eSAdam Denchfield         cg->sts = snorm * snorm;
6789566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
6799566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
6809566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->stepdirection, cg->d_work));
681c8bcdf1eSAdam Denchfield         /* Construct the constant ytDgkp1 */
6829566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, cg->g_work, &gkp1_yk));
6839566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, cg->y_work, &tau_k));
684c8bcdf1eSAdam Denchfield         tau_k = tau_k * gd / (dk_yk * dk_yk);
685c8bcdf1eSAdam Denchfield         tmp   = gd / dk_yk;
686c8bcdf1eSAdam Denchfield         /* beta is the constant which adds the dk contribution */
687484c7b14SAdam Denchfield         beta  = gkp1_yk / dk_yk - step * tmp - tau_k;
688c8bcdf1eSAdam Denchfield         /* Update this for the last term in beta */
6899566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->y_work, tao->stepdirection, &dk_yk));
690c8bcdf1eSAdam Denchfield         beta += tmp * dk_yk / (dnorm * dnorm); /* projection of y_work onto dk */
6919566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->stepdirection, cg->g_work, &gd));
6929566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->d_work, &gd_old));
69350b47da0SAdam Denchfield         beta = PetscMax(PetscMax(beta, cg->hz_eta * gd_old / (dnorm * dnorm)), cg->dk_eta * gd / (dnorm * dnorm));
694c8bcdf1eSAdam Denchfield         /* Do the update */
6959566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
69650b47da0SAdam Denchfield       }
697c8bcdf1eSAdam Denchfield       break;
698484c7b14SAdam Denchfield 
699c8bcdf1eSAdam Denchfield     case CG_KouDai:
700110fc3b0SBarry Smith       /* Kou, Cai-Xia, and Yu-Hong Dai. "A modified self-scaling memoryless Broyden-Fletcher-Goldfarb-Shanno method for unconstrained optimization."
701484c7b14SAdam Denchfield          Journal of Optimization Theory and Applications 165, no. 1 (2015): 209-224. */
7029566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
7039566063dSJacob Faibussowitsch       PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
7049566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
70550b47da0SAdam Denchfield       snorm   = step * dnorm;
70650b47da0SAdam Denchfield       cg->yts = dk_yk * step;
70748a46eb9SPierre Jolivet       if (cg->use_dynamic_restart) PetscCall(TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold));
70850b47da0SAdam Denchfield       if (cg->dynamic_restart) {
7099566063dSJacob Faibussowitsch         PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
710c8bcdf1eSAdam Denchfield       } else {
711c8bcdf1eSAdam Denchfield         if (!cg->diag_scaling) {
7129566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
7139566063dSJacob Faibussowitsch           PetscCall(TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm * snorm, &tau_k, cg->alpha));
714c8bcdf1eSAdam Denchfield           beta = tau_k * (gkp1_yk / dk_yk - ynorm2 * gd / (dk_yk * dk_yk)) - step * gd / dk_yk;
715c8bcdf1eSAdam Denchfield           if (beta < cg->zeta * tau_k * gd / (dnorm * dnorm)) /* 0.1 is KD's zeta parameter */
716c8bcdf1eSAdam Denchfield           {
717c8bcdf1eSAdam Denchfield             beta  = cg->zeta * tau_k * gd / (dnorm * dnorm);
718c8bcdf1eSAdam Denchfield             gamma = 0.0;
719c8bcdf1eSAdam Denchfield           } else {
720c8bcdf1eSAdam Denchfield             if (gkp1_yk < 0 && cg->neg_xi) gamma = -1.0 * gd / dk_yk;
721484c7b14SAdam Denchfield             /* This seems to be very effective when there's no tau_k scaling.
722484c7b14SAdam Denchfield                This guarantees a large descent step every iteration, going through DK 2015 Lemma 3.1's proof but allowing for negative xi */
723ad540459SPierre Jolivet             else gamma = cg->xi * gd / dk_yk;
724c8bcdf1eSAdam Denchfield           }
725c8bcdf1eSAdam Denchfield           /* d <- -t*g + beta*t*d + t*tmp*yk */
7269566063dSJacob Faibussowitsch           PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, gamma * tau_k, beta, tao->gradient, cg->yk));
727c8bcdf1eSAdam Denchfield         } else {
728c8bcdf1eSAdam Denchfield           /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
729c8bcdf1eSAdam Denchfield           cg->yty = ynorm2;
730c8bcdf1eSAdam Denchfield           cg->sts = snorm * snorm;
7319566063dSJacob Faibussowitsch           PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
7329566063dSJacob Faibussowitsch           PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
733c8bcdf1eSAdam Denchfield           /* Construct the constant ytDgkp1 */
7349566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, cg->g_work, &gkp1D_yk));
735c8bcdf1eSAdam Denchfield           /* Construct the constant for scaling Dkyk in the update */
736c8bcdf1eSAdam Denchfield           gamma = gd / dk_yk;
737c8bcdf1eSAdam Denchfield           /* tau_k = -ytDy/(ytd)^2 * gd */
7389566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, cg->y_work, &tau_k));
739c8bcdf1eSAdam Denchfield           tau_k = tau_k * gd / (dk_yk * dk_yk);
740c8bcdf1eSAdam Denchfield           /* beta is the constant which adds the d_k contribution */
741c8bcdf1eSAdam Denchfield           beta  = gkp1D_yk / dk_yk - step * gamma - tau_k;
742c8bcdf1eSAdam Denchfield           /* Here is the requisite check */
7439566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->stepdirection, cg->g_work, &tmp));
744c8bcdf1eSAdam Denchfield           if (cg->neg_xi) {
745c8bcdf1eSAdam Denchfield             /* modified KD implementation */
746c8bcdf1eSAdam Denchfield             if (gkp1D_yk / dk_yk < 0) gamma = -1.0 * gd / dk_yk;
747ad540459SPierre Jolivet             else gamma = cg->xi * gd / dk_yk;
748c8bcdf1eSAdam Denchfield             if (beta < cg->zeta * tmp / (dnorm * dnorm)) {
749c8bcdf1eSAdam Denchfield               beta  = cg->zeta * tmp / (dnorm * dnorm);
750c8bcdf1eSAdam Denchfield               gamma = 0.0;
751c8bcdf1eSAdam Denchfield             }
752c8bcdf1eSAdam Denchfield           } else { /* original KD 2015 implementation */
753c8bcdf1eSAdam Denchfield             if (beta < cg->zeta * tmp / (dnorm * dnorm)) {
754c8bcdf1eSAdam Denchfield               beta  = cg->zeta * tmp / (dnorm * dnorm);
755c8bcdf1eSAdam Denchfield               gamma = 0.0;
756ad540459SPierre Jolivet             } else gamma = cg->xi * gd / dk_yk;
757c8bcdf1eSAdam Denchfield           }
758c8bcdf1eSAdam Denchfield           /* Do the update in two steps */
7599566063dSJacob Faibussowitsch           PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
7609566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->stepdirection, gamma, cg->y_work));
76150b47da0SAdam Denchfield         }
76250b47da0SAdam Denchfield       }
763c8bcdf1eSAdam Denchfield       break;
764484c7b14SAdam Denchfield 
765484c7b14SAdam Denchfield     case CG_SSML_BFGS:
766484c7b14SAdam Denchfield       /* Perry, J. M. "A class of conjugate gradient algorithms with a two-step variable-metric memory."
767484c7b14SAdam Denchfield          Discussion Papers 269 (1977). */
7689566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
7699566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
770484c7b14SAdam Denchfield       snorm   = step * dnorm;
771484c7b14SAdam Denchfield       cg->yts = dk_yk * step;
772484c7b14SAdam Denchfield       cg->yty = ynorm2;
773484c7b14SAdam Denchfield       cg->sts = snorm * snorm;
774484c7b14SAdam Denchfield       if (!cg->diag_scaling) {
7759566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
7769566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha));
777484c7b14SAdam Denchfield         tmp  = gd / dk_yk;
778484c7b14SAdam Denchfield         beta = tau_k * (gkp1_yk / dk_yk - cg->yty * gd / (dk_yk * dk_yk)) - step * tmp;
779484c7b14SAdam Denchfield         /* d <- -t*g + beta*t*d + t*tmp*yk */
7809566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp * tau_k, beta, tao->gradient, cg->yk));
781484c7b14SAdam Denchfield       } else {
782484c7b14SAdam Denchfield         /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
7839566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
7849566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
785484c7b14SAdam Denchfield         /* compute scalar gamma */
7869566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
7879566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->y_work, cg->yk, &tmp));
788484c7b14SAdam Denchfield         gamma = gd / dk_yk;
789484c7b14SAdam Denchfield         /* Compute scalar beta */
790484c7b14SAdam Denchfield         beta  = (gkp1_yk / dk_yk - gd * tmp / (dk_yk * dk_yk)) - step * gd / dk_yk;
791484c7b14SAdam Denchfield         /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
7929566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work));
793484c7b14SAdam Denchfield       }
794484c7b14SAdam Denchfield       break;
795484c7b14SAdam Denchfield 
796484c7b14SAdam Denchfield     case CG_SSML_DFP:
7979566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
7989566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
799484c7b14SAdam Denchfield       snorm   = step * dnorm;
800484c7b14SAdam Denchfield       cg->yts = dk_yk * step;
801484c7b14SAdam Denchfield       cg->yty = ynorm2;
802484c7b14SAdam Denchfield       cg->sts = snorm * snorm;
803484c7b14SAdam Denchfield       if (!cg->diag_scaling) {
804484c7b14SAdam Denchfield         /* Instead of a regular convex combination, we will solve a quadratic formula. */
8059566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha));
8069566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
807484c7b14SAdam Denchfield         tau_k = cg->dfp_scale * tau_k;
808484c7b14SAdam Denchfield         tmp   = tau_k * gkp1_yk / cg->yty;
809484c7b14SAdam Denchfield         beta  = -step * gd / dk_yk;
810484c7b14SAdam Denchfield         /* d <- -t*g + beta*d + tmp*yk */
8119566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk));
812484c7b14SAdam Denchfield       } else {
813484c7b14SAdam Denchfield         /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless DFP step */
8149566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
8159566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
816484c7b14SAdam Denchfield         /* compute scalar gamma */
8179566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
8189566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->y_work, cg->yk, &tmp));
819484c7b14SAdam Denchfield         gamma = (gkp1_yk / tmp);
820484c7b14SAdam Denchfield         /* Compute scalar beta */
821484c7b14SAdam Denchfield         beta  = -step * gd / dk_yk;
822484c7b14SAdam Denchfield         /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
8239566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work));
824484c7b14SAdam Denchfield       }
825484c7b14SAdam Denchfield       break;
826484c7b14SAdam Denchfield 
827484c7b14SAdam Denchfield     case CG_SSML_BROYDEN:
8289566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
8299566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
830484c7b14SAdam Denchfield       snorm   = step * dnorm;
831484c7b14SAdam Denchfield       cg->yts = step * dk_yk;
832484c7b14SAdam Denchfield       cg->yty = ynorm2;
833484c7b14SAdam Denchfield       cg->sts = snorm * snorm;
834484c7b14SAdam Denchfield       if (!cg->diag_scaling) {
835484c7b14SAdam Denchfield         /* Instead of a regular convex combination, we will solve a quadratic formula. */
8369566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(cg->yty, step * dk_yk, snorm * snorm, &tau_bfgs, cg->bfgs_scale));
8379566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(cg->yty, step * dk_yk, snorm * snorm, &tau_dfp, cg->dfp_scale));
8389566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
839484c7b14SAdam Denchfield         tau_k = cg->theta * tau_bfgs + (1.0 - cg->theta) * tau_dfp;
840484c7b14SAdam Denchfield         /* If bfgs_scale = 1.0, it should reproduce the bfgs tau_bfgs. If bfgs_scale = 0.0,
841484c7b14SAdam Denchfield            it should reproduce the tau_dfp scaling. Same with dfp_scale.   */
842484c7b14SAdam Denchfield         tmp   = cg->theta * tau_bfgs * gd / dk_yk + (1 - cg->theta) * tau_dfp * gkp1_yk / cg->yty;
843484c7b14SAdam Denchfield         beta  = cg->theta * tau_bfgs * (gkp1_yk / dk_yk - cg->yty * gd / (dk_yk * dk_yk)) - step * gd / dk_yk;
844484c7b14SAdam Denchfield         /* d <- -t*g + beta*d + tmp*yk */
8459566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk));
846484c7b14SAdam Denchfield       } else {
847484c7b14SAdam Denchfield         /* We have diagonal scaling enabled */
8489566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
8499566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
850484c7b14SAdam Denchfield         /* compute scalar gamma */
8519566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
8529566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->y_work, cg->yk, &tmp));
853484c7b14SAdam Denchfield         gamma = cg->theta * gd / dk_yk + (1 - cg->theta) * (gkp1_yk / tmp);
854484c7b14SAdam Denchfield         /* Compute scalar beta */
855484c7b14SAdam Denchfield         beta  = cg->theta * (gkp1_yk / dk_yk - gd * tmp / (dk_yk * dk_yk)) - step * gd / dk_yk;
856484c7b14SAdam Denchfield         /* Compute stepdirection dkp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
8579566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work));
858484c7b14SAdam Denchfield       }
859484c7b14SAdam Denchfield       break;
860484c7b14SAdam Denchfield 
8619371c9d4SSatish Balay     default: break;
862c8bcdf1eSAdam Denchfield     }
863c8bcdf1eSAdam Denchfield   }
864c8bcdf1eSAdam Denchfield   PetscFunctionReturn(0);
865c8bcdf1eSAdam Denchfield }
866c8bcdf1eSAdam Denchfield 
8679371c9d4SSatish Balay PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao tao, PetscReal gnorm) {
868c8bcdf1eSAdam Denchfield   TAO_BNCG                    *cg        = (TAO_BNCG *)tao->data;
869c8bcdf1eSAdam Denchfield   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
8708ca2df50S   PetscReal                    step = 1.0, gnorm2, gd, dnorm = 0.0;
871c8bcdf1eSAdam Denchfield   PetscReal                    gnorm2_old, f_old, resnorm, gnorm_old;
872c624ebd3SAlp Dener   PetscBool                    pcgd_fallback = PETSC_FALSE;
873c8bcdf1eSAdam Denchfield 
874c8bcdf1eSAdam Denchfield   PetscFunctionBegin;
875c8bcdf1eSAdam Denchfield   /* We are now going to perform a line search along the direction. */
876c8bcdf1eSAdam Denchfield   /* Store solution and gradient info before it changes */
8779566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->solution, cg->X_old));
8789566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->gradient, cg->G_old));
8799566063dSJacob Faibussowitsch   PetscCall(VecCopy(cg->unprojected_gradient, cg->unprojected_gradient_old));
880c8bcdf1eSAdam Denchfield 
881c8bcdf1eSAdam Denchfield   gnorm_old  = gnorm;
882c8bcdf1eSAdam Denchfield   gnorm2_old = gnorm_old * gnorm_old;
883c8bcdf1eSAdam Denchfield   f_old      = cg->f;
884484c7b14SAdam Denchfield   /* Perform bounded line search. If we are recycling a solution from a previous */
885484c7b14SAdam Denchfield   /* TaoSolve, then we want to immediately skip to calculating a new direction rather than performing a linesearch */
886414d97d3SAlp Dener   if (!(tao->recycle && 0 == tao->niter)) {
887484c7b14SAdam Denchfield     /* Above logic: the below code happens every iteration, except for the first iteration of a recycled TaoSolve */
8889566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
8899566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status));
8909566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
891c8bcdf1eSAdam Denchfield 
892c8bcdf1eSAdam Denchfield     /*  Check linesearch failure */
893c8bcdf1eSAdam Denchfield     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
894c8bcdf1eSAdam Denchfield       ++cg->ls_fails;
895c624ebd3SAlp Dener       if (cg->cg_type == CG_GradientDescent) {
896c8bcdf1eSAdam Denchfield         /* Nothing left to do but fail out of the optimization */
897c8bcdf1eSAdam Denchfield         step        = 0.0;
898c8bcdf1eSAdam Denchfield         tao->reason = TAO_DIVERGED_LS_FAILURE;
899c8bcdf1eSAdam Denchfield       } else {
900484c7b14SAdam Denchfield         /* Restore previous point, perform preconditioned GD and regular GD steps at the last good point */
9019566063dSJacob Faibussowitsch         PetscCall(VecCopy(cg->X_old, tao->solution));
9029566063dSJacob Faibussowitsch         PetscCall(VecCopy(cg->G_old, tao->gradient));
9039566063dSJacob Faibussowitsch         PetscCall(VecCopy(cg->unprojected_gradient_old, cg->unprojected_gradient));
904c8bcdf1eSAdam Denchfield         gnorm  = gnorm_old;
905c8bcdf1eSAdam Denchfield         gnorm2 = gnorm2_old;
906c8bcdf1eSAdam Denchfield         cg->f  = f_old;
907c8bcdf1eSAdam Denchfield 
908484c7b14SAdam Denchfield         /* Fall back on preconditioned CG (so long as you're not already using it) */
909484c7b14SAdam Denchfield         if (cg->cg_type != CG_PCGradientDescent && cg->diag_scaling) {
910e2570530SAlp Dener           pcgd_fallback = PETSC_TRUE;
9119566063dSJacob Faibussowitsch           PetscCall(TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback));
912484c7b14SAdam Denchfield 
9139566063dSJacob Faibussowitsch           PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
9149566063dSJacob Faibussowitsch           PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
915c8bcdf1eSAdam Denchfield 
9169566063dSJacob Faibussowitsch           PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
9179566063dSJacob Faibussowitsch           PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status));
9189566063dSJacob Faibussowitsch           PetscCall(TaoAddLineSearchCounts(tao));
919c8bcdf1eSAdam Denchfield 
920484c7b14SAdam Denchfield           pcgd_fallback = PETSC_FALSE;
921484c7b14SAdam Denchfield           if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
922484c7b14SAdam Denchfield             /* Going to perform a regular gradient descent step. */
923484c7b14SAdam Denchfield             ++cg->ls_fails;
924484c7b14SAdam Denchfield             step = 0.0;
925484c7b14SAdam Denchfield           }
926484c7b14SAdam Denchfield         }
927484c7b14SAdam Denchfield         /* Fall back on the scaled gradient step */
928484c7b14SAdam Denchfield         if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
929484c7b14SAdam Denchfield           ++cg->ls_fails;
9309566063dSJacob Faibussowitsch           PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
9319566063dSJacob Faibussowitsch           PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
9329566063dSJacob Faibussowitsch           PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
9339566063dSJacob Faibussowitsch           PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status));
9349566063dSJacob Faibussowitsch           PetscCall(TaoAddLineSearchCounts(tao));
935484c7b14SAdam Denchfield         }
936484c7b14SAdam Denchfield 
937c8bcdf1eSAdam Denchfield         if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
938c8bcdf1eSAdam Denchfield           /* Nothing left to do but fail out of the optimization */
93950b47da0SAdam Denchfield           ++cg->ls_fails;
940c8bcdf1eSAdam Denchfield           step        = 0.0;
941c8bcdf1eSAdam Denchfield           tao->reason = TAO_DIVERGED_LS_FAILURE;
942484c7b14SAdam Denchfield         } else {
943484c7b14SAdam Denchfield           /* One of the fallbacks worked. Set them both back equal to false. */
944484c7b14SAdam Denchfield           pcgd_fallback = PETSC_FALSE;
945c8bcdf1eSAdam Denchfield         }
946c8bcdf1eSAdam Denchfield       }
947c8bcdf1eSAdam Denchfield     }
948c8bcdf1eSAdam Denchfield     /* Convergence test for line search failure */
949c8bcdf1eSAdam Denchfield     if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
950c8bcdf1eSAdam Denchfield 
951c8bcdf1eSAdam Denchfield     /* Standard convergence test */
9529566063dSJacob Faibussowitsch     PetscCall(VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W));
9539566063dSJacob Faibussowitsch     PetscCall(VecNorm(cg->W, NORM_2, &resnorm));
9543c859ba3SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
9559566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its));
9569566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step));
957dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
958c8bcdf1eSAdam Denchfield     if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
959484c7b14SAdam Denchfield   }
960c8bcdf1eSAdam Denchfield   /* Assert we have an updated step and we need at least one more iteration. */
961c8bcdf1eSAdam Denchfield   /* Calculate the next direction */
962c8bcdf1eSAdam Denchfield   /* Estimate the active set at the new solution */
9639566063dSJacob Faibussowitsch   PetscCall(TaoBNCGEstimateActiveSet(tao, cg->as_type));
964c8bcdf1eSAdam Denchfield   /* Compute the projected gradient and its norm */
9659566063dSJacob Faibussowitsch   PetscCall(VecCopy(cg->unprojected_gradient, tao->gradient));
9669566063dSJacob Faibussowitsch   PetscCall(VecISSet(tao->gradient, cg->active_idx, 0.0));
9679566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
968c8bcdf1eSAdam Denchfield   gnorm2 = gnorm * gnorm;
969c8bcdf1eSAdam Denchfield 
970484c7b14SAdam Denchfield   /* Calculate some quantities used in the StepDirectionUpdate. */
9719566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->stepdirection, NORM_2, &dnorm));
972484c7b14SAdam Denchfield   /* Update the step direction. */
9739566063dSJacob Faibussowitsch   PetscCall(TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback));
974484c7b14SAdam Denchfield   ++tao->niter;
9759566063dSJacob Faibussowitsch   PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
976c8bcdf1eSAdam Denchfield 
977c8bcdf1eSAdam Denchfield   if (cg->cg_type != CG_GradientDescent) {
978c8bcdf1eSAdam Denchfield     /* Figure out which previously active variables became inactive this iteration */
9799566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cg->new_inactives));
98048a46eb9SPierre Jolivet     if (cg->inactive_idx && cg->inactive_old) PetscCall(ISDifference(cg->inactive_idx, cg->inactive_old, &cg->new_inactives));
981c8bcdf1eSAdam Denchfield     /* Selectively reset the CG step those freshly inactive variables */
982c8bcdf1eSAdam Denchfield     if (cg->new_inactives) {
9839566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step));
9849566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad));
9859566063dSJacob Faibussowitsch       PetscCall(VecCopy(cg->inactive_grad, cg->inactive_step));
9869566063dSJacob Faibussowitsch       PetscCall(VecScale(cg->inactive_step, -1.0));
9879566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step));
9889566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad));
989c8bcdf1eSAdam Denchfield     }
990c8bcdf1eSAdam Denchfield     /* Verify that this is a descent direction */
9919566063dSJacob Faibussowitsch     PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
9929566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->stepdirection, NORM_2, &dnorm));
99350b47da0SAdam Denchfield     if (PetscIsInfOrNanReal(gd) || (gd / (dnorm * dnorm) <= -1e10 || gd / (dnorm * dnorm) >= -1e-10)) {
994c8bcdf1eSAdam Denchfield       /* Not a descent direction, so we reset back to projected gradient descent */
9959566063dSJacob Faibussowitsch       PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
9969566063dSJacob Faibussowitsch       PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
997c8bcdf1eSAdam Denchfield       ++cg->descent_error;
998c8bcdf1eSAdam Denchfield     } else {
999c8bcdf1eSAdam Denchfield     }
1000c8bcdf1eSAdam Denchfield   }
1001ac9112b8SAlp Dener   PetscFunctionReturn(0);
1002ac9112b8SAlp Dener }
1003484c7b14SAdam Denchfield 
10049371c9d4SSatish Balay PetscErrorCode TaoBNCGSetH0(Tao tao, Mat H0) {
1005484c7b14SAdam Denchfield   TAO_BNCG *cg = (TAO_BNCG *)tao->data;
1006484c7b14SAdam Denchfield 
1007484c7b14SAdam Denchfield   PetscFunctionBegin;
10089566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)H0));
1009484c7b14SAdam Denchfield   cg->pc = H0;
1010484c7b14SAdam Denchfield   PetscFunctionReturn(0);
1011484c7b14SAdam Denchfield }
1012