xref: /petsc/src/tao/bound/impls/bncg/bncg.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 
2861be54a6SAlp Dener PetscErrorCode TaoBNCGEstimateActiveSet(Tao tao, PetscInt asType)
2961be54a6SAlp Dener {
3061be54a6SAlp Dener   TAO_BNCG                     *cg = (TAO_BNCG *)tao->data;
3161be54a6SAlp Dener 
3261be54a6SAlp Dener   PetscFunctionBegin;
339566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->inactive_old));
3461be54a6SAlp Dener   if (cg->inactive_idx) {
359566063dSJacob Faibussowitsch     PetscCall(ISDuplicate(cg->inactive_idx, &cg->inactive_old));
369566063dSJacob Faibussowitsch     PetscCall(ISCopy(cg->inactive_idx, cg->inactive_old));
3761be54a6SAlp Dener   }
3861be54a6SAlp Dener   switch (asType) {
3961be54a6SAlp Dener   case CG_AS_NONE:
409566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cg->inactive_idx));
419566063dSJacob Faibussowitsch     PetscCall(VecWhichInactive(tao->XL, tao->solution, cg->unprojected_gradient, tao->XU, PETSC_TRUE, &cg->inactive_idx));
429566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cg->active_idx));
439566063dSJacob Faibussowitsch     PetscCall(ISComplementVec(cg->inactive_idx, tao->solution, &cg->active_idx));
4461be54a6SAlp Dener     break;
4561be54a6SAlp Dener 
4661be54a6SAlp Dener   case CG_AS_BERTSEKAS:
4761be54a6SAlp Dener     /* Use gradient descent to estimate the active set */
489566063dSJacob Faibussowitsch     PetscCall(VecCopy(cg->unprojected_gradient, cg->W));
499566063dSJacob Faibussowitsch     PetscCall(VecScale(cg->W, -1.0));
50d0609cedSBarry Smith     PetscCall(TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, cg->unprojected_gradient, cg->W, cg->work, cg->as_step, &cg->as_tol,
51d0609cedSBarry Smith                                       &cg->active_lower, &cg->active_upper, &cg->active_fixed, &cg->active_idx, &cg->inactive_idx));
52c4b75bccSAlp Dener     break;
5361be54a6SAlp Dener   default:
5461be54a6SAlp Dener     break;
5561be54a6SAlp Dener   }
5661be54a6SAlp Dener   PetscFunctionReturn(0);
5761be54a6SAlp Dener }
5861be54a6SAlp Dener 
59a1318120SAlp Dener PetscErrorCode TaoBNCGBoundStep(Tao tao, PetscInt asType, Vec step)
6061be54a6SAlp Dener {
6161be54a6SAlp Dener   TAO_BNCG                     *cg = (TAO_BNCG *)tao->data;
6261be54a6SAlp Dener 
6361be54a6SAlp Dener   PetscFunctionBegin;
64a1318120SAlp Dener   switch (asType) {
6561be54a6SAlp Dener   case CG_AS_NONE:
669566063dSJacob Faibussowitsch     PetscCall(VecISSet(step, cg->active_idx, 0.0));
6761be54a6SAlp Dener     break;
6861be54a6SAlp Dener 
6961be54a6SAlp Dener   case CG_AS_BERTSEKAS:
709566063dSJacob Faibussowitsch     PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, cg->active_lower, cg->active_upper, cg->active_fixed, 1.0, step));
7161be54a6SAlp Dener     break;
7261be54a6SAlp Dener 
7361be54a6SAlp Dener   default:
7461be54a6SAlp Dener     break;
7561be54a6SAlp Dener   }
7661be54a6SAlp Dener   PetscFunctionReturn(0);
7761be54a6SAlp Dener }
7861be54a6SAlp Dener 
79ac9112b8SAlp Dener static PetscErrorCode TaoSolve_BNCG(Tao tao)
80ac9112b8SAlp Dener {
81ac9112b8SAlp Dener   TAO_BNCG                     *cg = (TAO_BNCG*)tao->data;
82484c7b14SAdam Denchfield   PetscReal                    step=1.0,gnorm,gnorm2, resnorm;
83c4b75bccSAlp Dener   PetscInt                     nDiff;
84ac9112b8SAlp Dener 
85ac9112b8SAlp Dener   PetscFunctionBegin;
86ac9112b8SAlp Dener   /*   Project the current point onto the feasible set */
879566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
889566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU));
89ac9112b8SAlp Dener 
90c8bcdf1eSAdam Denchfield   /* Project the initial point onto the feasible region */
919566063dSJacob Faibussowitsch   PetscCall(TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution));
92484c7b14SAdam Denchfield 
93414d97d3SAlp Dener   if (nDiff > 0 || !tao->recycle) {
949566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient));
95484c7b14SAdam Denchfield   }
969566063dSJacob Faibussowitsch   PetscCall(VecNorm(cg->unprojected_gradient,NORM_2,&gnorm));
973c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(cg->f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
98ac9112b8SAlp Dener 
9961be54a6SAlp Dener   /* Estimate the active set and compute the projected gradient */
1009566063dSJacob Faibussowitsch   PetscCall(TaoBNCGEstimateActiveSet(tao, cg->as_type));
10161be54a6SAlp Dener 
102ac9112b8SAlp Dener   /* Project the gradient and calculate the norm */
1039566063dSJacob Faibussowitsch   PetscCall(VecCopy(cg->unprojected_gradient, tao->gradient));
1049566063dSJacob Faibussowitsch   PetscCall(VecISSet(tao->gradient, cg->active_idx, 0.0));
1059566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->gradient,NORM_2,&gnorm));
106ac9112b8SAlp Dener   gnorm2 = gnorm*gnorm;
107ac9112b8SAlp Dener 
108c8bcdf1eSAdam Denchfield   /* Initialize counters */
109e031d6f5SAlp Dener   tao->niter = 0;
11050b47da0SAdam Denchfield   cg->ls_fails = cg->descent_error = 0;
111c8bcdf1eSAdam Denchfield   cg->resets = -1;
112484c7b14SAdam Denchfield   cg->skipped_updates = cg->pure_gd_steps = 0;
113c8bcdf1eSAdam Denchfield   cg->iter_quad = 0;
114c8bcdf1eSAdam Denchfield 
115c8bcdf1eSAdam Denchfield   /* Convergence test at the starting point. */
116ac9112b8SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
117484c7b14SAdam Denchfield 
1189566063dSJacob Faibussowitsch   PetscCall(VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W));
1199566063dSJacob Faibussowitsch   PetscCall(VecNorm(cg->W, NORM_2, &resnorm));
1203c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(resnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
1219566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its));
1229566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step));
123*dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao,convergencetest ,tao->cnvP);
124ac9112b8SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
125484c7b14SAdam Denchfield   /* Calculate initial direction. */
126414d97d3SAlp Dener   if (!tao->recycle) {
127484c7b14SAdam Denchfield     /* We are not recycling a solution/history from a past TaoSolve */
1289566063dSJacob Faibussowitsch     PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
129ac9112b8SAlp Dener   }
130c8bcdf1eSAdam Denchfield   /* Initial gradient descent step. Scaling by 1.0 also does a decent job for some problems. */
131c8bcdf1eSAdam Denchfield   while (1) {
132e1e80dc8SAlp Dener     /* Call general purpose update function */
133e1e80dc8SAlp Dener     if (tao->ops->update) {
134*dbbe0bcdSBarry Smith       PetscUseTypeMethod(tao,update , tao->niter, tao->user_update);
1357494f0b1SStefano Zampini       PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient));
136e1e80dc8SAlp Dener     }
1379566063dSJacob Faibussowitsch     PetscCall(TaoBNCGConductIteration(tao, gnorm));
138c8bcdf1eSAdam Denchfield     if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
139ac9112b8SAlp Dener   }
140ac9112b8SAlp Dener }
141ac9112b8SAlp Dener 
142ac9112b8SAlp Dener static PetscErrorCode TaoSetUp_BNCG(Tao tao)
143ac9112b8SAlp Dener {
144ac9112b8SAlp Dener   TAO_BNCG         *cg = (TAO_BNCG*)tao->data;
145ac9112b8SAlp Dener 
146ac9112b8SAlp Dener   PetscFunctionBegin;
147c4b75bccSAlp Dener   if (!tao->gradient) {
1489566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&tao->gradient));
149c4b75bccSAlp Dener   }
150c4b75bccSAlp Dener   if (!tao->stepdirection) {
1519566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&tao->stepdirection));
152c4b75bccSAlp Dener   }
153c4b75bccSAlp Dener   if (!cg->W) {
1549566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&cg->W));
155c4b75bccSAlp Dener   }
156c4b75bccSAlp Dener   if (!cg->work) {
1579566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&cg->work));
158c4b75bccSAlp Dener   }
159c8bcdf1eSAdam Denchfield   if (!cg->sk) {
1609566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&cg->sk));
161c8bcdf1eSAdam Denchfield   }
162c8bcdf1eSAdam Denchfield   if (!cg->yk) {
1639566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->gradient,&cg->yk));
164c8bcdf1eSAdam Denchfield   }
165c4b75bccSAlp Dener   if (!cg->X_old) {
1669566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&cg->X_old));
167c4b75bccSAlp Dener   }
168c4b75bccSAlp Dener   if (!cg->G_old) {
1699566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->gradient,&cg->G_old));
170c8bcdf1eSAdam Denchfield   }
171c8bcdf1eSAdam Denchfield   if (cg->diag_scaling) {
1729566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&cg->d_work));
1739566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&cg->y_work));
1749566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&cg->g_work));
175c4b75bccSAlp Dener   }
176c4b75bccSAlp Dener   if (!cg->unprojected_gradient) {
1779566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->gradient,&cg->unprojected_gradient));
178c4b75bccSAlp Dener   }
179c4b75bccSAlp Dener   if (!cg->unprojected_gradient_old) {
1809566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->gradient,&cg->unprojected_gradient_old));
181c4b75bccSAlp Dener   }
1829566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(cg->B, cg->sk, cg->yk));
1831baa6e33SBarry Smith   if (cg->pc) PetscCall(MatLMVMSetJ0(cg->B, cg->pc));
184ac9112b8SAlp Dener   PetscFunctionReturn(0);
185ac9112b8SAlp Dener }
186ac9112b8SAlp Dener 
187ac9112b8SAlp Dener static PetscErrorCode TaoDestroy_BNCG(Tao tao)
188ac9112b8SAlp Dener {
189ac9112b8SAlp Dener   TAO_BNCG       *cg = (TAO_BNCG*) tao->data;
190ac9112b8SAlp Dener 
191ac9112b8SAlp Dener   PetscFunctionBegin;
192ac9112b8SAlp Dener   if (tao->setupcalled) {
1939566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->W));
1949566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->work));
1959566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->X_old));
1969566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->G_old));
1979566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->unprojected_gradient));
1989566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->unprojected_gradient_old));
1999566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->g_work));
2009566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->d_work));
2019566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->y_work));
2029566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->sk));
2039566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cg->yk));
204ac9112b8SAlp Dener   }
2059566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->active_lower));
2069566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->active_upper));
2079566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->active_fixed));
2089566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->active_idx));
2099566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->inactive_idx));
2109566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->inactive_old));
2119566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cg->new_inactives));
2129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&cg->B));
213484c7b14SAdam Denchfield   if (cg->pc) {
2149566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&cg->pc));
215484c7b14SAdam Denchfield   }
2169566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
217ac9112b8SAlp Dener   PetscFunctionReturn(0);
218ac9112b8SAlp Dener }
219ac9112b8SAlp Dener 
220*dbbe0bcdSBarry Smith static PetscErrorCode TaoSetFromOptions_BNCG(Tao tao,PetscOptionItems *PetscOptionsObject)
221ac9112b8SAlp Dener {
222ac9112b8SAlp Dener   TAO_BNCG       *cg = (TAO_BNCG*)tao->data;
223ac9112b8SAlp Dener 
224ac9112b8SAlp Dener   PetscFunctionBegin;
225d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Nonlinear Conjugate Gradient method for unconstrained optimization");
2269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_bncg_type","cg formula", "", CG_Table, CGTypes, CG_Table[cg->cg_type], &cg->cg_type,NULL));
2278ebe3e4eSStefano Zampini   if (cg->cg_type != CG_SSML_BFGS) cg->alpha = -1.0; /* Setting defaults for non-BFGS methods. User can change it below. */
228484c7b14SAdam Denchfield   if (CG_GradientDescent == cg->cg_type) {
229484c7b14SAdam Denchfield     cg->cg_type = CG_PCGradientDescent;
230484c7b14SAdam Denchfield     /* Set scaling equal to none or, at best, scalar scaling. */
231484c7b14SAdam Denchfield     cg->unscaled_restart = PETSC_TRUE;
232484c7b14SAdam Denchfield     cg->diag_scaling = PETSC_FALSE;
233484c7b14SAdam Denchfield   }
2349566063dSJacob 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));
2359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_hz_eta","(developer) cutoff tolerance for HZ", "", cg->hz_eta,&cg->hz_eta,NULL));
2369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_eps","(developer) cutoff value for restarts", "", cg->epsilon,&cg->epsilon,NULL));
2379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_dk_eta","(developer) cutoff tolerance for DK", "", cg->dk_eta,&cg->dk_eta,NULL));
2389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_xi","(developer) Parameter in the KD method", "", cg->xi,&cg->xi,NULL));
2399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_theta", "(developer) update parameter for the Broyden method", "", cg->theta, &cg->theta, NULL));
2409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_hz_theta", "(developer) parameter for the HZ (2006) method", "", cg->hz_theta, &cg->hz_theta, NULL));
2419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_alpha","(developer) parameter for the scalar scaling","",cg->alpha,&cg->alpha,NULL));
2429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_bfgs_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->bfgs_scale, &cg->bfgs_scale, NULL));
2439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_dfp_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->dfp_scale, &cg->dfp_scale, NULL));
2449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_bncg_diag_scaling","Enable diagonal Broyden-like preconditioning","",cg->diag_scaling,&cg->diag_scaling,NULL));
2459566063dSJacob 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));
2469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_bncg_unscaled_restart","(developer) use unscaled gradient restarts","",cg->unscaled_restart,&cg->unscaled_restart,NULL));
2479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_zeta", "(developer) Free parameter for the Kou-Dai method", "", cg->zeta, &cg->zeta, NULL));
2489566063dSJacob 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));
2499566063dSJacob 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));
2509566063dSJacob 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));
2519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_bncg_no_scaling","Disable all scaling except in restarts","",cg->no_scaling,&cg->no_scaling,NULL));
252484c7b14SAdam Denchfield   if (cg->no_scaling) {
253484c7b14SAdam Denchfield     cg->diag_scaling = PETSC_FALSE;
254484c7b14SAdam Denchfield     cg->alpha = -1.0;
255484c7b14SAdam Denchfield   }
256b474139fSKarl Rupp   if (cg->alpha == -1.0 && cg->cg_type == CG_KouDai && !cg->diag_scaling) { /* Some more default options that appear to be good. */
257484c7b14SAdam Denchfield     cg->neg_xi = PETSC_TRUE;
258484c7b14SAdam Denchfield   }
2599566063dSJacob 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));
2609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_as_tol", "(developer) initial tolerance used when estimating actively bounded variables","",cg->as_tol,&cg->as_tol,NULL));
2619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_as_step", "(developer) step length used when estimating actively bounded variables","",cg->as_step,&cg->as_step,NULL));
2629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_delta_min", "(developer) minimum scaling factor used for scaled gradient restarts","",cg->delta_min,&cg->delta_min,NULL));
2639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bncg_delta_max", "(developer) maximum scaling factor used for scaled gradient restarts","",cg->delta_max,&cg->delta_max,NULL));
26450b47da0SAdam Denchfield 
265d0609cedSBarry Smith   PetscOptionsHeadEnd();
2669566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(cg->B, ((PetscObject)tao)->prefix));
2679566063dSJacob Faibussowitsch   PetscCall(MatAppendOptionsPrefix(cg->B, "tao_bncg_"));
2689566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(cg->B));
269ac9112b8SAlp Dener   PetscFunctionReturn(0);
270ac9112b8SAlp Dener }
271ac9112b8SAlp Dener 
272ac9112b8SAlp Dener static PetscErrorCode TaoView_BNCG(Tao tao, PetscViewer viewer)
273ac9112b8SAlp Dener {
274ac9112b8SAlp Dener   PetscBool      isascii;
275ac9112b8SAlp Dener   TAO_BNCG       *cg = (TAO_BNCG*)tao->data;
276ac9112b8SAlp Dener 
277ac9112b8SAlp Dener   PetscFunctionBegin;
2789566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
279ac9112b8SAlp Dener   if (isascii) {
2809566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
2819566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cg->cg_type]));
28263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Skipped Stepdirection Updates: %" PetscInt_FMT "\n", cg->skipped_updates));
28363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", cg->resets));
28463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Pure gradient steps: %" PetscInt_FMT "\n", cg->pure_gd_steps));
28563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Not a descent direction: %" PetscInt_FMT "\n", cg->descent_error));
28663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Line search fails: %" PetscInt_FMT "\n", cg->ls_fails));
287484c7b14SAdam Denchfield     if (cg->diag_scaling) {
2889566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
289484c7b14SAdam Denchfield       if (isascii) {
2909566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
2919566063dSJacob Faibussowitsch         PetscCall(MatView(cg->B, viewer));
2929566063dSJacob Faibussowitsch         PetscCall(PetscViewerPopFormat(viewer));
293484c7b14SAdam Denchfield       }
294484c7b14SAdam Denchfield     }
2959566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
296ac9112b8SAlp Dener   }
297ac9112b8SAlp Dener   PetscFunctionReturn(0);
298ac9112b8SAlp Dener }
299ac9112b8SAlp Dener 
300c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal yty, PetscReal yts, PetscReal sts, PetscReal *scale, PetscReal alpha)
301c8bcdf1eSAdam Denchfield {
302c8bcdf1eSAdam Denchfield   PetscReal a, b, c, sig1, sig2;
303c8bcdf1eSAdam Denchfield 
304c8bcdf1eSAdam Denchfield   PetscFunctionBegin;
305c8bcdf1eSAdam Denchfield   *scale = 0.0;
3068ebe3e4eSStefano Zampini   if (1.0 == alpha) *scale = yts/yty;
3078ebe3e4eSStefano Zampini   else if (0.0 == alpha) *scale = sts/yts;
30850b47da0SAdam Denchfield   else if (-1.0 == alpha) *scale = 1.0;
309c8bcdf1eSAdam Denchfield   else {
310c8bcdf1eSAdam Denchfield     a = yty;
311c8bcdf1eSAdam Denchfield     b = yts;
312c8bcdf1eSAdam Denchfield     c = sts;
313c8bcdf1eSAdam Denchfield     a *= alpha;
314c8bcdf1eSAdam Denchfield     b *= -(2.0*alpha - 1.0);
315c8bcdf1eSAdam Denchfield     c *= alpha - 1.0;
316c8bcdf1eSAdam Denchfield     sig1 = (-b + PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a);
317c8bcdf1eSAdam Denchfield     sig2 = (-b - PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a);
318c8bcdf1eSAdam Denchfield     /* accept the positive root as the scalar */
3198ebe3e4eSStefano Zampini     if (sig1 > 0.0) *scale = sig1;
3208ebe3e4eSStefano Zampini     else if (sig2 > 0.0) *scale = sig2;
3218ebe3e4eSStefano Zampini     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Cannot find positive scalar");
322c8bcdf1eSAdam Denchfield   }
323c8bcdf1eSAdam Denchfield   PetscFunctionReturn(0);
324c8bcdf1eSAdam Denchfield }
325c8bcdf1eSAdam Denchfield 
326ac9112b8SAlp Dener /*MC
327ac9112b8SAlp Dener   TAOBNCG - Bound-constrained Nonlinear Conjugate Gradient method.
328ac9112b8SAlp Dener 
329ac9112b8SAlp Dener   Options Database Keys:
33050b47da0SAdam Denchfield + -tao_bncg_recycle - enable recycling the latest calculated gradient vector in subsequent TaoSolve() calls (currently disabled)
331c4b75bccSAlp Dener . -tao_bncg_eta <r> - restart tolerance
33261be54a6SAlp Dener . -tao_bncg_type <taocg_type> - cg formula
333c4b75bccSAlp Dener . -tao_bncg_as_type <none,bertsekas> - active set estimation method
334c4b75bccSAlp Dener . -tao_bncg_as_tol <r> - tolerance used in Bertsekas active-set estimation
335c4b75bccSAlp Dener . -tao_bncg_as_step <r> - trial step length used in Bertsekas active-set estimation
33650b47da0SAdam 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.
33750b47da0SAdam Denchfield . -tao_bncg_theta <r> - convex combination parameter for the Broyden method
33850b47da0SAdam Denchfield . -tao_bncg_hz_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods
33950b47da0SAdam Denchfield . -tao_bncg_dk_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods
34050b47da0SAdam Denchfield . -tao_bncg_xi <r> - Multiplicative constant of the gamma term in the KD method
34150b47da0SAdam Denchfield . -tao_bncg_hz_theta <r> - Multiplicative constant of the theta term for the HZ method
34250b47da0SAdam Denchfield . -tao_bncg_bfgs_scale <r> - Scaling parameter of the bfgs contribution to the scalar Broyden method
34350b47da0SAdam Denchfield . -tao_bncg_dfp_scale <r> - Scaling parameter of the dfp contribution to the scalar Broyden method
34450b47da0SAdam Denchfield . -tao_bncg_diag_scaling <b> - Whether or not to use diagonal initialization/preconditioning for the CG methods. Default True.
34550b47da0SAdam Denchfield . -tao_bncg_dynamic_restart <b> - use dynamic restart strategy in the HZ, DK, KD methods
34650b47da0SAdam Denchfield . -tao_bncg_unscaled_restart <b> - whether or not to scale the gradient when doing gradient descent restarts
34750b47da0SAdam Denchfield . -tao_bncg_zeta <r> - Scaling parameter in the KD method
348484c7b14SAdam Denchfield . -tao_bncg_delta_min <r> - Minimum bound for rescaling during restarted gradient descent steps
349484c7b14SAdam Denchfield . -tao_bncg_delta_max <r> - Maximum bound for rescaling during restarted gradient descent steps
35050b47da0SAdam Denchfield . -tao_bncg_min_quad <i> - Number of quadratic-like steps in a row necessary to do a dynamic restart
35150b47da0SAdam 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
35250b47da0SAdam Denchfield . -tao_bncg_spaced_restart <b> - whether or not to do gradient descent steps every x*n iterations
353484c7b14SAdam Denchfield . -tao_bncg_no_scaling <b> - If true, eliminates all scaling, including defaults.
3543850be85SAlp Dener - -tao_bncg_neg_xi <b> - Whether or not to use negative xi in the KD method under certain conditions
355ac9112b8SAlp Dener 
356ac9112b8SAlp Dener   Notes:
357ac9112b8SAlp Dener     CG formulas are:
3583850be85SAlp Dener + "gd" - Gradient Descent
3593850be85SAlp Dener . "fr" - Fletcher-Reeves
3603850be85SAlp Dener . "pr" - Polak-Ribiere-Polyak
3613850be85SAlp Dener . "prp" - Polak-Ribiere-Plus
3623850be85SAlp Dener . "hs" - Hestenes-Steifel
3633850be85SAlp Dener . "dy" - Dai-Yuan
3643850be85SAlp Dener . "ssml_bfgs" - Self-Scaling Memoryless BFGS
3653850be85SAlp Dener . "ssml_dfp"  - Self-Scaling Memoryless DFP
3663850be85SAlp Dener . "ssml_brdn" - Self-Scaling Memoryless Broyden
3673850be85SAlp Dener . "hz" - Hager-Zhang (CG_DESCENT 5.3)
3683850be85SAlp Dener . "dk" - Dai-Kou (2013)
3693850be85SAlp Dener - "kd" - Kou-Dai (2015)
3709abc5736SPatrick Sanan 
371ac9112b8SAlp Dener   Level: beginner
372ac9112b8SAlp Dener M*/
373ac9112b8SAlp Dener 
374ac9112b8SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNCG(Tao tao)
375ac9112b8SAlp Dener {
376ac9112b8SAlp Dener   TAO_BNCG       *cg;
377ac9112b8SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
378ac9112b8SAlp Dener 
379ac9112b8SAlp Dener   PetscFunctionBegin;
380ac9112b8SAlp Dener   tao->ops->setup = TaoSetUp_BNCG;
381ac9112b8SAlp Dener   tao->ops->solve = TaoSolve_BNCG;
382ac9112b8SAlp Dener   tao->ops->view = TaoView_BNCG;
383ac9112b8SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNCG;
384ac9112b8SAlp Dener   tao->ops->destroy = TaoDestroy_BNCG;
385ac9112b8SAlp Dener 
386ac9112b8SAlp Dener   /* Override default settings (unless already changed) */
387ac9112b8SAlp Dener   if (!tao->max_it_changed) tao->max_it = 2000;
388ac9112b8SAlp Dener   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
389ac9112b8SAlp Dener 
390ac9112b8SAlp Dener   /*  Note: nondefault values should be used for nonlinear conjugate gradient  */
391ac9112b8SAlp Dener   /*  method.  In particular, gtol should be less that 0.5; the value used in  */
392ac9112b8SAlp Dener   /*  Nocedal and Wright is 0.10.  We use the default values for the  */
393ac9112b8SAlp Dener   /*  linesearch because it seems to work better. */
3949566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
3959566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
3969566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
3979566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
398ac9112b8SAlp Dener 
3999566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao,&cg));
400ac9112b8SAlp Dener   tao->data = (void*)cg;
4019566063dSJacob Faibussowitsch   PetscCall(KSPInitializePackage());
4029566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &cg->B));
4039566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)cg->B, (PetscObject)tao, 1));
4049566063dSJacob Faibussowitsch   PetscCall(MatSetType(cg->B, MATLMVMDIAGBROYDEN));
40550b47da0SAdam Denchfield 
406484c7b14SAdam Denchfield   cg->pc = NULL;
407484c7b14SAdam Denchfield 
40850b47da0SAdam Denchfield   cg->dk_eta = 0.5;
40950b47da0SAdam Denchfield   cg->hz_eta = 0.4;
410c8bcdf1eSAdam Denchfield   cg->dynamic_restart = PETSC_FALSE;
411c8bcdf1eSAdam Denchfield   cg->unscaled_restart = PETSC_FALSE;
412484c7b14SAdam Denchfield   cg->no_scaling = PETSC_FALSE;
413484c7b14SAdam Denchfield   cg->delta_min = 1e-7;
414484c7b14SAdam Denchfield   cg->delta_max = 100;
415c8bcdf1eSAdam Denchfield   cg->theta = 1.0;
416c8bcdf1eSAdam Denchfield   cg->hz_theta = 1.0;
417c8bcdf1eSAdam Denchfield   cg->dfp_scale = 1.0;
418c8bcdf1eSAdam Denchfield   cg->bfgs_scale = 1.0;
41950b47da0SAdam Denchfield   cg->zeta = 0.1;
42050b47da0SAdam Denchfield   cg->min_quad = 6;
421c8bcdf1eSAdam Denchfield   cg->min_restart_num = 6; /* As in CG_DESCENT and KD2015*/
422c8bcdf1eSAdam Denchfield   cg->xi = 1.0;
42350b47da0SAdam Denchfield   cg->neg_xi = PETSC_TRUE;
424c8bcdf1eSAdam Denchfield   cg->spaced_restart = PETSC_FALSE;
425c8bcdf1eSAdam Denchfield   cg->tol_quad = 1e-8;
42661be54a6SAlp Dener   cg->as_step = 0.001;
42761be54a6SAlp Dener   cg->as_tol = 0.001;
42850b47da0SAdam Denchfield   cg->eps_23 = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); /* Just a little tighter*/
42961be54a6SAlp Dener   cg->as_type = CG_AS_BERTSEKAS;
430c8bcdf1eSAdam Denchfield   cg->cg_type = CG_SSML_BFGS;
431c8bcdf1eSAdam Denchfield   cg->alpha = 1.0;
432c8bcdf1eSAdam Denchfield   cg->diag_scaling = PETSC_TRUE;
433c8bcdf1eSAdam Denchfield   PetscFunctionReturn(0);
434c8bcdf1eSAdam Denchfield }
435c8bcdf1eSAdam Denchfield 
436c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGResetUpdate(Tao tao, PetscReal gnormsq)
437c8bcdf1eSAdam Denchfield {
438c8bcdf1eSAdam Denchfield    TAO_BNCG          *cg = (TAO_BNCG*)tao->data;
439c8bcdf1eSAdam Denchfield    PetscReal         scaling;
440c8bcdf1eSAdam Denchfield 
441c8bcdf1eSAdam Denchfield    PetscFunctionBegin;
442c8bcdf1eSAdam Denchfield    ++cg->resets;
443c8bcdf1eSAdam Denchfield    scaling = 2.0 * PetscMax(1.0, PetscAbsScalar(cg->f)) / PetscMax(gnormsq, cg->eps_23);
444484c7b14SAdam Denchfield    scaling = PetscMin(cg->delta_max, PetscMax(cg->delta_min, scaling));
445484c7b14SAdam Denchfield    if (cg->unscaled_restart) {
446484c7b14SAdam Denchfield      scaling = 1.0;
447484c7b14SAdam Denchfield      ++cg->pure_gd_steps;
448484c7b14SAdam Denchfield    }
4499566063dSJacob Faibussowitsch    PetscCall(VecAXPBY(tao->stepdirection, -scaling, 0.0, tao->gradient));
450c8bcdf1eSAdam Denchfield    /* Also want to reset our diagonal scaling with each restart */
4511baa6e33SBarry Smith    if (cg->diag_scaling) PetscCall(MatLMVMReset(cg->B, PETSC_FALSE));
452c8bcdf1eSAdam Denchfield    PetscFunctionReturn(0);
453c8bcdf1eSAdam Denchfield  }
454c8bcdf1eSAdam Denchfield 
455c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGCheckDynamicRestart(Tao tao, PetscReal stepsize, PetscReal gd, PetscReal gd_old, PetscBool *dynrestart, PetscReal fold)
456c8bcdf1eSAdam Denchfield {
457c8bcdf1eSAdam Denchfield    TAO_BNCG          *cg = (TAO_BNCG*)tao->data;
458c8bcdf1eSAdam Denchfield    PetscReal         quadinterp;
459c8bcdf1eSAdam Denchfield 
460c8bcdf1eSAdam Denchfield    PetscFunctionBegin;
46150b47da0SAdam Denchfield    if (cg->f < cg->min_quad/10) {
46250b47da0SAdam Denchfield      *dynrestart = PETSC_FALSE;
46350b47da0SAdam Denchfield      PetscFunctionReturn(0);
46450b47da0SAdam Denchfield    } /* just skip this since this strategy doesn't work well for functions near zero */
465484c7b14SAdam Denchfield    quadinterp = 2.0*(cg->f - fold)/(stepsize*(gd + gd_old));
46650b47da0SAdam Denchfield    if (PetscAbs(quadinterp - 1.0) < cg->tol_quad) ++cg->iter_quad;
467c8bcdf1eSAdam Denchfield    else {
468c8bcdf1eSAdam Denchfield      cg->iter_quad = 0;
469c8bcdf1eSAdam Denchfield      *dynrestart = PETSC_FALSE;
470c8bcdf1eSAdam Denchfield    }
471c8bcdf1eSAdam Denchfield    if (cg->iter_quad >= cg->min_quad) {
472c8bcdf1eSAdam Denchfield      cg->iter_quad = 0;
473c8bcdf1eSAdam Denchfield      *dynrestart = PETSC_TRUE;
474c8bcdf1eSAdam Denchfield    }
475c8bcdf1eSAdam Denchfield    PetscFunctionReturn(0);
476c8bcdf1eSAdam Denchfield  }
477c8bcdf1eSAdam Denchfield 
4788ca2df50S PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao tao, PetscReal gnorm2, PetscReal step, PetscReal fold, PetscReal gnorm2_old, PetscReal dnorm, PetscBool pcgd_fallback)
47950b47da0SAdam Denchfield {
480c8bcdf1eSAdam Denchfield   TAO_BNCG          *cg = (TAO_BNCG*)tao->data;
48150b47da0SAdam Denchfield   PetscReal         gamma = 1.0, tau_k, beta;
482484c7b14SAdam Denchfield   PetscReal         tmp = 1.0, ynorm, ynorm2 = 1.0, snorm = 1.0, dk_yk=1.0, gd;
48350b47da0SAdam Denchfield   PetscReal         gkp1_yk, gd_old, tau_bfgs, tau_dfp, gkp1D_yk, gtDg;
484c8bcdf1eSAdam Denchfield   PetscInt          dim;
485484c7b14SAdam Denchfield   PetscBool         cg_restart = PETSC_FALSE;
486c8bcdf1eSAdam Denchfield   PetscFunctionBegin;
487c8bcdf1eSAdam Denchfield 
48850b47da0SAdam Denchfield   /* Local curvature check to see if we need to restart */
489414d97d3SAlp Dener   if (tao->niter >= 1 || tao->recycle) {
4909566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient));
4919566063dSJacob Faibussowitsch     PetscCall(VecNorm(cg->yk, NORM_2, &ynorm));
492c8bcdf1eSAdam Denchfield     ynorm2 = ynorm*ynorm;
4939566063dSJacob Faibussowitsch     PetscCall(VecDot(cg->yk, tao->stepdirection, &dk_yk));
494484c7b14SAdam Denchfield     if (step*dnorm < PETSC_MACHINE_EPSILON || step*dk_yk < PETSC_MACHINE_EPSILON) {
495e2570530SAlp Dener       cg_restart = PETSC_TRUE;
496484c7b14SAdam Denchfield       ++cg->skipped_updates;
497484c7b14SAdam Denchfield     }
49850b47da0SAdam Denchfield     if (cg->spaced_restart) {
4999566063dSJacob Faibussowitsch       PetscCall(VecGetSize(tao->gradient, &dim));
500e2570530SAlp Dener       if (tao->niter % (dim*cg->min_restart_num)) cg_restart = PETSC_TRUE;
50150b47da0SAdam Denchfield     }
50250b47da0SAdam Denchfield   }
50350b47da0SAdam Denchfield   /* If the user wants regular restarts, do it every 6n iterations, where n=dimension */
50450b47da0SAdam Denchfield   if (cg->spaced_restart) {
5059566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->gradient, &dim));
506e2570530SAlp Dener     if (0 == tao->niter % (6*dim)) cg_restart = PETSC_TRUE;
50750b47da0SAdam Denchfield   }
50850b47da0SAdam Denchfield   /* Compute the diagonal scaling vector if applicable */
5091baa6e33SBarry Smith   if (cg->diag_scaling) PetscCall(MatLMVMUpdate(cg->B, tao->solution, tao->gradient));
51050b47da0SAdam Denchfield 
511484c7b14SAdam Denchfield   /* A note on diagonal scaling (to be added to paper):
512484c7b14SAdam Denchfield    For the FR, PR, PRP, and DY methods, the diagonally scaled versions
513484c7b14SAdam Denchfield    must be derived as a preconditioned CG method rather than as
514484c7b14SAdam Denchfield    a Hessian initialization like in the Broyden methods. */
51550b47da0SAdam Denchfield 
516484c7b14SAdam Denchfield   /* In that case, one writes the objective function as
517484c7b14SAdam Denchfield    f(x) \equiv f(Ay). Gradient evaluations yield g(x_k) = A g(Ay_k) = A g(x_k).
518484c7b14SAdam Denchfield    Furthermore, the direction d_k \equiv (x_k - x_{k-1})/step according to
519484c7b14SAdam Denchfield    HZ (2006) becomes A^{-1} d_k, such that d_k^T g_k remains the
520484c7b14SAdam Denchfield    same under preconditioning. Note that A is diagonal, such that A^T = A. */
52150b47da0SAdam Denchfield 
522484c7b14SAdam Denchfield   /* This yields questions like what the dot product d_k^T y_k
523484c7b14SAdam Denchfield    should look like. HZ mistakenly treats that as the same under
524484c7b14SAdam Denchfield    preconditioning, but that is not necessarily true. */
52550b47da0SAdam Denchfield 
526484c7b14SAdam Denchfield   /* Observe y_k \equiv g_k - g_{k-1}, and under the P.C. transformation,
527484c7b14SAdam 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}),
528484c7b14SAdam 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
529484c7b14SAdam Denchfield    NOT the same if our preconditioning matrix is updated between iterations.
530484c7b14SAdam Denchfield    This same issue is found when considering dot products of the form g_{k+1}^T y_k. */
53150b47da0SAdam Denchfield 
53250b47da0SAdam Denchfield   /* Compute CG step direction */
53350b47da0SAdam Denchfield   if (cg_restart) {
5349566063dSJacob Faibussowitsch     PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
535484c7b14SAdam Denchfield   } else if (pcgd_fallback) {
536484c7b14SAdam Denchfield     /* Just like preconditioned CG */
5379566063dSJacob Faibussowitsch     PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5389566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work));
53950b47da0SAdam Denchfield   } else if (ynorm2 > PETSC_MACHINE_EPSILON) {
54050b47da0SAdam Denchfield     switch (cg->cg_type) {
541484c7b14SAdam Denchfield     case CG_PCGradientDescent:
54250b47da0SAdam Denchfield       if (!cg->diag_scaling) {
543484c7b14SAdam Denchfield         if (!cg->no_scaling) {
54450b47da0SAdam Denchfield         cg->sts = step*step*dnorm*dnorm;
5459566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->sts, &tau_k, cg->alpha));
546484c7b14SAdam Denchfield         } else {
547484c7b14SAdam Denchfield           tau_k = 1.0;
548484c7b14SAdam Denchfield           ++cg->pure_gd_steps;
549484c7b14SAdam Denchfield         }
5509566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, 0.0, tao->gradient));
55150b47da0SAdam Denchfield       } else {
5529566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5539566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work));
55450b47da0SAdam Denchfield       }
55550b47da0SAdam Denchfield       break;
556484c7b14SAdam Denchfield 
55750b47da0SAdam Denchfield     case CG_HestenesStiefel:
55850b47da0SAdam Denchfield       /* Classic Hestenes-Stiefel method, modified with scalar and diagonal preconditioning. */
55950b47da0SAdam Denchfield       if (!cg->diag_scaling) {
56050b47da0SAdam Denchfield         cg->sts = step*step*dnorm*dnorm;
5619566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
5629566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->sts, &tau_k, cg->alpha));
56350b47da0SAdam Denchfield         beta = tau_k*gkp1_yk/dk_yk;
5649566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
56550b47da0SAdam Denchfield       } else {
5669566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5679566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, cg->g_work, &gkp1_yk));
56850b47da0SAdam Denchfield         beta = gkp1_yk/dk_yk;
5699566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
57050b47da0SAdam Denchfield       }
571c8bcdf1eSAdam Denchfield       break;
572484c7b14SAdam Denchfield 
573c8bcdf1eSAdam Denchfield     case CG_FletcherReeves:
5749566063dSJacob Faibussowitsch       PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old));
5759566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient));
5769566063dSJacob Faibussowitsch       PetscCall(VecNorm(cg->yk, NORM_2, &ynorm));
57750b47da0SAdam Denchfield       ynorm2 = ynorm*ynorm;
5789566063dSJacob Faibussowitsch       PetscCall(VecDot(cg->yk, tao->stepdirection, &dk_yk));
57950b47da0SAdam Denchfield       if (!cg->diag_scaling) {
5809566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, step*step*dnorm*dnorm, &tau_k, cg->alpha));
58150b47da0SAdam Denchfield         beta = tau_k*gnorm2/gnorm2_old;
5829566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
58350b47da0SAdam Denchfield       } else {
5849566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old)); /* Before it's updated */
5859566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5869566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->gradient, cg->g_work, &tmp));
58750b47da0SAdam Denchfield         beta = tmp/gnorm2_old;
5889566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
58950b47da0SAdam Denchfield       }
590c8bcdf1eSAdam Denchfield       break;
591484c7b14SAdam Denchfield 
59250b47da0SAdam Denchfield     case CG_PolakRibierePolyak:
59350b47da0SAdam Denchfield       snorm = step*dnorm;
59450b47da0SAdam Denchfield       if (!cg->diag_scaling) {
5959566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old));
5969566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
5979566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_k, cg->alpha));
59850b47da0SAdam Denchfield         beta = tau_k*gkp1_yk/gnorm2_old;
5999566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
60050b47da0SAdam Denchfield       } else {
6019566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old));
6029566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
6039566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
60450b47da0SAdam Denchfield         beta = gkp1_yk/gnorm2_old;
6059566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
60650b47da0SAdam Denchfield       }
607c8bcdf1eSAdam Denchfield       break;
608484c7b14SAdam Denchfield 
609c8bcdf1eSAdam Denchfield     case CG_PolakRibierePlus:
6109566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient));
6119566063dSJacob Faibussowitsch       PetscCall(VecNorm(cg->yk, NORM_2, &ynorm));
61250b47da0SAdam Denchfield       ynorm2 = ynorm*ynorm;
61350b47da0SAdam Denchfield       if (!cg->diag_scaling) {
6149566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old));
6159566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
6169566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_k, cg->alpha));
61750b47da0SAdam Denchfield         beta = tau_k*gkp1_yk/gnorm2_old;
61850b47da0SAdam Denchfield         beta = PetscMax(beta, 0.0);
6199566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
62050b47da0SAdam Denchfield       } else {
6219566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old)); /* Old gtDg */
6229566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
6239566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
62450b47da0SAdam Denchfield         beta = gkp1_yk/gnorm2_old;
62550b47da0SAdam Denchfield         beta = PetscMax(beta, 0.0);
6269566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
62750b47da0SAdam Denchfield       }
628c8bcdf1eSAdam Denchfield       break;
629484c7b14SAdam Denchfield 
630484c7b14SAdam Denchfield     case CG_DaiYuan:
631484c7b14SAdam Denchfield       /* Dai, Yu-Hong, and Yaxiang Yuan. "A nonlinear conjugate gradient method with a strong global convergence property."
632484c7b14SAdam Denchfield          SIAM Journal on optimization 10, no. 1 (1999): 177-182. */
63350b47da0SAdam Denchfield       if (!cg->diag_scaling) {
6349566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->stepdirection, tao->gradient, &gd));
6359566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
6369566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->yts, &tau_k, cg->alpha));
63750b47da0SAdam Denchfield         beta = tau_k*gnorm2/(gd - gd_old);
6389566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
63950b47da0SAdam Denchfield       } else {
6409566063dSJacob Faibussowitsch         PetscCall(MatMult(cg->B, tao->stepdirection, cg->d_work));
6419566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
6429566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, tao->gradient, &gtDg));
6439566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->stepdirection, cg->G_old, &gd_old));
6449566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->d_work, cg->g_work, &dk_yk));
64550b47da0SAdam Denchfield         dk_yk = dk_yk - gd_old;
64650b47da0SAdam Denchfield         beta = gtDg/dk_yk;
6479566063dSJacob Faibussowitsch         PetscCall(VecScale(cg->d_work, beta));
6489566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(tao->stepdirection, -1.0, cg->g_work, cg->d_work));
64950b47da0SAdam Denchfield       }
650c8bcdf1eSAdam Denchfield       break;
651484c7b14SAdam Denchfield 
652c8bcdf1eSAdam Denchfield     case CG_HagerZhang:
653484c7b14SAdam Denchfield       /* Hager, William W., and Hongchao Zhang. "Algorithm 851: CG_DESCENT, a conjugate gradient method with guaranteed descent."
654484c7b14SAdam Denchfield          ACM Transactions on Mathematical Software (TOMS) 32, no. 1 (2006): 113-137. */
6559566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
6569566063dSJacob Faibussowitsch       PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
6579566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
65850b47da0SAdam Denchfield       snorm = dnorm*step;
65950b47da0SAdam Denchfield       cg->yts = step*dk_yk;
660c8bcdf1eSAdam Denchfield       if (cg->use_dynamic_restart) {
6619566063dSJacob Faibussowitsch         PetscCall(TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold));
662c8bcdf1eSAdam Denchfield       }
66350b47da0SAdam Denchfield       if (cg->dynamic_restart) {
6649566063dSJacob Faibussowitsch         PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
665c8bcdf1eSAdam Denchfield       } else {
666c8bcdf1eSAdam Denchfield         if (!cg->diag_scaling) {
6679566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
6689566063dSJacob Faibussowitsch           PetscCall(TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha));
669c8bcdf1eSAdam Denchfield           /* Supplying cg->alpha = -1.0 will give the CG_DESCENT 5.3 special case of tau_k = 1.0 */
670c8bcdf1eSAdam Denchfield           tmp = gd/dk_yk;
671c8bcdf1eSAdam Denchfield           beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk));
672c8bcdf1eSAdam Denchfield           /* Bound beta as in CG_DESCENT 5.3, as implemented, with the third comparison from DK 2013 */
67350b47da0SAdam Denchfield           beta = PetscMax(PetscMax(beta, cg->hz_eta*tau_k*gd_old/(dnorm*dnorm)), cg->dk_eta*tau_k*gd/(dnorm*dnorm));
674c8bcdf1eSAdam Denchfield           /* d <- -t*g + beta*t*d */
6759566063dSJacob Faibussowitsch           PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
676c8bcdf1eSAdam Denchfield         } else {
677c8bcdf1eSAdam Denchfield           /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
678c8bcdf1eSAdam Denchfield           cg->yty = ynorm2;
679c8bcdf1eSAdam Denchfield           cg->sts = snorm*snorm;
68050b47da0SAdam Denchfield           /* Apply the diagonal scaling to all my vectors */
6819566063dSJacob Faibussowitsch           PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
6829566063dSJacob Faibussowitsch           PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
6839566063dSJacob Faibussowitsch           PetscCall(MatSolve(cg->B, tao->stepdirection, cg->d_work));
684c8bcdf1eSAdam Denchfield           /* Construct the constant ytDgkp1 */
6859566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, cg->g_work, &gkp1_yk));
686c8bcdf1eSAdam Denchfield           /* Construct the constant for scaling Dkyk in the update */
687c8bcdf1eSAdam Denchfield           tmp = gd/dk_yk;
6889566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, cg->y_work, &tau_k));
689c8bcdf1eSAdam Denchfield           tau_k = -tau_k*gd/(dk_yk*dk_yk);
690c8bcdf1eSAdam Denchfield           /* beta is the constant which adds the dk contribution */
691484c7b14SAdam Denchfield           beta = gkp1_yk/dk_yk + cg->hz_theta*tau_k; /* HZ; (1.15) from DK 2013 */
692c8bcdf1eSAdam Denchfield           /* From HZ2013, modified to account for diagonal scaling*/
6939566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->G_old, cg->d_work, &gd_old));
6949566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->stepdirection, cg->g_work, &gd));
69550b47da0SAdam Denchfield           beta = PetscMax(PetscMax(beta, cg->hz_eta*gd_old/(dnorm*dnorm)), cg->dk_eta*gd/(dnorm*dnorm));
696c8bcdf1eSAdam Denchfield           /* Do the update */
6979566063dSJacob Faibussowitsch           PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
69850b47da0SAdam Denchfield         }
69950b47da0SAdam Denchfield       }
700c8bcdf1eSAdam Denchfield       break;
701484c7b14SAdam Denchfield 
702c8bcdf1eSAdam Denchfield     case CG_DaiKou:
703484c7b14SAdam Denchfield       /* Dai, Yu-Hong, and Cai-Xia Kou. "A nonlinear conjugate gradient algorithm with an optimal property and an improved Wolfe line search."
704484c7b14SAdam Denchfield          SIAM Journal on Optimization 23, no. 1 (2013): 296-320. */
7059566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
7069566063dSJacob Faibussowitsch       PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
7079566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
70850b47da0SAdam Denchfield       snorm = step*dnorm;
70950b47da0SAdam Denchfield       cg->yts = dk_yk*step;
710c8bcdf1eSAdam Denchfield       if (!cg->diag_scaling) {
7119566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
7129566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha));
713c8bcdf1eSAdam Denchfield         /* Use cg->alpha = -1.0 to get tau_k = 1.0 as in CG_DESCENT 5.3 */
714c8bcdf1eSAdam Denchfield         tmp = gd/dk_yk;
71550b47da0SAdam Denchfield         beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk) + gd/(dnorm*dnorm)) - step*gd/dk_yk;
71650b47da0SAdam Denchfield         beta = PetscMax(PetscMax(beta, cg->hz_eta*tau_k*gd_old/(dnorm*dnorm)), cg->dk_eta*tau_k*gd/(dnorm*dnorm));
717c8bcdf1eSAdam Denchfield         /* d <- -t*g + beta*t*d */
7189566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, 0.0, beta, tao->gradient, cg->yk));
719c8bcdf1eSAdam Denchfield       } else {
720c8bcdf1eSAdam Denchfield         /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
721c8bcdf1eSAdam Denchfield         cg->yty = ynorm2;
722c8bcdf1eSAdam Denchfield         cg->sts = snorm*snorm;
7239566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
7249566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
7259566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->stepdirection, cg->d_work));
726c8bcdf1eSAdam Denchfield         /* Construct the constant ytDgkp1 */
7279566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, cg->g_work, &gkp1_yk));
7289566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, cg->y_work, &tau_k));
729c8bcdf1eSAdam Denchfield         tau_k = tau_k*gd/(dk_yk*dk_yk);
730c8bcdf1eSAdam Denchfield         tmp = gd/dk_yk;
731c8bcdf1eSAdam Denchfield         /* beta is the constant which adds the dk contribution */
732484c7b14SAdam Denchfield         beta = gkp1_yk/dk_yk - step*tmp - tau_k;
733c8bcdf1eSAdam Denchfield         /* Update this for the last term in beta */
7349566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->y_work, tao->stepdirection, &dk_yk));
735c8bcdf1eSAdam Denchfield         beta += tmp*dk_yk/(dnorm*dnorm); /* projection of y_work onto dk */
7369566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->stepdirection, cg->g_work, &gd));
7379566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->G_old, cg->d_work, &gd_old));
73850b47da0SAdam Denchfield         beta = PetscMax(PetscMax(beta, cg->hz_eta*gd_old/(dnorm*dnorm)), cg->dk_eta*gd/(dnorm*dnorm));
739c8bcdf1eSAdam Denchfield         /* Do the update */
7409566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
74150b47da0SAdam Denchfield       }
742c8bcdf1eSAdam Denchfield       break;
743484c7b14SAdam Denchfield 
744c8bcdf1eSAdam Denchfield     case CG_KouDai:
745110fc3b0SBarry Smith       /* Kou, Cai-Xia, and Yu-Hong Dai. "A modified self-scaling memoryless Broyden-Fletcher-Goldfarb-Shanno method for unconstrained optimization."
746484c7b14SAdam Denchfield          Journal of Optimization Theory and Applications 165, no. 1 (2015): 209-224. */
7479566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
7489566063dSJacob Faibussowitsch       PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
7499566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
75050b47da0SAdam Denchfield       snorm = step*dnorm;
75150b47da0SAdam Denchfield       cg->yts = dk_yk*step;
752c8bcdf1eSAdam Denchfield       if (cg->use_dynamic_restart) {
7539566063dSJacob Faibussowitsch         PetscCall(TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold));
754c8bcdf1eSAdam Denchfield       }
75550b47da0SAdam Denchfield       if (cg->dynamic_restart) {
7569566063dSJacob Faibussowitsch         PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
757c8bcdf1eSAdam Denchfield       } else {
758c8bcdf1eSAdam Denchfield         if (!cg->diag_scaling) {
7599566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
7609566063dSJacob Faibussowitsch           PetscCall(TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha));
761c8bcdf1eSAdam Denchfield           beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk)) - step*gd/dk_yk;
762c8bcdf1eSAdam Denchfield           if (beta < cg->zeta*tau_k*gd/(dnorm*dnorm)) /* 0.1 is KD's zeta parameter */
763c8bcdf1eSAdam Denchfield           {
764c8bcdf1eSAdam Denchfield             beta = cg->zeta*tau_k*gd/(dnorm*dnorm);
765c8bcdf1eSAdam Denchfield             gamma = 0.0;
766c8bcdf1eSAdam Denchfield           } else {
767c8bcdf1eSAdam Denchfield             if (gkp1_yk < 0 && cg->neg_xi) gamma = -1.0*gd/dk_yk;
768484c7b14SAdam Denchfield             /* This seems to be very effective when there's no tau_k scaling.
769484c7b14SAdam Denchfield                This guarantees a large descent step every iteration, going through DK 2015 Lemma 3.1's proof but allowing for negative xi */
77050b47da0SAdam Denchfield             else {
77150b47da0SAdam Denchfield               gamma = cg->xi*gd/dk_yk;
77250b47da0SAdam Denchfield             }
773c8bcdf1eSAdam Denchfield           }
774c8bcdf1eSAdam Denchfield           /* d <- -t*g + beta*t*d + t*tmp*yk */
7759566063dSJacob Faibussowitsch           PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, gamma*tau_k, beta, tao->gradient, cg->yk));
776c8bcdf1eSAdam Denchfield         } else {
777c8bcdf1eSAdam Denchfield           /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
778c8bcdf1eSAdam Denchfield           cg->yty = ynorm2;
779c8bcdf1eSAdam Denchfield           cg->sts = snorm*snorm;
7809566063dSJacob Faibussowitsch           PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
7819566063dSJacob Faibussowitsch           PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
782c8bcdf1eSAdam Denchfield           /* Construct the constant ytDgkp1 */
7839566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, cg->g_work, &gkp1D_yk));
784c8bcdf1eSAdam Denchfield           /* Construct the constant for scaling Dkyk in the update */
785c8bcdf1eSAdam Denchfield           gamma = gd/dk_yk;
786c8bcdf1eSAdam Denchfield           /* tau_k = -ytDy/(ytd)^2 * gd */
7879566063dSJacob Faibussowitsch           PetscCall(VecDot(cg->yk, cg->y_work, &tau_k));
788c8bcdf1eSAdam Denchfield           tau_k = tau_k*gd/(dk_yk*dk_yk);
789c8bcdf1eSAdam Denchfield           /* beta is the constant which adds the d_k contribution */
790c8bcdf1eSAdam Denchfield           beta = gkp1D_yk/dk_yk - step*gamma - tau_k;
791c8bcdf1eSAdam Denchfield           /* Here is the requisite check */
7929566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->stepdirection, cg->g_work, &tmp));
793c8bcdf1eSAdam Denchfield           if (cg->neg_xi) {
794c8bcdf1eSAdam Denchfield             /* modified KD implementation */
795c8bcdf1eSAdam Denchfield             if (gkp1D_yk/dk_yk < 0) gamma = -1.0*gd/dk_yk;
79650b47da0SAdam Denchfield             else {
79750b47da0SAdam Denchfield               gamma = cg->xi*gd/dk_yk;
79850b47da0SAdam Denchfield             }
799c8bcdf1eSAdam Denchfield             if (beta < cg->zeta*tmp/(dnorm*dnorm)) {
800c8bcdf1eSAdam Denchfield               beta = cg->zeta*tmp/(dnorm*dnorm);
801c8bcdf1eSAdam Denchfield               gamma = 0.0;
802c8bcdf1eSAdam Denchfield             }
803c8bcdf1eSAdam Denchfield           } else { /* original KD 2015 implementation */
804c8bcdf1eSAdam Denchfield             if (beta < cg->zeta*tmp/(dnorm*dnorm)) {
805c8bcdf1eSAdam Denchfield               beta = cg->zeta*tmp/(dnorm*dnorm);
806c8bcdf1eSAdam Denchfield               gamma = 0.0;
807c8bcdf1eSAdam Denchfield             } else {
808c8bcdf1eSAdam Denchfield               gamma = cg->xi*gd/dk_yk;
809c8bcdf1eSAdam Denchfield             }
810c8bcdf1eSAdam Denchfield           }
811c8bcdf1eSAdam Denchfield           /* Do the update in two steps */
8129566063dSJacob Faibussowitsch           PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
8139566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->stepdirection, gamma, cg->y_work));
81450b47da0SAdam Denchfield         }
81550b47da0SAdam Denchfield       }
816c8bcdf1eSAdam Denchfield       break;
817484c7b14SAdam Denchfield 
818484c7b14SAdam Denchfield     case CG_SSML_BFGS:
819484c7b14SAdam Denchfield       /* Perry, J. M. "A class of conjugate gradient algorithms with a two-step variable-metric memory."
820484c7b14SAdam Denchfield          Discussion Papers 269 (1977). */
8219566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
8229566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
823484c7b14SAdam Denchfield       snorm = step*dnorm;
824484c7b14SAdam Denchfield       cg->yts = dk_yk*step;
825484c7b14SAdam Denchfield       cg->yty = ynorm2;
826484c7b14SAdam Denchfield       cg->sts = snorm*snorm;
827484c7b14SAdam Denchfield       if (!cg->diag_scaling) {
8289566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
8299566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha));
830484c7b14SAdam Denchfield         tmp = gd/dk_yk;
831484c7b14SAdam Denchfield         beta = tau_k*(gkp1_yk/dk_yk - cg->yty*gd/(dk_yk*dk_yk)) - step*tmp;
832484c7b14SAdam Denchfield         /* d <- -t*g + beta*t*d + t*tmp*yk */
8339566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp*tau_k, beta, tao->gradient, cg->yk));
834484c7b14SAdam Denchfield       } else {
835484c7b14SAdam Denchfield         /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
8369566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
8379566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
838484c7b14SAdam Denchfield         /* compute scalar gamma */
8399566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
8409566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->y_work, cg->yk, &tmp));
841484c7b14SAdam Denchfield         gamma = gd/dk_yk;
842484c7b14SAdam Denchfield         /* Compute scalar beta */
843484c7b14SAdam Denchfield         beta = (gkp1_yk/dk_yk - gd*tmp/(dk_yk*dk_yk)) - step*gd/dk_yk;
844484c7b14SAdam Denchfield         /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
8459566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work));
846484c7b14SAdam Denchfield       }
847484c7b14SAdam Denchfield       break;
848484c7b14SAdam Denchfield 
849484c7b14SAdam Denchfield     case CG_SSML_DFP:
8509566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
8519566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
852484c7b14SAdam Denchfield       snorm = step*dnorm;
853484c7b14SAdam Denchfield       cg->yts = dk_yk*step;
854484c7b14SAdam Denchfield       cg->yty = ynorm2;
855484c7b14SAdam Denchfield       cg->sts = snorm*snorm;
856484c7b14SAdam Denchfield       if (!cg->diag_scaling) {
857484c7b14SAdam Denchfield         /* Instead of a regular convex combination, we will solve a quadratic formula. */
8589566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha));
8599566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
860484c7b14SAdam Denchfield         tau_k = cg->dfp_scale*tau_k;
861484c7b14SAdam Denchfield         tmp = tau_k*gkp1_yk/cg->yty;
862484c7b14SAdam Denchfield         beta = -step*gd/dk_yk;
863484c7b14SAdam Denchfield         /* d <- -t*g + beta*d + tmp*yk */
8649566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk));
865484c7b14SAdam Denchfield       } else {
866484c7b14SAdam Denchfield         /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless DFP step */
8679566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
8689566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
869484c7b14SAdam Denchfield         /* compute scalar gamma */
8709566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
8719566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->y_work, cg->yk, &tmp));
872484c7b14SAdam Denchfield         gamma = (gkp1_yk/tmp);
873484c7b14SAdam Denchfield         /* Compute scalar beta */
874484c7b14SAdam Denchfield         beta = -step*gd/dk_yk;
875484c7b14SAdam Denchfield         /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
8769566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work));
877484c7b14SAdam Denchfield       }
878484c7b14SAdam Denchfield       break;
879484c7b14SAdam Denchfield 
880484c7b14SAdam Denchfield     case CG_SSML_BROYDEN:
8819566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
8829566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
883484c7b14SAdam Denchfield       snorm = step*dnorm;
884484c7b14SAdam Denchfield       cg->yts = step*dk_yk;
885484c7b14SAdam Denchfield       cg->yty = ynorm2;
886484c7b14SAdam Denchfield       cg->sts = snorm*snorm;
887484c7b14SAdam Denchfield       if (!cg->diag_scaling) {
888484c7b14SAdam Denchfield         /* Instead of a regular convex combination, we will solve a quadratic formula. */
8899566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(cg->yty, step*dk_yk, snorm*snorm, &tau_bfgs, cg->bfgs_scale));
8909566063dSJacob Faibussowitsch         PetscCall(TaoBNCGComputeScalarScaling(cg->yty, step*dk_yk, snorm*snorm, &tau_dfp, cg->dfp_scale));
8919566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
892484c7b14SAdam Denchfield         tau_k = cg->theta*tau_bfgs + (1.0-cg->theta)*tau_dfp;
893484c7b14SAdam Denchfield         /* If bfgs_scale = 1.0, it should reproduce the bfgs tau_bfgs. If bfgs_scale = 0.0,
894484c7b14SAdam Denchfield            it should reproduce the tau_dfp scaling. Same with dfp_scale.   */
895484c7b14SAdam Denchfield         tmp = cg->theta*tau_bfgs*gd/dk_yk + (1-cg->theta)*tau_dfp*gkp1_yk/cg->yty;
896484c7b14SAdam Denchfield         beta = cg->theta*tau_bfgs*(gkp1_yk/dk_yk - cg->yty*gd/(dk_yk*dk_yk)) - step*gd/dk_yk;
897484c7b14SAdam Denchfield         /* d <- -t*g + beta*d + tmp*yk */
8989566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk));
899484c7b14SAdam Denchfield       } else {
900484c7b14SAdam Denchfield         /* We have diagonal scaling enabled */
9019566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
9029566063dSJacob Faibussowitsch         PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
903484c7b14SAdam Denchfield         /* compute scalar gamma */
9049566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
9059566063dSJacob Faibussowitsch         PetscCall(VecDot(cg->y_work, cg->yk, &tmp));
906484c7b14SAdam Denchfield         gamma = cg->theta*gd/dk_yk + (1-cg->theta)*(gkp1_yk/tmp);
907484c7b14SAdam Denchfield         /* Compute scalar beta */
908484c7b14SAdam Denchfield         beta = cg->theta*(gkp1_yk/dk_yk - gd*tmp/(dk_yk*dk_yk)) - step*gd/dk_yk;
909484c7b14SAdam Denchfield         /* Compute stepdirection dkp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
9109566063dSJacob Faibussowitsch         PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work));
911484c7b14SAdam Denchfield       }
912484c7b14SAdam Denchfield       break;
913484c7b14SAdam Denchfield 
914c8bcdf1eSAdam Denchfield     default:
915c8bcdf1eSAdam Denchfield       break;
916484c7b14SAdam Denchfield 
917c8bcdf1eSAdam Denchfield     }
918c8bcdf1eSAdam Denchfield   }
919c8bcdf1eSAdam Denchfield   PetscFunctionReturn(0);
920c8bcdf1eSAdam Denchfield }
921c8bcdf1eSAdam Denchfield 
922c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao tao, PetscReal gnorm)
923c8bcdf1eSAdam Denchfield {
924c8bcdf1eSAdam Denchfield   TAO_BNCG                     *cg = (TAO_BNCG*)tao->data;
925c8bcdf1eSAdam Denchfield   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
9268ca2df50S   PetscReal                    step=1.0,gnorm2,gd,dnorm=0.0;
927c8bcdf1eSAdam Denchfield   PetscReal                    gnorm2_old,f_old,resnorm, gnorm_old;
928c624ebd3SAlp Dener   PetscBool                    pcgd_fallback = PETSC_FALSE;
929c8bcdf1eSAdam Denchfield 
930c8bcdf1eSAdam Denchfield   PetscFunctionBegin;
931c8bcdf1eSAdam Denchfield   /* We are now going to perform a line search along the direction. */
932c8bcdf1eSAdam Denchfield   /* Store solution and gradient info before it changes */
9339566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->solution, cg->X_old));
9349566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->gradient, cg->G_old));
9359566063dSJacob Faibussowitsch   PetscCall(VecCopy(cg->unprojected_gradient, cg->unprojected_gradient_old));
936c8bcdf1eSAdam Denchfield 
937c8bcdf1eSAdam Denchfield   gnorm_old = gnorm;
938c8bcdf1eSAdam Denchfield   gnorm2_old = gnorm_old*gnorm_old;
939c8bcdf1eSAdam Denchfield   f_old = cg->f;
940484c7b14SAdam Denchfield   /* Perform bounded line search. If we are recycling a solution from a previous */
941484c7b14SAdam Denchfield   /* TaoSolve, then we want to immediately skip to calculating a new direction rather than performing a linesearch */
942414d97d3SAlp Dener   if (!(tao->recycle && 0 == tao->niter)) {
943484c7b14SAdam Denchfield     /* Above logic: the below code happens every iteration, except for the first iteration of a recycled TaoSolve */
9449566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
9459566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status));
9469566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
947c8bcdf1eSAdam Denchfield 
948c8bcdf1eSAdam Denchfield     /*  Check linesearch failure */
949c8bcdf1eSAdam Denchfield     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
950c8bcdf1eSAdam Denchfield       ++cg->ls_fails;
951c624ebd3SAlp Dener       if (cg->cg_type == CG_GradientDescent) {
952c8bcdf1eSAdam Denchfield         /* Nothing left to do but fail out of the optimization */
953c8bcdf1eSAdam Denchfield         step = 0.0;
954c8bcdf1eSAdam Denchfield         tao->reason = TAO_DIVERGED_LS_FAILURE;
955c8bcdf1eSAdam Denchfield       } else {
956484c7b14SAdam Denchfield         /* Restore previous point, perform preconditioned GD and regular GD steps at the last good point */
9579566063dSJacob Faibussowitsch         PetscCall(VecCopy(cg->X_old, tao->solution));
9589566063dSJacob Faibussowitsch         PetscCall(VecCopy(cg->G_old, tao->gradient));
9599566063dSJacob Faibussowitsch         PetscCall(VecCopy(cg->unprojected_gradient_old, cg->unprojected_gradient));
960c8bcdf1eSAdam Denchfield         gnorm = gnorm_old;
961c8bcdf1eSAdam Denchfield         gnorm2 = gnorm2_old;
962c8bcdf1eSAdam Denchfield         cg->f = f_old;
963c8bcdf1eSAdam Denchfield 
964484c7b14SAdam Denchfield         /* Fall back on preconditioned CG (so long as you're not already using it) */
965484c7b14SAdam Denchfield         if (cg->cg_type != CG_PCGradientDescent && cg->diag_scaling) {
966e2570530SAlp Dener           pcgd_fallback = PETSC_TRUE;
9679566063dSJacob Faibussowitsch           PetscCall(TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback));
968484c7b14SAdam Denchfield 
9699566063dSJacob Faibussowitsch           PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
9709566063dSJacob Faibussowitsch           PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
971c8bcdf1eSAdam Denchfield 
9729566063dSJacob Faibussowitsch           PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
9739566063dSJacob Faibussowitsch           PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status));
9749566063dSJacob Faibussowitsch           PetscCall(TaoAddLineSearchCounts(tao));
975c8bcdf1eSAdam Denchfield 
976484c7b14SAdam Denchfield           pcgd_fallback = PETSC_FALSE;
977484c7b14SAdam Denchfield           if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
978484c7b14SAdam Denchfield             /* Going to perform a regular gradient descent step. */
979484c7b14SAdam Denchfield             ++cg->ls_fails;
980484c7b14SAdam Denchfield             step = 0.0;
981484c7b14SAdam Denchfield           }
982484c7b14SAdam Denchfield         }
983484c7b14SAdam Denchfield         /* Fall back on the scaled gradient step */
984484c7b14SAdam Denchfield         if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
985484c7b14SAdam Denchfield           ++cg->ls_fails;
9869566063dSJacob Faibussowitsch           PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
9879566063dSJacob Faibussowitsch           PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
9889566063dSJacob Faibussowitsch           PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
9899566063dSJacob Faibussowitsch           PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status));
9909566063dSJacob Faibussowitsch           PetscCall(TaoAddLineSearchCounts(tao));
991484c7b14SAdam Denchfield         }
992484c7b14SAdam Denchfield 
993c8bcdf1eSAdam Denchfield         if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
994c8bcdf1eSAdam Denchfield           /* Nothing left to do but fail out of the optimization */
99550b47da0SAdam Denchfield           ++cg->ls_fails;
996c8bcdf1eSAdam Denchfield           step = 0.0;
997c8bcdf1eSAdam Denchfield           tao->reason = TAO_DIVERGED_LS_FAILURE;
998484c7b14SAdam Denchfield         } else {
999484c7b14SAdam Denchfield           /* One of the fallbacks worked. Set them both back equal to false. */
1000484c7b14SAdam Denchfield           pcgd_fallback = PETSC_FALSE;
1001c8bcdf1eSAdam Denchfield         }
1002c8bcdf1eSAdam Denchfield       }
1003c8bcdf1eSAdam Denchfield     }
1004c8bcdf1eSAdam Denchfield     /* Convergence test for line search failure */
1005c8bcdf1eSAdam Denchfield     if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
1006c8bcdf1eSAdam Denchfield 
1007c8bcdf1eSAdam Denchfield     /* Standard convergence test */
10089566063dSJacob Faibussowitsch     PetscCall(VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W));
10099566063dSJacob Faibussowitsch     PetscCall(VecNorm(cg->W, NORM_2, &resnorm));
10103c859ba3SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(resnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
10119566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its));
10129566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step));
1013*dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao,convergencetest ,tao->cnvP);
1014c8bcdf1eSAdam Denchfield     if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
1015484c7b14SAdam Denchfield   }
1016c8bcdf1eSAdam Denchfield   /* Assert we have an updated step and we need at least one more iteration. */
1017c8bcdf1eSAdam Denchfield   /* Calculate the next direction */
1018c8bcdf1eSAdam Denchfield   /* Estimate the active set at the new solution */
10199566063dSJacob Faibussowitsch   PetscCall(TaoBNCGEstimateActiveSet(tao, cg->as_type));
1020c8bcdf1eSAdam Denchfield   /* Compute the projected gradient and its norm */
10219566063dSJacob Faibussowitsch   PetscCall(VecCopy(cg->unprojected_gradient, tao->gradient));
10229566063dSJacob Faibussowitsch   PetscCall(VecISSet(tao->gradient, cg->active_idx, 0.0));
10239566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->gradient,NORM_2,&gnorm));
1024c8bcdf1eSAdam Denchfield   gnorm2 = gnorm*gnorm;
1025c8bcdf1eSAdam Denchfield 
1026484c7b14SAdam Denchfield   /* Calculate some quantities used in the StepDirectionUpdate. */
10279566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->stepdirection, NORM_2, &dnorm));
1028484c7b14SAdam Denchfield   /* Update the step direction. */
10299566063dSJacob Faibussowitsch   PetscCall(TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback));
1030484c7b14SAdam Denchfield   ++tao->niter;
10319566063dSJacob Faibussowitsch   PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
1032c8bcdf1eSAdam Denchfield 
1033c8bcdf1eSAdam Denchfield   if (cg->cg_type != CG_GradientDescent) {
1034c8bcdf1eSAdam Denchfield     /* Figure out which previously active variables became inactive this iteration */
10359566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cg->new_inactives));
1036c8bcdf1eSAdam Denchfield     if (cg->inactive_idx && cg->inactive_old) {
10379566063dSJacob Faibussowitsch       PetscCall(ISDifference(cg->inactive_idx, cg->inactive_old, &cg->new_inactives));
1038c8bcdf1eSAdam Denchfield     }
1039c8bcdf1eSAdam Denchfield     /* Selectively reset the CG step those freshly inactive variables */
1040c8bcdf1eSAdam Denchfield     if (cg->new_inactives) {
10419566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step));
10429566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad));
10439566063dSJacob Faibussowitsch       PetscCall(VecCopy(cg->inactive_grad, cg->inactive_step));
10449566063dSJacob Faibussowitsch       PetscCall(VecScale(cg->inactive_step, -1.0));
10459566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step));
10469566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad));
1047c8bcdf1eSAdam Denchfield     }
1048c8bcdf1eSAdam Denchfield     /* Verify that this is a descent direction */
10499566063dSJacob Faibussowitsch     PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
10509566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->stepdirection, NORM_2, &dnorm));
105150b47da0SAdam Denchfield     if (PetscIsInfOrNanReal(gd) || (gd/(dnorm*dnorm) <= -1e10 || gd/(dnorm*dnorm) >= -1e-10)) {
1052c8bcdf1eSAdam Denchfield       /* Not a descent direction, so we reset back to projected gradient descent */
10539566063dSJacob Faibussowitsch       PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
10549566063dSJacob Faibussowitsch       PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
1055c8bcdf1eSAdam Denchfield       ++cg->descent_error;
1056c8bcdf1eSAdam Denchfield     } else {
1057c8bcdf1eSAdam Denchfield     }
1058c8bcdf1eSAdam Denchfield   }
1059ac9112b8SAlp Dener   PetscFunctionReturn(0);
1060ac9112b8SAlp Dener }
1061484c7b14SAdam Denchfield 
1062484c7b14SAdam Denchfield PetscErrorCode TaoBNCGSetH0(Tao tao, Mat H0)
1063484c7b14SAdam Denchfield {
1064484c7b14SAdam Denchfield   TAO_BNCG                     *cg = (TAO_BNCG*)tao->data;
1065484c7b14SAdam Denchfield 
1066484c7b14SAdam Denchfield   PetscFunctionBegin;
10679566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)H0));
1068484c7b14SAdam Denchfield   cg->pc = H0;
1069484c7b14SAdam Denchfield   PetscFunctionReturn(0);
1070484c7b14SAdam Denchfield }
1071