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
5d6e07cdcSHong Zhang const char *const TaoBNCGTypes[] = {"gd", "pcgd", "hs", "fr", "prp", "prp_plus", "dy", "hz", "dk", "kd", "ssml_bfgs", "ssml_dfp", "ssml_brdn", "TAOBNCGType", "TAO_BNCG_", NULL};
6ac9112b8SAlp Dener
761be54a6SAlp Dener #define CG_AS_NONE 0
861be54a6SAlp Dener #define CG_AS_BERTSEKAS 1
961be54a6SAlp Dener #define CG_AS_SIZE 2
10ac9112b8SAlp Dener
1161be54a6SAlp Dener static const char *CG_AS_TYPE[64] = {"none", "bertsekas"};
12ac9112b8SAlp Dener
TaoBNCGEstimateActiveSet(Tao tao,PetscInt asType)13d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNCGEstimateActiveSet(Tao tao, PetscInt asType)
14d71ae5a4SJacob Faibussowitsch {
1561be54a6SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data;
1661be54a6SAlp Dener
1761be54a6SAlp Dener PetscFunctionBegin;
189566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->inactive_old));
1961be54a6SAlp Dener if (cg->inactive_idx) {
209566063dSJacob Faibussowitsch PetscCall(ISDuplicate(cg->inactive_idx, &cg->inactive_old));
219566063dSJacob Faibussowitsch PetscCall(ISCopy(cg->inactive_idx, cg->inactive_old));
2261be54a6SAlp Dener }
2361be54a6SAlp Dener switch (asType) {
2461be54a6SAlp Dener case CG_AS_NONE:
259566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->inactive_idx));
269566063dSJacob Faibussowitsch PetscCall(VecWhichInactive(tao->XL, tao->solution, cg->unprojected_gradient, tao->XU, PETSC_TRUE, &cg->inactive_idx));
279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->active_idx));
289566063dSJacob Faibussowitsch PetscCall(ISComplementVec(cg->inactive_idx, tao->solution, &cg->active_idx));
2961be54a6SAlp Dener break;
3061be54a6SAlp Dener case CG_AS_BERTSEKAS:
3161be54a6SAlp Dener /* Use gradient descent to estimate the active set */
329566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->unprojected_gradient, cg->W));
339566063dSJacob Faibussowitsch PetscCall(VecScale(cg->W, -1.0));
349371c9d4SSatish 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));
35c4b75bccSAlp Dener break;
36d71ae5a4SJacob Faibussowitsch default:
37d71ae5a4SJacob Faibussowitsch break;
3861be54a6SAlp Dener }
393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4061be54a6SAlp Dener }
4161be54a6SAlp Dener
TaoBNCGBoundStep(Tao tao,PetscInt asType,Vec step)42d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNCGBoundStep(Tao tao, PetscInt asType, Vec step)
43d71ae5a4SJacob Faibussowitsch {
4461be54a6SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data;
4561be54a6SAlp Dener
4661be54a6SAlp Dener PetscFunctionBegin;
47a1318120SAlp Dener switch (asType) {
48d71ae5a4SJacob Faibussowitsch case CG_AS_NONE:
49976ed0a4SStefano Zampini if (cg->active_idx) PetscCall(VecISSet(step, cg->active_idx, 0.0));
50d71ae5a4SJacob Faibussowitsch break;
51d71ae5a4SJacob Faibussowitsch case CG_AS_BERTSEKAS:
52d71ae5a4SJacob Faibussowitsch PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, cg->active_lower, cg->active_upper, cg->active_fixed, 1.0, step));
53d71ae5a4SJacob Faibussowitsch break;
54d71ae5a4SJacob Faibussowitsch default:
55d71ae5a4SJacob Faibussowitsch break;
5661be54a6SAlp Dener }
573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5861be54a6SAlp Dener }
5961be54a6SAlp Dener
TaoSolve_BNCG(Tao tao)60d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_BNCG(Tao tao)
61d71ae5a4SJacob Faibussowitsch {
62ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data;
63484c7b14SAdam Denchfield PetscReal step = 1.0, gnorm, gnorm2, resnorm;
64c4b75bccSAlp Dener PetscInt nDiff;
65ac9112b8SAlp Dener
66ac9112b8SAlp Dener PetscFunctionBegin;
67ac9112b8SAlp Dener /* Project the current point onto the feasible set */
689566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao));
699566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
70ac9112b8SAlp Dener
71c8bcdf1eSAdam Denchfield /* Project the initial point onto the feasible region */
729566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
73484c7b14SAdam Denchfield
7448a46eb9SPierre Jolivet if (nDiff > 0 || !tao->recycle) PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient));
759566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->unprojected_gradient, NORM_2, &gnorm));
76*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(cg->f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
77ac9112b8SAlp Dener
7861be54a6SAlp Dener /* Estimate the active set and compute the projected gradient */
799566063dSJacob Faibussowitsch PetscCall(TaoBNCGEstimateActiveSet(tao, cg->as_type));
8061be54a6SAlp Dener
81ac9112b8SAlp Dener /* Project the gradient and calculate the norm */
829566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->unprojected_gradient, tao->gradient));
83976ed0a4SStefano Zampini if (cg->active_idx) PetscCall(VecISSet(tao->gradient, cg->active_idx, 0.0));
849566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
85ac9112b8SAlp Dener gnorm2 = gnorm * gnorm;
86ac9112b8SAlp Dener
87c8bcdf1eSAdam Denchfield /* Initialize counters */
88e031d6f5SAlp Dener tao->niter = 0;
8950b47da0SAdam Denchfield cg->ls_fails = cg->descent_error = 0;
90c8bcdf1eSAdam Denchfield cg->resets = -1;
91484c7b14SAdam Denchfield cg->skipped_updates = cg->pure_gd_steps = 0;
92c8bcdf1eSAdam Denchfield cg->iter_quad = 0;
93c8bcdf1eSAdam Denchfield
94c8bcdf1eSAdam Denchfield /* Convergence test at the starting point. */
95ac9112b8SAlp Dener tao->reason = TAO_CONTINUE_ITERATING;
96484c7b14SAdam Denchfield
979566063dSJacob Faibussowitsch PetscCall(VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W));
989566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->W, NORM_2, &resnorm));
99*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
1009566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its));
1019566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step));
102dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
1033ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
104484c7b14SAdam Denchfield /* Calculate initial direction. */
105414d97d3SAlp Dener if (!tao->recycle) {
106484c7b14SAdam Denchfield /* We are not recycling a solution/history from a past TaoSolve */
1079566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
108ac9112b8SAlp Dener }
109c8bcdf1eSAdam Denchfield /* Initial gradient descent step. Scaling by 1.0 also does a decent job for some problems. */
110c8bcdf1eSAdam Denchfield while (1) {
111e1e80dc8SAlp Dener /* Call general purpose update function */
112e1e80dc8SAlp Dener if (tao->ops->update) {
113dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
114270bebe6SStefano Zampini PetscCall(TaoComputeObjective(tao, tao->solution, &cg->f));
115e1e80dc8SAlp Dener }
1169566063dSJacob Faibussowitsch PetscCall(TaoBNCGConductIteration(tao, gnorm));
1173ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
118ac9112b8SAlp Dener }
119ac9112b8SAlp Dener }
120ac9112b8SAlp Dener
TaoSetUp_BNCG(Tao tao)121d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_BNCG(Tao tao)
122d71ae5a4SJacob Faibussowitsch {
123ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data;
124ac9112b8SAlp Dener
125ac9112b8SAlp Dener PetscFunctionBegin;
12648a46eb9SPierre Jolivet if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
12748a46eb9SPierre Jolivet if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
12848a46eb9SPierre Jolivet if (!cg->W) PetscCall(VecDuplicate(tao->solution, &cg->W));
12948a46eb9SPierre Jolivet if (!cg->work) PetscCall(VecDuplicate(tao->solution, &cg->work));
13048a46eb9SPierre Jolivet if (!cg->sk) PetscCall(VecDuplicate(tao->solution, &cg->sk));
13148a46eb9SPierre Jolivet if (!cg->yk) PetscCall(VecDuplicate(tao->gradient, &cg->yk));
13248a46eb9SPierre Jolivet if (!cg->X_old) PetscCall(VecDuplicate(tao->solution, &cg->X_old));
13348a46eb9SPierre Jolivet if (!cg->G_old) PetscCall(VecDuplicate(tao->gradient, &cg->G_old));
134c8bcdf1eSAdam Denchfield if (cg->diag_scaling) {
1359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &cg->d_work));
1369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &cg->y_work));
1379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &cg->g_work));
138c4b75bccSAlp Dener }
13948a46eb9SPierre Jolivet if (!cg->unprojected_gradient) PetscCall(VecDuplicate(tao->gradient, &cg->unprojected_gradient));
14048a46eb9SPierre Jolivet if (!cg->unprojected_gradient_old) PetscCall(VecDuplicate(tao->gradient, &cg->unprojected_gradient_old));
1419566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(cg->B, cg->sk, cg->yk));
1421baa6e33SBarry Smith if (cg->pc) PetscCall(MatLMVMSetJ0(cg->B, cg->pc));
1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
144ac9112b8SAlp Dener }
145ac9112b8SAlp Dener
TaoDestroy_BNCG(Tao tao)146d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_BNCG(Tao tao)
147d71ae5a4SJacob Faibussowitsch {
148ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data;
149ac9112b8SAlp Dener
150ac9112b8SAlp Dener PetscFunctionBegin;
151ac9112b8SAlp Dener if (tao->setupcalled) {
1529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->W));
1539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->work));
1549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->X_old));
1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->G_old));
1569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->unprojected_gradient));
1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->unprojected_gradient_old));
1589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->g_work));
1599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->d_work));
1609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->y_work));
1619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->sk));
1629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->yk));
163ac9112b8SAlp Dener }
1649566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->active_lower));
1659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->active_upper));
1669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->active_fixed));
1679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->active_idx));
1689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->inactive_idx));
1699566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->inactive_old));
1709566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->new_inactives));
1719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&cg->B));
17248a46eb9SPierre Jolivet if (cg->pc) PetscCall(MatDestroy(&cg->pc));
1739566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data));
1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
175ac9112b8SAlp Dener }
176ac9112b8SAlp Dener
TaoSetFromOptions_BNCG(Tao tao,PetscOptionItems PetscOptionsObject)177ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_BNCG(Tao tao, PetscOptionItems PetscOptionsObject)
178d71ae5a4SJacob Faibussowitsch {
179ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data;
180ac9112b8SAlp Dener
181ac9112b8SAlp Dener PetscFunctionBegin;
182d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Nonlinear Conjugate Gradient method for unconstrained optimization");
183d6e07cdcSHong Zhang PetscCall(PetscOptionsEnum("-tao_bncg_type", "CG update formula", "TaoBNCGTypes", TaoBNCGTypes, (PetscEnum)cg->cg_type, (PetscEnum *)&cg->cg_type, NULL));
184d6e07cdcSHong Zhang if (cg->cg_type != TAO_BNCG_SSML_BFGS) cg->alpha = -1.0; /* Setting defaults for non-BFGS methods. User can change it below. */
185d6e07cdcSHong Zhang if (TAO_BNCG_GD == cg->cg_type) {
186d6e07cdcSHong Zhang cg->cg_type = TAO_BNCG_PCGD;
187484c7b14SAdam Denchfield /* Set scaling equal to none or, at best, scalar scaling. */
188484c7b14SAdam Denchfield cg->unscaled_restart = PETSC_TRUE;
189484c7b14SAdam Denchfield cg->diag_scaling = PETSC_FALSE;
190484c7b14SAdam Denchfield }
1919566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_hz_eta", "(developer) cutoff tolerance for HZ", "", cg->hz_eta, &cg->hz_eta, NULL));
1929566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_eps", "(developer) cutoff value for restarts", "", cg->epsilon, &cg->epsilon, NULL));
1939566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_dk_eta", "(developer) cutoff tolerance for DK", "", cg->dk_eta, &cg->dk_eta, NULL));
1949566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_xi", "(developer) Parameter in the KD method", "", cg->xi, &cg->xi, NULL));
1959566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_theta", "(developer) update parameter for the Broyden method", "", cg->theta, &cg->theta, NULL));
1969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_hz_theta", "(developer) parameter for the HZ (2006) method", "", cg->hz_theta, &cg->hz_theta, NULL));
1979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_alpha", "(developer) parameter for the scalar scaling", "", cg->alpha, &cg->alpha, NULL));
1989566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_bfgs_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->bfgs_scale, &cg->bfgs_scale, NULL));
1999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_dfp_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->dfp_scale, &cg->dfp_scale, NULL));
2009566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_bncg_diag_scaling", "Enable diagonal Broyden-like preconditioning", "", cg->diag_scaling, &cg->diag_scaling, NULL));
2019566063dSJacob 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));
2029566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_bncg_unscaled_restart", "(developer) use unscaled gradient restarts", "", cg->unscaled_restart, &cg->unscaled_restart, NULL));
2039566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_zeta", "(developer) Free parameter for the Kou-Dai method", "", cg->zeta, &cg->zeta, NULL));
2049566063dSJacob 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));
2059566063dSJacob 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));
2069566063dSJacob 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));
2079566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_bncg_no_scaling", "Disable all scaling except in restarts", "", cg->no_scaling, &cg->no_scaling, NULL));
208484c7b14SAdam Denchfield if (cg->no_scaling) {
209484c7b14SAdam Denchfield cg->diag_scaling = PETSC_FALSE;
210484c7b14SAdam Denchfield cg->alpha = -1.0;
211484c7b14SAdam Denchfield }
212d6e07cdcSHong Zhang if (cg->alpha == -1.0 && cg->cg_type == TAO_BNCG_KD && !cg->diag_scaling) { /* Some more default options that appear to be good. */
213484c7b14SAdam Denchfield cg->neg_xi = PETSC_TRUE;
214484c7b14SAdam Denchfield }
2159566063dSJacob 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));
216d6e07cdcSHong Zhang PetscCall(PetscOptionsEList("-tao_bncg_as_type", "active set estimation method", "", CG_AS_TYPE, CG_AS_SIZE, CG_AS_TYPE[cg->as_type], &cg->as_type, NULL));
2179566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", cg->as_tol, &cg->as_tol, NULL));
2189566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_as_step", "(developer) step length used when estimating actively bounded variables", "", cg->as_step, &cg->as_step, NULL));
2199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_delta_min", "(developer) minimum scaling factor used for scaled gradient restarts", "", cg->delta_min, &cg->delta_min, NULL));
2209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_delta_max", "(developer) maximum scaling factor used for scaled gradient restarts", "", cg->delta_max, &cg->delta_max, NULL));
22150b47da0SAdam Denchfield
222d0609cedSBarry Smith PetscOptionsHeadEnd();
2239566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(cg->B, ((PetscObject)tao)->prefix));
2249566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(cg->B, "tao_bncg_"));
2259566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(cg->B));
2263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
227ac9112b8SAlp Dener }
228ac9112b8SAlp Dener
TaoView_BNCG(Tao tao,PetscViewer viewer)229d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_BNCG(Tao tao, PetscViewer viewer)
230d71ae5a4SJacob Faibussowitsch {
231ac9112b8SAlp Dener PetscBool isascii;
232ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data;
233ac9112b8SAlp Dener
234ac9112b8SAlp Dener PetscFunctionBegin;
2359566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
236ac9112b8SAlp Dener if (isascii) {
2379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
238d6e07cdcSHong Zhang PetscCall(PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", TaoBNCGTypes[cg->cg_type]));
23963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Skipped Stepdirection Updates: %" PetscInt_FMT "\n", cg->skipped_updates));
24063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", cg->resets));
24163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Pure gradient steps: %" PetscInt_FMT "\n", cg->pure_gd_steps));
24263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Not a descent direction: %" PetscInt_FMT "\n", cg->descent_error));
24363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Line search fails: %" PetscInt_FMT "\n", cg->ls_fails));
244484c7b14SAdam Denchfield if (cg->diag_scaling) {
2459566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
246484c7b14SAdam Denchfield if (isascii) {
2479566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
2489566063dSJacob Faibussowitsch PetscCall(MatView(cg->B, viewer));
2499566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer));
250484c7b14SAdam Denchfield }
251484c7b14SAdam Denchfield }
2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
253ac9112b8SAlp Dener }
2543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
255ac9112b8SAlp Dener }
256ac9112b8SAlp Dener
TaoBNCGComputeScalarScaling(PetscReal yty,PetscReal yts,PetscReal sts,PetscReal * scale,PetscReal alpha)257d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal yty, PetscReal yts, PetscReal sts, PetscReal *scale, PetscReal alpha)
258d71ae5a4SJacob Faibussowitsch {
259c8bcdf1eSAdam Denchfield PetscReal a, b, c, sig1, sig2;
260c8bcdf1eSAdam Denchfield
261c8bcdf1eSAdam Denchfield PetscFunctionBegin;
262c8bcdf1eSAdam Denchfield *scale = 0.0;
2638ebe3e4eSStefano Zampini if (1.0 == alpha) *scale = yts / yty;
2648ebe3e4eSStefano Zampini else if (0.0 == alpha) *scale = sts / yts;
26550b47da0SAdam Denchfield else if (-1.0 == alpha) *scale = 1.0;
266c8bcdf1eSAdam Denchfield else {
267c8bcdf1eSAdam Denchfield a = yty;
268c8bcdf1eSAdam Denchfield b = yts;
269c8bcdf1eSAdam Denchfield c = sts;
270c8bcdf1eSAdam Denchfield a *= alpha;
271c8bcdf1eSAdam Denchfield b *= -(2.0 * alpha - 1.0);
272c8bcdf1eSAdam Denchfield c *= alpha - 1.0;
273c8bcdf1eSAdam Denchfield sig1 = (-b + PetscSqrtReal(b * b - 4.0 * a * c)) / (2.0 * a);
274c8bcdf1eSAdam Denchfield sig2 = (-b - PetscSqrtReal(b * b - 4.0 * a * c)) / (2.0 * a);
275c8bcdf1eSAdam Denchfield /* accept the positive root as the scalar */
2768ebe3e4eSStefano Zampini if (sig1 > 0.0) *scale = sig1;
2778ebe3e4eSStefano Zampini else if (sig2 > 0.0) *scale = sig2;
2788ebe3e4eSStefano Zampini else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Cannot find positive scalar");
279c8bcdf1eSAdam Denchfield }
2803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
281c8bcdf1eSAdam Denchfield }
282c8bcdf1eSAdam Denchfield
283ac9112b8SAlp Dener /*MC
284ac9112b8SAlp Dener TAOBNCG - Bound-constrained Nonlinear Conjugate Gradient method.
285ac9112b8SAlp Dener
286ac9112b8SAlp Dener Options Database Keys:
28750b47da0SAdam Denchfield + -tao_bncg_recycle - enable recycling the latest calculated gradient vector in subsequent TaoSolve() calls (currently disabled)
288c4b75bccSAlp Dener . -tao_bncg_eta <r> - restart tolerance
28961be54a6SAlp Dener . -tao_bncg_type <taocg_type> - cg formula
290c4b75bccSAlp Dener . -tao_bncg_as_type <none,bertsekas> - active set estimation method
291c4b75bccSAlp Dener . -tao_bncg_as_tol <r> - tolerance used in Bertsekas active-set estimation
292c4b75bccSAlp Dener . -tao_bncg_as_step <r> - trial step length used in Bertsekas active-set estimation
29350b47da0SAdam 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.
29450b47da0SAdam Denchfield . -tao_bncg_theta <r> - convex combination parameter for the Broyden method
29550b47da0SAdam Denchfield . -tao_bncg_hz_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods
29650b47da0SAdam Denchfield . -tao_bncg_dk_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods
29750b47da0SAdam Denchfield . -tao_bncg_xi <r> - Multiplicative constant of the gamma term in the KD method
29850b47da0SAdam Denchfield . -tao_bncg_hz_theta <r> - Multiplicative constant of the theta term for the HZ method
29950b47da0SAdam Denchfield . -tao_bncg_bfgs_scale <r> - Scaling parameter of the bfgs contribution to the scalar Broyden method
30050b47da0SAdam Denchfield . -tao_bncg_dfp_scale <r> - Scaling parameter of the dfp contribution to the scalar Broyden method
30150b47da0SAdam Denchfield . -tao_bncg_diag_scaling <b> - Whether or not to use diagonal initialization/preconditioning for the CG methods. Default True.
30250b47da0SAdam Denchfield . -tao_bncg_dynamic_restart <b> - use dynamic restart strategy in the HZ, DK, KD methods
30350b47da0SAdam Denchfield . -tao_bncg_unscaled_restart <b> - whether or not to scale the gradient when doing gradient descent restarts
30450b47da0SAdam Denchfield . -tao_bncg_zeta <r> - Scaling parameter in the KD method
305484c7b14SAdam Denchfield . -tao_bncg_delta_min <r> - Minimum bound for rescaling during restarted gradient descent steps
306484c7b14SAdam Denchfield . -tao_bncg_delta_max <r> - Maximum bound for rescaling during restarted gradient descent steps
30750b47da0SAdam Denchfield . -tao_bncg_min_quad <i> - Number of quadratic-like steps in a row necessary to do a dynamic restart
30850b47da0SAdam 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
30950b47da0SAdam Denchfield . -tao_bncg_spaced_restart <b> - whether or not to do gradient descent steps every x*n iterations
310484c7b14SAdam Denchfield . -tao_bncg_no_scaling <b> - If true, eliminates all scaling, including defaults.
3113850be85SAlp Dener - -tao_bncg_neg_xi <b> - Whether or not to use negative xi in the KD method under certain conditions
312ac9112b8SAlp Dener
313ac9112b8SAlp Dener Notes:
314ac9112b8SAlp Dener CG formulas are:
3153850be85SAlp Dener + "gd" - Gradient Descent
3163850be85SAlp Dener . "fr" - Fletcher-Reeves
3173850be85SAlp Dener . "pr" - Polak-Ribiere-Polyak
3183850be85SAlp Dener . "prp" - Polak-Ribiere-Plus
3193850be85SAlp Dener . "hs" - Hestenes-Steifel
3203850be85SAlp Dener . "dy" - Dai-Yuan
3213850be85SAlp Dener . "ssml_bfgs" - Self-Scaling Memoryless BFGS
3223850be85SAlp Dener . "ssml_dfp" - Self-Scaling Memoryless DFP
3233850be85SAlp Dener . "ssml_brdn" - Self-Scaling Memoryless Broyden
3243850be85SAlp Dener . "hz" - Hager-Zhang (CG_DESCENT 5.3)
3253850be85SAlp Dener . "dk" - Dai-Kou (2013)
3263850be85SAlp Dener - "kd" - Kou-Dai (2015)
3279abc5736SPatrick Sanan
328ac9112b8SAlp Dener Level: beginner
329ac9112b8SAlp Dener M*/
330ac9112b8SAlp Dener
TaoCreate_BNCG(Tao tao)331d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BNCG(Tao tao)
332d71ae5a4SJacob Faibussowitsch {
333ac9112b8SAlp Dener TAO_BNCG *cg;
334ac9112b8SAlp Dener const char *morethuente_type = TAOLINESEARCHMT;
335ac9112b8SAlp Dener
336ac9112b8SAlp Dener PetscFunctionBegin;
337ac9112b8SAlp Dener tao->ops->setup = TaoSetUp_BNCG;
338ac9112b8SAlp Dener tao->ops->solve = TaoSolve_BNCG;
339ac9112b8SAlp Dener tao->ops->view = TaoView_BNCG;
340ac9112b8SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNCG;
341ac9112b8SAlp Dener tao->ops->destroy = TaoDestroy_BNCG;
342ac9112b8SAlp Dener
343ac9112b8SAlp Dener /* Override default settings (unless already changed) */
344606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao));
345606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_it, 2000);
346606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_funcs, 4000);
347ac9112b8SAlp Dener
348ac9112b8SAlp Dener /* Note: nondefault values should be used for nonlinear conjugate gradient */
349606f75f6SBarry Smith /* method. In particular, gtol should be less than 0.5; the value used in */
350ac9112b8SAlp Dener /* Nocedal and Wright is 0.10. We use the default values for the */
351ac9112b8SAlp Dener /* linesearch because it seems to work better. */
3529566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
3539566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
3549566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
3559566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
356ac9112b8SAlp Dener
3574dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&cg));
358ac9112b8SAlp Dener tao->data = (void *)cg;
3599566063dSJacob Faibussowitsch PetscCall(KSPInitializePackage());
3609566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &cg->B));
3619566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)cg->B, (PetscObject)tao, 1));
3629566063dSJacob Faibussowitsch PetscCall(MatSetType(cg->B, MATLMVMDIAGBROYDEN));
36350b47da0SAdam Denchfield
364484c7b14SAdam Denchfield cg->pc = NULL;
365484c7b14SAdam Denchfield
36650b47da0SAdam Denchfield cg->dk_eta = 0.5;
36750b47da0SAdam Denchfield cg->hz_eta = 0.4;
368c8bcdf1eSAdam Denchfield cg->dynamic_restart = PETSC_FALSE;
369c8bcdf1eSAdam Denchfield cg->unscaled_restart = PETSC_FALSE;
370484c7b14SAdam Denchfield cg->no_scaling = PETSC_FALSE;
371484c7b14SAdam Denchfield cg->delta_min = 1e-7;
372484c7b14SAdam Denchfield cg->delta_max = 100;
373c8bcdf1eSAdam Denchfield cg->theta = 1.0;
374c8bcdf1eSAdam Denchfield cg->hz_theta = 1.0;
375c8bcdf1eSAdam Denchfield cg->dfp_scale = 1.0;
376c8bcdf1eSAdam Denchfield cg->bfgs_scale = 1.0;
37750b47da0SAdam Denchfield cg->zeta = 0.1;
37850b47da0SAdam Denchfield cg->min_quad = 6;
379c8bcdf1eSAdam Denchfield cg->min_restart_num = 6; /* As in CG_DESCENT and KD2015*/
380c8bcdf1eSAdam Denchfield cg->xi = 1.0;
38150b47da0SAdam Denchfield cg->neg_xi = PETSC_TRUE;
382c8bcdf1eSAdam Denchfield cg->spaced_restart = PETSC_FALSE;
383c8bcdf1eSAdam Denchfield cg->tol_quad = 1e-8;
38461be54a6SAlp Dener cg->as_step = 0.001;
38561be54a6SAlp Dener cg->as_tol = 0.001;
38650b47da0SAdam Denchfield cg->eps_23 = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0); /* Just a little tighter*/
38761be54a6SAlp Dener cg->as_type = CG_AS_BERTSEKAS;
388d6e07cdcSHong Zhang cg->cg_type = TAO_BNCG_SSML_BFGS;
389c8bcdf1eSAdam Denchfield cg->alpha = 1.0;
390c8bcdf1eSAdam Denchfield cg->diag_scaling = PETSC_TRUE;
3913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
392c8bcdf1eSAdam Denchfield }
393c8bcdf1eSAdam Denchfield
TaoBNCGResetUpdate(Tao tao,PetscReal gnormsq)394d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNCGResetUpdate(Tao tao, PetscReal gnormsq)
395d71ae5a4SJacob Faibussowitsch {
396c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG *)tao->data;
397c8bcdf1eSAdam Denchfield PetscReal scaling;
398c8bcdf1eSAdam Denchfield
399c8bcdf1eSAdam Denchfield PetscFunctionBegin;
400c8bcdf1eSAdam Denchfield ++cg->resets;
401c8bcdf1eSAdam Denchfield scaling = 2.0 * PetscMax(1.0, PetscAbsScalar(cg->f)) / PetscMax(gnormsq, cg->eps_23);
402484c7b14SAdam Denchfield scaling = PetscMin(cg->delta_max, PetscMax(cg->delta_min, scaling));
403484c7b14SAdam Denchfield if (cg->unscaled_restart) {
404484c7b14SAdam Denchfield scaling = 1.0;
405484c7b14SAdam Denchfield ++cg->pure_gd_steps;
406484c7b14SAdam Denchfield }
4079566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -scaling, 0.0, tao->gradient));
408c8bcdf1eSAdam Denchfield /* Also want to reset our diagonal scaling with each restart */
4091baa6e33SBarry Smith if (cg->diag_scaling) PetscCall(MatLMVMReset(cg->B, PETSC_FALSE));
4103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
411c8bcdf1eSAdam Denchfield }
412c8bcdf1eSAdam Denchfield
TaoBNCGCheckDynamicRestart(Tao tao,PetscReal stepsize,PetscReal gd,PetscReal gd_old,PetscBool * dynrestart,PetscReal fold)413d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNCGCheckDynamicRestart(Tao tao, PetscReal stepsize, PetscReal gd, PetscReal gd_old, PetscBool *dynrestart, PetscReal fold)
414d71ae5a4SJacob Faibussowitsch {
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;
4213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
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 }
4333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
434c8bcdf1eSAdam Denchfield }
435c8bcdf1eSAdam Denchfield
TaoBNCGStepDirectionUpdate(Tao tao,PetscReal gnorm2,PetscReal step,PetscReal fold,PetscReal gnorm2_old,PetscReal dnorm,PetscBool pcgd_fallback)436d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao tao, PetscReal gnorm2, PetscReal step, PetscReal fold, PetscReal gnorm2_old, PetscReal dnorm, PetscBool pcgd_fallback)
437d71ae5a4SJacob Faibussowitsch {
438c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG *)tao->data;
43950b47da0SAdam Denchfield PetscReal gamma = 1.0, tau_k, beta;
440484c7b14SAdam Denchfield PetscReal tmp = 1.0, ynorm, ynorm2 = 1.0, snorm = 1.0, dk_yk = 1.0, gd;
44150b47da0SAdam Denchfield PetscReal gkp1_yk, gd_old, tau_bfgs, tau_dfp, gkp1D_yk, gtDg;
442c8bcdf1eSAdam Denchfield PetscInt dim;
443484c7b14SAdam Denchfield PetscBool cg_restart = PETSC_FALSE;
444c8bcdf1eSAdam Denchfield
4454d86920dSPierre Jolivet PetscFunctionBegin;
44650b47da0SAdam Denchfield /* Local curvature check to see if we need to restart */
447414d97d3SAlp Dener if (tao->niter >= 1 || tao->recycle) {
4489566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient));
4499566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->yk, NORM_2, &ynorm));
450c8bcdf1eSAdam Denchfield ynorm2 = ynorm * ynorm;
4519566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->stepdirection, &dk_yk));
452484c7b14SAdam Denchfield if (step * dnorm < PETSC_MACHINE_EPSILON || step * dk_yk < PETSC_MACHINE_EPSILON) {
453e2570530SAlp Dener cg_restart = PETSC_TRUE;
454484c7b14SAdam Denchfield ++cg->skipped_updates;
455484c7b14SAdam Denchfield }
45650b47da0SAdam Denchfield if (cg->spaced_restart) {
4579566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->gradient, &dim));
458e2570530SAlp Dener if (tao->niter % (dim * cg->min_restart_num)) cg_restart = PETSC_TRUE;
45950b47da0SAdam Denchfield }
46050b47da0SAdam Denchfield }
46150b47da0SAdam Denchfield /* If the user wants regular restarts, do it every 6n iterations, where n=dimension */
46250b47da0SAdam Denchfield if (cg->spaced_restart) {
4639566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->gradient, &dim));
464e2570530SAlp Dener if (0 == tao->niter % (6 * dim)) cg_restart = PETSC_TRUE;
46550b47da0SAdam Denchfield }
46650b47da0SAdam Denchfield /* Compute the diagonal scaling vector if applicable */
4671baa6e33SBarry Smith if (cg->diag_scaling) PetscCall(MatLMVMUpdate(cg->B, tao->solution, tao->gradient));
46850b47da0SAdam Denchfield
469484c7b14SAdam Denchfield /* A note on diagonal scaling (to be added to paper):
470484c7b14SAdam Denchfield For the FR, PR, PRP, and DY methods, the diagonally scaled versions
471484c7b14SAdam Denchfield must be derived as a preconditioned CG method rather than as
472484c7b14SAdam Denchfield a Hessian initialization like in the Broyden methods. */
47350b47da0SAdam Denchfield
474484c7b14SAdam Denchfield /* In that case, one writes the objective function as
475484c7b14SAdam Denchfield f(x) \equiv f(Ay). Gradient evaluations yield g(x_k) = A g(Ay_k) = A g(x_k).
476484c7b14SAdam Denchfield Furthermore, the direction d_k \equiv (x_k - x_{k-1})/step according to
477484c7b14SAdam Denchfield HZ (2006) becomes A^{-1} d_k, such that d_k^T g_k remains the
478484c7b14SAdam Denchfield same under preconditioning. Note that A is diagonal, such that A^T = A. */
47950b47da0SAdam Denchfield
480484c7b14SAdam Denchfield /* This yields questions like what the dot product d_k^T y_k
481484c7b14SAdam Denchfield should look like. HZ mistakenly treats that as the same under
482484c7b14SAdam Denchfield preconditioning, but that is not necessarily true. */
48350b47da0SAdam Denchfield
484484c7b14SAdam Denchfield /* Observe y_k \equiv g_k - g_{k-1}, and under the P.C. transformation,
485484c7b14SAdam 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}),
486484c7b14SAdam 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
4877addb90fSBarry Smith NOT the same if our matrix used to construct the preconditioner is updated between iterations.
488484c7b14SAdam Denchfield This same issue is found when considering dot products of the form g_{k+1}^T y_k. */
48950b47da0SAdam Denchfield
49050b47da0SAdam Denchfield /* Compute CG step direction */
49150b47da0SAdam Denchfield if (cg_restart) {
4929566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
493484c7b14SAdam Denchfield } else if (pcgd_fallback) {
494484c7b14SAdam Denchfield /* Just like preconditioned CG */
4959566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
4969566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work));
49750b47da0SAdam Denchfield } else if (ynorm2 > PETSC_MACHINE_EPSILON) {
49850b47da0SAdam Denchfield switch (cg->cg_type) {
499d6e07cdcSHong Zhang case TAO_BNCG_PCGD:
50050b47da0SAdam Denchfield if (!cg->diag_scaling) {
501484c7b14SAdam Denchfield if (!cg->no_scaling) {
50250b47da0SAdam Denchfield cg->sts = step * step * dnorm * dnorm;
5039566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, cg->sts, &tau_k, cg->alpha));
504484c7b14SAdam Denchfield } else {
505484c7b14SAdam Denchfield tau_k = 1.0;
506484c7b14SAdam Denchfield ++cg->pure_gd_steps;
507484c7b14SAdam Denchfield }
5089566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, 0.0, tao->gradient));
50950b47da0SAdam Denchfield } else {
5109566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5119566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work));
51250b47da0SAdam Denchfield }
51350b47da0SAdam Denchfield break;
514484c7b14SAdam Denchfield
515d6e07cdcSHong Zhang case TAO_BNCG_HS:
51650b47da0SAdam Denchfield /* Classic Hestenes-Stiefel method, modified with scalar and diagonal preconditioning. */
51750b47da0SAdam Denchfield if (!cg->diag_scaling) {
51850b47da0SAdam Denchfield cg->sts = step * step * dnorm * dnorm;
5199566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
5209566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, cg->sts, &tau_k, cg->alpha));
52150b47da0SAdam Denchfield beta = tau_k * gkp1_yk / dk_yk;
5229566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
52350b47da0SAdam Denchfield } else {
5249566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5259566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->g_work, &gkp1_yk));
52650b47da0SAdam Denchfield beta = gkp1_yk / dk_yk;
5279566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
52850b47da0SAdam Denchfield }
529c8bcdf1eSAdam Denchfield break;
530484c7b14SAdam Denchfield
531d6e07cdcSHong Zhang case TAO_BNCG_FR:
5329566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old));
5339566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient));
5349566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->yk, NORM_2, &ynorm));
53550b47da0SAdam Denchfield ynorm2 = ynorm * ynorm;
5369566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->stepdirection, &dk_yk));
53750b47da0SAdam Denchfield if (!cg->diag_scaling) {
5389566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, step * step * dnorm * dnorm, &tau_k, cg->alpha));
53950b47da0SAdam Denchfield beta = tau_k * gnorm2 / gnorm2_old;
5409566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
54150b47da0SAdam Denchfield } else {
5429566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old)); /* Before it's updated */
5439566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5449566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, cg->g_work, &tmp));
54550b47da0SAdam Denchfield beta = tmp / gnorm2_old;
5469566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
54750b47da0SAdam Denchfield }
548c8bcdf1eSAdam Denchfield break;
549484c7b14SAdam Denchfield
550d6e07cdcSHong Zhang case TAO_BNCG_PRP:
55150b47da0SAdam Denchfield snorm = step * dnorm;
55250b47da0SAdam Denchfield if (!cg->diag_scaling) {
5539566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old));
5549566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
5559566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, snorm * snorm, &tau_k, cg->alpha));
55650b47da0SAdam Denchfield beta = tau_k * gkp1_yk / gnorm2_old;
5579566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
55850b47da0SAdam Denchfield } else {
5599566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old));
5609566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5619566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
56250b47da0SAdam Denchfield beta = gkp1_yk / gnorm2_old;
5639566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
56450b47da0SAdam Denchfield }
565c8bcdf1eSAdam Denchfield break;
566484c7b14SAdam Denchfield
567d6e07cdcSHong Zhang case TAO_BNCG_PRP_PLUS:
5689566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient));
5699566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->yk, NORM_2, &ynorm));
57050b47da0SAdam Denchfield ynorm2 = ynorm * ynorm;
57150b47da0SAdam Denchfield if (!cg->diag_scaling) {
5729566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old));
5739566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
5749566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, snorm * snorm, &tau_k, cg->alpha));
57550b47da0SAdam Denchfield beta = tau_k * gkp1_yk / gnorm2_old;
57650b47da0SAdam Denchfield beta = PetscMax(beta, 0.0);
5779566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
57850b47da0SAdam Denchfield } else {
5799566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old)); /* Old gtDg */
5809566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
5819566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
58250b47da0SAdam Denchfield beta = gkp1_yk / gnorm2_old;
58350b47da0SAdam Denchfield beta = PetscMax(beta, 0.0);
5849566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
58550b47da0SAdam Denchfield }
586c8bcdf1eSAdam Denchfield break;
587484c7b14SAdam Denchfield
588d6e07cdcSHong Zhang case TAO_BNCG_DY:
589484c7b14SAdam Denchfield /* Dai, Yu-Hong, and Yaxiang Yuan. "A nonlinear conjugate gradient method with a strong global convergence property."
590484c7b14SAdam Denchfield SIAM Journal on optimization 10, no. 1 (1999): 177-182. */
59150b47da0SAdam Denchfield if (!cg->diag_scaling) {
5929566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, tao->gradient, &gd));
5939566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
5949566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, cg->yts, &tau_k, cg->alpha));
59550b47da0SAdam Denchfield beta = tau_k * gnorm2 / (gd - gd_old);
5969566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
59750b47da0SAdam Denchfield } else {
5989566063dSJacob Faibussowitsch PetscCall(MatMult(cg->B, tao->stepdirection, cg->d_work));
5999566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
6009566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, tao->gradient, >Dg));
6019566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, cg->G_old, &gd_old));
6029566063dSJacob Faibussowitsch PetscCall(VecDot(cg->d_work, cg->g_work, &dk_yk));
60350b47da0SAdam Denchfield dk_yk = dk_yk - gd_old;
60450b47da0SAdam Denchfield beta = gtDg / dk_yk;
6059566063dSJacob Faibussowitsch PetscCall(VecScale(cg->d_work, beta));
6069566063dSJacob Faibussowitsch PetscCall(VecWAXPY(tao->stepdirection, -1.0, cg->g_work, cg->d_work));
60750b47da0SAdam Denchfield }
608c8bcdf1eSAdam Denchfield break;
609484c7b14SAdam Denchfield
610d6e07cdcSHong Zhang case TAO_BNCG_HZ:
611484c7b14SAdam Denchfield /* Hager, William W., and Hongchao Zhang. "Algorithm 851: CG_DESCENT, a conjugate gradient method with guaranteed descent."
612484c7b14SAdam Denchfield ACM Transactions on Mathematical Software (TOMS) 32, no. 1 (2006): 113-137. */
6139566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
6149566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
6159566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
61650b47da0SAdam Denchfield snorm = dnorm * step;
61750b47da0SAdam Denchfield cg->yts = step * dk_yk;
61848a46eb9SPierre Jolivet if (cg->use_dynamic_restart) PetscCall(TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold));
61950b47da0SAdam Denchfield if (cg->dynamic_restart) {
6209566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
621c8bcdf1eSAdam Denchfield } else {
622c8bcdf1eSAdam Denchfield if (!cg->diag_scaling) {
6239566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
6249566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm * snorm, &tau_k, cg->alpha));
625c8bcdf1eSAdam Denchfield /* Supplying cg->alpha = -1.0 will give the CG_DESCENT 5.3 special case of tau_k = 1.0 */
626c8bcdf1eSAdam Denchfield tmp = gd / dk_yk;
627c8bcdf1eSAdam Denchfield beta = tau_k * (gkp1_yk / dk_yk - ynorm2 * gd / (dk_yk * dk_yk));
628c8bcdf1eSAdam Denchfield /* Bound beta as in CG_DESCENT 5.3, as implemented, with the third comparison from DK 2013 */
62950b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta * tau_k * gd_old / (dnorm * dnorm)), cg->dk_eta * tau_k * gd / (dnorm * dnorm));
630c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d */
6319566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient));
632c8bcdf1eSAdam Denchfield } else {
633c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
634c8bcdf1eSAdam Denchfield cg->yty = ynorm2;
635c8bcdf1eSAdam Denchfield cg->sts = snorm * snorm;
63650b47da0SAdam Denchfield /* Apply the diagonal scaling to all my vectors */
6379566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
6389566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
6399566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->stepdirection, cg->d_work));
640c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */
6419566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->g_work, &gkp1_yk));
642c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */
643c8bcdf1eSAdam Denchfield tmp = gd / dk_yk;
6449566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->y_work, &tau_k));
645c8bcdf1eSAdam Denchfield tau_k = -tau_k * gd / (dk_yk * dk_yk);
646c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */
647484c7b14SAdam Denchfield beta = gkp1_yk / dk_yk + cg->hz_theta * tau_k; /* HZ; (1.15) from DK 2013 */
648c8bcdf1eSAdam Denchfield /* From HZ2013, modified to account for diagonal scaling*/
6499566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->d_work, &gd_old));
6509566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, cg->g_work, &gd));
65150b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta * gd_old / (dnorm * dnorm)), cg->dk_eta * gd / (dnorm * dnorm));
652c8bcdf1eSAdam Denchfield /* Do the update */
6539566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
65450b47da0SAdam Denchfield }
65550b47da0SAdam Denchfield }
656c8bcdf1eSAdam Denchfield break;
657484c7b14SAdam Denchfield
658d6e07cdcSHong Zhang case TAO_BNCG_DK:
659484c7b14SAdam Denchfield /* Dai, Yu-Hong, and Cai-Xia Kou. "A nonlinear conjugate gradient algorithm with an optimal property and an improved Wolfe line search."
660484c7b14SAdam Denchfield SIAM Journal on Optimization 23, no. 1 (2013): 296-320. */
6619566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
6629566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
6639566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
66450b47da0SAdam Denchfield snorm = step * dnorm;
66550b47da0SAdam Denchfield cg->yts = dk_yk * step;
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 /* Use cg->alpha = -1.0 to get tau_k = 1.0 as in CG_DESCENT 5.3 */
670c8bcdf1eSAdam Denchfield tmp = gd / dk_yk;
67150b47da0SAdam Denchfield beta = tau_k * (gkp1_yk / dk_yk - ynorm2 * gd / (dk_yk * dk_yk) + gd / (dnorm * dnorm)) - step * gd / dk_yk;
67250b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta * tau_k * gd_old / (dnorm * dnorm)), cg->dk_eta * tau_k * gd / (dnorm * dnorm));
673c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d */
6749566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, 0.0, beta, tao->gradient, cg->yk));
675c8bcdf1eSAdam Denchfield } else {
676c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
677c8bcdf1eSAdam Denchfield cg->yty = ynorm2;
678c8bcdf1eSAdam Denchfield cg->sts = snorm * snorm;
6799566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
6809566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
6819566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->stepdirection, cg->d_work));
682c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */
6839566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->g_work, &gkp1_yk));
6849566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->y_work, &tau_k));
685c8bcdf1eSAdam Denchfield tau_k = tau_k * gd / (dk_yk * dk_yk);
686c8bcdf1eSAdam Denchfield tmp = gd / dk_yk;
687c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */
688484c7b14SAdam Denchfield beta = gkp1_yk / dk_yk - step * tmp - tau_k;
689c8bcdf1eSAdam Denchfield /* Update this for the last term in beta */
6909566063dSJacob Faibussowitsch PetscCall(VecDot(cg->y_work, tao->stepdirection, &dk_yk));
691c8bcdf1eSAdam Denchfield beta += tmp * dk_yk / (dnorm * dnorm); /* projection of y_work onto dk */
6929566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, cg->g_work, &gd));
6939566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->d_work, &gd_old));
69450b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta * gd_old / (dnorm * dnorm)), cg->dk_eta * gd / (dnorm * dnorm));
695c8bcdf1eSAdam Denchfield /* Do the update */
6969566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
69750b47da0SAdam Denchfield }
698c8bcdf1eSAdam Denchfield break;
699484c7b14SAdam Denchfield
700d6e07cdcSHong Zhang case TAO_BNCG_KD:
701110fc3b0SBarry Smith /* Kou, Cai-Xia, and Yu-Hong Dai. "A modified self-scaling memoryless Broyden-Fletcher-Goldfarb-Shanno method for unconstrained optimization."
702484c7b14SAdam Denchfield Journal of Optimization Theory and Applications 165, no. 1 (2015): 209-224. */
7039566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
7049566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old));
7059566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
70650b47da0SAdam Denchfield snorm = step * dnorm;
70750b47da0SAdam Denchfield cg->yts = dk_yk * step;
70848a46eb9SPierre Jolivet if (cg->use_dynamic_restart) PetscCall(TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold));
70950b47da0SAdam Denchfield if (cg->dynamic_restart) {
7109566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
711c8bcdf1eSAdam Denchfield } else {
712c8bcdf1eSAdam Denchfield if (!cg->diag_scaling) {
7139566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
7149566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm * snorm, &tau_k, cg->alpha));
715c8bcdf1eSAdam Denchfield beta = tau_k * (gkp1_yk / dk_yk - ynorm2 * gd / (dk_yk * dk_yk)) - step * gd / dk_yk;
716c8bcdf1eSAdam Denchfield if (beta < cg->zeta * tau_k * gd / (dnorm * dnorm)) /* 0.1 is KD's zeta parameter */
717c8bcdf1eSAdam Denchfield {
718c8bcdf1eSAdam Denchfield beta = cg->zeta * tau_k * gd / (dnorm * dnorm);
719c8bcdf1eSAdam Denchfield gamma = 0.0;
720c8bcdf1eSAdam Denchfield } else {
721c8bcdf1eSAdam Denchfield if (gkp1_yk < 0 && cg->neg_xi) gamma = -1.0 * gd / dk_yk;
722484c7b14SAdam Denchfield /* This seems to be very effective when there's no tau_k scaling.
723484c7b14SAdam Denchfield This guarantees a large descent step every iteration, going through DK 2015 Lemma 3.1's proof but allowing for negative xi */
724ad540459SPierre Jolivet else gamma = cg->xi * gd / dk_yk;
725c8bcdf1eSAdam Denchfield }
726c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d + t*tmp*yk */
7279566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, gamma * tau_k, beta, tao->gradient, cg->yk));
728c8bcdf1eSAdam Denchfield } else {
729c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
730c8bcdf1eSAdam Denchfield cg->yty = ynorm2;
731c8bcdf1eSAdam Denchfield cg->sts = snorm * snorm;
7329566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
7339566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
734c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */
7359566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->g_work, &gkp1D_yk));
736c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */
737c8bcdf1eSAdam Denchfield gamma = gd / dk_yk;
738c8bcdf1eSAdam Denchfield /* tau_k = -ytDy/(ytd)^2 * gd */
7399566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->y_work, &tau_k));
740c8bcdf1eSAdam Denchfield tau_k = tau_k * gd / (dk_yk * dk_yk);
741c8bcdf1eSAdam Denchfield /* beta is the constant which adds the d_k contribution */
742c8bcdf1eSAdam Denchfield beta = gkp1D_yk / dk_yk - step * gamma - tau_k;
743c8bcdf1eSAdam Denchfield /* Here is the requisite check */
7449566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, cg->g_work, &tmp));
745c8bcdf1eSAdam Denchfield if (cg->neg_xi) {
746c8bcdf1eSAdam Denchfield /* modified KD implementation */
747c8bcdf1eSAdam Denchfield if (gkp1D_yk / dk_yk < 0) gamma = -1.0 * gd / dk_yk;
748ad540459SPierre Jolivet else gamma = cg->xi * gd / dk_yk;
749c8bcdf1eSAdam Denchfield if (beta < cg->zeta * tmp / (dnorm * dnorm)) {
750c8bcdf1eSAdam Denchfield beta = cg->zeta * tmp / (dnorm * dnorm);
751c8bcdf1eSAdam Denchfield gamma = 0.0;
752c8bcdf1eSAdam Denchfield }
753c8bcdf1eSAdam Denchfield } else { /* original KD 2015 implementation */
754c8bcdf1eSAdam Denchfield if (beta < cg->zeta * tmp / (dnorm * dnorm)) {
755c8bcdf1eSAdam Denchfield beta = cg->zeta * tmp / (dnorm * dnorm);
756c8bcdf1eSAdam Denchfield gamma = 0.0;
757ad540459SPierre Jolivet } else gamma = cg->xi * gd / dk_yk;
758c8bcdf1eSAdam Denchfield }
759c8bcdf1eSAdam Denchfield /* Do the update in two steps */
7609566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work));
7619566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->stepdirection, gamma, cg->y_work));
76250b47da0SAdam Denchfield }
76350b47da0SAdam Denchfield }
764c8bcdf1eSAdam Denchfield break;
765484c7b14SAdam Denchfield
766d6e07cdcSHong Zhang case TAO_BNCG_SSML_BFGS:
767484c7b14SAdam Denchfield /* Perry, J. M. "A class of conjugate gradient algorithms with a two-step variable-metric memory."
768484c7b14SAdam Denchfield Discussion Papers 269 (1977). */
7699566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
7709566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
771484c7b14SAdam Denchfield snorm = step * dnorm;
772484c7b14SAdam Denchfield cg->yts = dk_yk * step;
773484c7b14SAdam Denchfield cg->yty = ynorm2;
774484c7b14SAdam Denchfield cg->sts = snorm * snorm;
775484c7b14SAdam Denchfield if (!cg->diag_scaling) {
7769566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
7779566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha));
778484c7b14SAdam Denchfield tmp = gd / dk_yk;
779484c7b14SAdam Denchfield beta = tau_k * (gkp1_yk / dk_yk - cg->yty * gd / (dk_yk * dk_yk)) - step * tmp;
780484c7b14SAdam Denchfield /* d <- -t*g + beta*t*d + t*tmp*yk */
7819566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp * tau_k, beta, tao->gradient, cg->yk));
782484c7b14SAdam Denchfield } else {
783484c7b14SAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
7849566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
7859566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
786484c7b14SAdam Denchfield /* compute scalar gamma */
7879566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
7889566063dSJacob Faibussowitsch PetscCall(VecDot(cg->y_work, cg->yk, &tmp));
789484c7b14SAdam Denchfield gamma = gd / dk_yk;
790484c7b14SAdam Denchfield /* Compute scalar beta */
791484c7b14SAdam Denchfield beta = (gkp1_yk / dk_yk - gd * tmp / (dk_yk * dk_yk)) - step * gd / dk_yk;
792484c7b14SAdam Denchfield /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
7939566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work));
794484c7b14SAdam Denchfield }
795484c7b14SAdam Denchfield break;
796484c7b14SAdam Denchfield
797d6e07cdcSHong Zhang case TAO_BNCG_SSML_DFP:
7989566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
7999566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
800484c7b14SAdam Denchfield snorm = step * dnorm;
801484c7b14SAdam Denchfield cg->yts = dk_yk * step;
802484c7b14SAdam Denchfield cg->yty = ynorm2;
803484c7b14SAdam Denchfield cg->sts = snorm * snorm;
804484c7b14SAdam Denchfield if (!cg->diag_scaling) {
805484c7b14SAdam Denchfield /* Instead of a regular convex combination, we will solve a quadratic formula. */
8069566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha));
8079566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
808484c7b14SAdam Denchfield tau_k = cg->dfp_scale * tau_k;
809484c7b14SAdam Denchfield tmp = tau_k * gkp1_yk / cg->yty;
810484c7b14SAdam Denchfield beta = -step * gd / dk_yk;
811484c7b14SAdam Denchfield /* d <- -t*g + beta*d + tmp*yk */
8129566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk));
813484c7b14SAdam Denchfield } else {
814484c7b14SAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless DFP step */
8159566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
8169566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
817484c7b14SAdam Denchfield /* compute scalar gamma */
8189566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
8199566063dSJacob Faibussowitsch PetscCall(VecDot(cg->y_work, cg->yk, &tmp));
820484c7b14SAdam Denchfield gamma = (gkp1_yk / tmp);
821484c7b14SAdam Denchfield /* Compute scalar beta */
822484c7b14SAdam Denchfield beta = -step * gd / dk_yk;
823484c7b14SAdam Denchfield /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
8249566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work));
825484c7b14SAdam Denchfield }
826484c7b14SAdam Denchfield break;
827484c7b14SAdam Denchfield
828d6e07cdcSHong Zhang case TAO_BNCG_SSML_BRDN:
8299566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
8309566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution));
831484c7b14SAdam Denchfield snorm = step * dnorm;
832484c7b14SAdam Denchfield cg->yts = step * dk_yk;
833484c7b14SAdam Denchfield cg->yty = ynorm2;
834484c7b14SAdam Denchfield cg->sts = snorm * snorm;
835484c7b14SAdam Denchfield if (!cg->diag_scaling) {
836484c7b14SAdam Denchfield /* Instead of a regular convex combination, we will solve a quadratic formula. */
8379566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(cg->yty, step * dk_yk, snorm * snorm, &tau_bfgs, cg->bfgs_scale));
8389566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(cg->yty, step * dk_yk, snorm * snorm, &tau_dfp, cg->dfp_scale));
8399566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk));
840484c7b14SAdam Denchfield tau_k = cg->theta * tau_bfgs + (1.0 - cg->theta) * tau_dfp;
841484c7b14SAdam Denchfield /* If bfgs_scale = 1.0, it should reproduce the bfgs tau_bfgs. If bfgs_scale = 0.0,
842484c7b14SAdam Denchfield it should reproduce the tau_dfp scaling. Same with dfp_scale. */
843484c7b14SAdam Denchfield tmp = cg->theta * tau_bfgs * gd / dk_yk + (1 - cg->theta) * tau_dfp * gkp1_yk / cg->yty;
844484c7b14SAdam Denchfield beta = cg->theta * tau_bfgs * (gkp1_yk / dk_yk - cg->yty * gd / (dk_yk * dk_yk)) - step * gd / dk_yk;
845484c7b14SAdam Denchfield /* d <- -t*g + beta*d + tmp*yk */
8469566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk));
847484c7b14SAdam Denchfield } else {
848484c7b14SAdam Denchfield /* We have diagonal scaling enabled */
8499566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work));
8509566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, cg->yk, cg->y_work));
851484c7b14SAdam Denchfield /* compute scalar gamma */
8529566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk));
8539566063dSJacob Faibussowitsch PetscCall(VecDot(cg->y_work, cg->yk, &tmp));
854484c7b14SAdam Denchfield gamma = cg->theta * gd / dk_yk + (1 - cg->theta) * (gkp1_yk / tmp);
855484c7b14SAdam Denchfield /* Compute scalar beta */
856484c7b14SAdam Denchfield beta = cg->theta * (gkp1_yk / dk_yk - gd * tmp / (dk_yk * dk_yk)) - step * gd / dk_yk;
857484c7b14SAdam Denchfield /* Compute stepdirection dkp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
8589566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work));
859484c7b14SAdam Denchfield }
860484c7b14SAdam Denchfield break;
861484c7b14SAdam Denchfield
862d71ae5a4SJacob Faibussowitsch default:
863d71ae5a4SJacob Faibussowitsch break;
864c8bcdf1eSAdam Denchfield }
865c8bcdf1eSAdam Denchfield }
8663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
867c8bcdf1eSAdam Denchfield }
868c8bcdf1eSAdam Denchfield
TaoBNCGConductIteration(Tao tao,PetscReal gnorm)869d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao tao, PetscReal gnorm)
870d71ae5a4SJacob Faibussowitsch {
871c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG *)tao->data;
872c8bcdf1eSAdam Denchfield TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
8738ca2df50S PetscReal step = 1.0, gnorm2, gd, dnorm = 0.0;
874c8bcdf1eSAdam Denchfield PetscReal gnorm2_old, f_old, resnorm, gnorm_old;
875c624ebd3SAlp Dener PetscBool pcgd_fallback = PETSC_FALSE;
876c8bcdf1eSAdam Denchfield
877c8bcdf1eSAdam Denchfield PetscFunctionBegin;
878c8bcdf1eSAdam Denchfield /* We are now going to perform a line search along the direction. */
879c8bcdf1eSAdam Denchfield /* Store solution and gradient info before it changes */
8809566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, cg->X_old));
8819566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, cg->G_old));
8829566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->unprojected_gradient, cg->unprojected_gradient_old));
883c8bcdf1eSAdam Denchfield
884c8bcdf1eSAdam Denchfield gnorm_old = gnorm;
885c8bcdf1eSAdam Denchfield gnorm2_old = gnorm_old * gnorm_old;
886c8bcdf1eSAdam Denchfield f_old = cg->f;
887484c7b14SAdam Denchfield /* Perform bounded line search. If we are recycling a solution from a previous */
888484c7b14SAdam Denchfield /* TaoSolve, then we want to immediately skip to calculating a new direction rather than performing a linesearch */
889414d97d3SAlp Dener if (!(tao->recycle && 0 == tao->niter)) {
890484c7b14SAdam Denchfield /* Above logic: the below code happens every iteration, except for the first iteration of a recycled TaoSolve */
8919566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
8929566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status));
8939566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao));
894c8bcdf1eSAdam Denchfield
895c8bcdf1eSAdam Denchfield /* Check linesearch failure */
896c8bcdf1eSAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
897c8bcdf1eSAdam Denchfield ++cg->ls_fails;
898d6e07cdcSHong Zhang if (cg->cg_type == TAO_BNCG_GD) {
899c8bcdf1eSAdam Denchfield /* Nothing left to do but fail out of the optimization */
900c8bcdf1eSAdam Denchfield step = 0.0;
901c8bcdf1eSAdam Denchfield tao->reason = TAO_DIVERGED_LS_FAILURE;
902c8bcdf1eSAdam Denchfield } else {
903484c7b14SAdam Denchfield /* Restore previous point, perform preconditioned GD and regular GD steps at the last good point */
9049566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->X_old, tao->solution));
9059566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->G_old, tao->gradient));
9069566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->unprojected_gradient_old, cg->unprojected_gradient));
907c8bcdf1eSAdam Denchfield gnorm = gnorm_old;
908c8bcdf1eSAdam Denchfield gnorm2 = gnorm2_old;
909c8bcdf1eSAdam Denchfield cg->f = f_old;
910c8bcdf1eSAdam Denchfield
911484c7b14SAdam Denchfield /* Fall back on preconditioned CG (so long as you're not already using it) */
912d6e07cdcSHong Zhang if (cg->cg_type != TAO_BNCG_PCGD && cg->diag_scaling) {
913e2570530SAlp Dener pcgd_fallback = PETSC_TRUE;
9149566063dSJacob Faibussowitsch PetscCall(TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback));
915484c7b14SAdam Denchfield
9169566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
9179566063dSJacob Faibussowitsch PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
918c8bcdf1eSAdam Denchfield
9199566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
9209566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status));
9219566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao));
922c8bcdf1eSAdam Denchfield
923484c7b14SAdam Denchfield pcgd_fallback = PETSC_FALSE;
924484c7b14SAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
925484c7b14SAdam Denchfield /* Going to perform a regular gradient descent step. */
926484c7b14SAdam Denchfield ++cg->ls_fails;
927484c7b14SAdam Denchfield step = 0.0;
928484c7b14SAdam Denchfield }
929484c7b14SAdam Denchfield }
930484c7b14SAdam Denchfield /* Fall back on the scaled gradient step */
931484c7b14SAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
932484c7b14SAdam Denchfield ++cg->ls_fails;
9339566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
9349566063dSJacob Faibussowitsch PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
9359566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
9369566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status));
9379566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao));
938484c7b14SAdam Denchfield }
939484c7b14SAdam Denchfield
940c8bcdf1eSAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
941c8bcdf1eSAdam Denchfield /* Nothing left to do but fail out of the optimization */
94250b47da0SAdam Denchfield ++cg->ls_fails;
943c8bcdf1eSAdam Denchfield step = 0.0;
944c8bcdf1eSAdam Denchfield tao->reason = TAO_DIVERGED_LS_FAILURE;
945484c7b14SAdam Denchfield } else {
946484c7b14SAdam Denchfield /* One of the fallbacks worked. Set them both back equal to false. */
947484c7b14SAdam Denchfield pcgd_fallback = PETSC_FALSE;
948c8bcdf1eSAdam Denchfield }
949c8bcdf1eSAdam Denchfield }
950c8bcdf1eSAdam Denchfield }
951c8bcdf1eSAdam Denchfield /* Convergence test for line search failure */
9523ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
953c8bcdf1eSAdam Denchfield
954c8bcdf1eSAdam Denchfield /* Standard convergence test */
9559566063dSJacob Faibussowitsch PetscCall(VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W));
9569566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->W, NORM_2, &resnorm));
957*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
9589566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its));
9599566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step));
960dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
9613ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
962484c7b14SAdam Denchfield }
963c8bcdf1eSAdam Denchfield /* Assert we have an updated step and we need at least one more iteration. */
964c8bcdf1eSAdam Denchfield /* Calculate the next direction */
965c8bcdf1eSAdam Denchfield /* Estimate the active set at the new solution */
9669566063dSJacob Faibussowitsch PetscCall(TaoBNCGEstimateActiveSet(tao, cg->as_type));
967c8bcdf1eSAdam Denchfield /* Compute the projected gradient and its norm */
9689566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->unprojected_gradient, tao->gradient));
969976ed0a4SStefano Zampini if (cg->active_idx) PetscCall(VecISSet(tao->gradient, cg->active_idx, 0.0));
9709566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
971c8bcdf1eSAdam Denchfield gnorm2 = gnorm * gnorm;
972c8bcdf1eSAdam Denchfield
973484c7b14SAdam Denchfield /* Calculate some quantities used in the StepDirectionUpdate. */
9749566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection, NORM_2, &dnorm));
975484c7b14SAdam Denchfield /* Update the step direction. */
9769566063dSJacob Faibussowitsch PetscCall(TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback));
977484c7b14SAdam Denchfield ++tao->niter;
9789566063dSJacob Faibussowitsch PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
979c8bcdf1eSAdam Denchfield
980d6e07cdcSHong Zhang if (cg->cg_type != TAO_BNCG_GD) {
981c8bcdf1eSAdam Denchfield /* Figure out which previously active variables became inactive this iteration */
9829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->new_inactives));
98348a46eb9SPierre Jolivet if (cg->inactive_idx && cg->inactive_old) PetscCall(ISDifference(cg->inactive_idx, cg->inactive_old, &cg->new_inactives));
984c8bcdf1eSAdam Denchfield /* Selectively reset the CG step those freshly inactive variables */
985c8bcdf1eSAdam Denchfield if (cg->new_inactives) {
9869566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step));
9879566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad));
9889566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->inactive_grad, cg->inactive_step));
9899566063dSJacob Faibussowitsch PetscCall(VecScale(cg->inactive_step, -1.0));
9909566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step));
9919566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad));
992c8bcdf1eSAdam Denchfield }
993c8bcdf1eSAdam Denchfield /* Verify that this is a descent direction */
9949566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
9959566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection, NORM_2, &dnorm));
99650b47da0SAdam Denchfield if (PetscIsInfOrNanReal(gd) || (gd / (dnorm * dnorm) <= -1e10 || gd / (dnorm * dnorm) >= -1e-10)) {
997c8bcdf1eSAdam Denchfield /* Not a descent direction, so we reset back to projected gradient descent */
9989566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2));
9999566063dSJacob Faibussowitsch PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection));
1000c8bcdf1eSAdam Denchfield ++cg->descent_error;
1001c8bcdf1eSAdam Denchfield } else {
1002c8bcdf1eSAdam Denchfield }
1003c8bcdf1eSAdam Denchfield }
10043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1005ac9112b8SAlp Dener }
1006484c7b14SAdam Denchfield
TaoBNCGSetH0(Tao tao,Mat H0)1007d6e07cdcSHong Zhang PETSC_INTERN PetscErrorCode TaoBNCGSetH0(Tao tao, Mat H0)
1008d71ae5a4SJacob Faibussowitsch {
1009484c7b14SAdam Denchfield TAO_BNCG *cg = (TAO_BNCG *)tao->data;
1010d6e07cdcSHong Zhang PetscBool same;
1011484c7b14SAdam Denchfield
1012484c7b14SAdam Denchfield PetscFunctionBegin;
1013d6e07cdcSHong Zhang PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBNCG, &same));
1014d6e07cdcSHong Zhang if (same) {
10159566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H0));
1016484c7b14SAdam Denchfield cg->pc = H0;
1017d6e07cdcSHong Zhang }
1018d6e07cdcSHong Zhang PetscFunctionReturn(PETSC_SUCCESS);
1019d6e07cdcSHong Zhang }
1020d6e07cdcSHong Zhang
1021d6e07cdcSHong Zhang /*@
1022d6e07cdcSHong Zhang TaoBNCGGetType - Return the type for the `TAOBNCG` solver
1023d6e07cdcSHong Zhang
1024d6e07cdcSHong Zhang Input Parameter:
1025d6e07cdcSHong Zhang . tao - the `Tao` solver context
1026d6e07cdcSHong Zhang
1027d6e07cdcSHong Zhang Output Parameter:
1028d6e07cdcSHong Zhang . type - `TAOBNCG` type
1029d6e07cdcSHong Zhang
1030d6e07cdcSHong Zhang Level: advanced
1031d6e07cdcSHong Zhang
1032d6e07cdcSHong Zhang .seealso: `Tao`, `TAOBNCG`, `TaoBNCGSetType()`, `TaoBNCGType`
1033d6e07cdcSHong Zhang @*/
TaoBNCGGetType(Tao tao,TaoBNCGType * type)1034d6e07cdcSHong Zhang PetscErrorCode TaoBNCGGetType(Tao tao, TaoBNCGType *type)
1035d6e07cdcSHong Zhang {
1036d6e07cdcSHong Zhang TAO_BNCG *cg = (TAO_BNCG *)tao->data;
1037d6e07cdcSHong Zhang PetscBool same;
1038d6e07cdcSHong Zhang
1039d6e07cdcSHong Zhang PetscFunctionBegin;
1040d6e07cdcSHong Zhang PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBNCG, &same));
1041d6e07cdcSHong Zhang PetscCheck(same, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "TAO solver is not BNCG type");
1042d6e07cdcSHong Zhang *type = cg->cg_type;
1043d6e07cdcSHong Zhang PetscFunctionReturn(PETSC_SUCCESS);
1044d6e07cdcSHong Zhang }
1045d6e07cdcSHong Zhang
1046d6e07cdcSHong Zhang /*@
1047d6e07cdcSHong Zhang TaoBNCGSetType - Set the type for the `TAOBNCG` solver
1048d6e07cdcSHong Zhang
1049d6e07cdcSHong Zhang Input Parameters:
1050d6e07cdcSHong Zhang + tao - the `Tao` solver context
1051d6e07cdcSHong Zhang - type - `TAOBNCG` type
1052d6e07cdcSHong Zhang
1053d6e07cdcSHong Zhang Level: advanced
1054d6e07cdcSHong Zhang
1055d6e07cdcSHong Zhang .seealso: `Tao`, `TAOBNCG`, `TaoBNCGGetType()`, `TaoBNCGType`
1056d6e07cdcSHong Zhang @*/
TaoBNCGSetType(Tao tao,TaoBNCGType type)1057d6e07cdcSHong Zhang PetscErrorCode TaoBNCGSetType(Tao tao, TaoBNCGType type)
1058d6e07cdcSHong Zhang {
1059d6e07cdcSHong Zhang TAO_BNCG *cg = (TAO_BNCG *)tao->data;
1060d6e07cdcSHong Zhang PetscBool same;
1061d6e07cdcSHong Zhang
1062d6e07cdcSHong Zhang PetscFunctionBegin;
1063d6e07cdcSHong Zhang PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBNCG, &same));
1064d6e07cdcSHong Zhang if (same) cg->cg_type = type;
10653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1066484c7b14SAdam Denchfield }
1067