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 5*d6e07cdcSHong 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 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 3161be54a6SAlp Dener case CG_AS_BERTSEKAS: 3261be54a6SAlp Dener /* Use gradient descent to estimate the active set */ 339566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->unprojected_gradient, cg->W)); 349566063dSJacob Faibussowitsch PetscCall(VecScale(cg->W, -1.0)); 359371c9d4SSatish 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)); 36c4b75bccSAlp Dener break; 37d71ae5a4SJacob Faibussowitsch default: 38d71ae5a4SJacob Faibussowitsch break; 3961be54a6SAlp Dener } 403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4161be54a6SAlp Dener } 4261be54a6SAlp Dener 43d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNCGBoundStep(Tao tao, PetscInt asType, Vec step) 44d71ae5a4SJacob Faibussowitsch { 4561be54a6SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 4661be54a6SAlp Dener 4761be54a6SAlp Dener PetscFunctionBegin; 48a1318120SAlp Dener switch (asType) { 49d71ae5a4SJacob Faibussowitsch case CG_AS_NONE: 50d71ae5a4SJacob Faibussowitsch PetscCall(VecISSet(step, cg->active_idx, 0.0)); 51d71ae5a4SJacob Faibussowitsch break; 5261be54a6SAlp Dener 53d71ae5a4SJacob Faibussowitsch case CG_AS_BERTSEKAS: 54d71ae5a4SJacob Faibussowitsch PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, cg->active_lower, cg->active_upper, cg->active_fixed, 1.0, step)); 55d71ae5a4SJacob Faibussowitsch break; 5661be54a6SAlp Dener 57d71ae5a4SJacob Faibussowitsch default: 58d71ae5a4SJacob Faibussowitsch break; 5961be54a6SAlp Dener } 603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6161be54a6SAlp Dener } 6261be54a6SAlp Dener 63d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_BNCG(Tao tao) 64d71ae5a4SJacob Faibussowitsch { 65ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 66484c7b14SAdam Denchfield PetscReal step = 1.0, gnorm, gnorm2, resnorm; 67c4b75bccSAlp Dener PetscInt nDiff; 68ac9112b8SAlp Dener 69ac9112b8SAlp Dener PetscFunctionBegin; 70ac9112b8SAlp Dener /* Project the current point onto the feasible set */ 719566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 729566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU)); 73ac9112b8SAlp Dener 74c8bcdf1eSAdam Denchfield /* Project the initial point onto the feasible region */ 759566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution)); 76484c7b14SAdam Denchfield 7748a46eb9SPierre Jolivet if (nDiff > 0 || !tao->recycle) PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient)); 789566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->unprojected_gradient, NORM_2, &gnorm)); 793c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(cg->f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 80ac9112b8SAlp Dener 8161be54a6SAlp Dener /* Estimate the active set and compute the projected gradient */ 829566063dSJacob Faibussowitsch PetscCall(TaoBNCGEstimateActiveSet(tao, cg->as_type)); 8361be54a6SAlp Dener 84ac9112b8SAlp Dener /* Project the gradient and calculate the norm */ 859566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->unprojected_gradient, tao->gradient)); 869566063dSJacob Faibussowitsch PetscCall(VecISSet(tao->gradient, cg->active_idx, 0.0)); 879566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm)); 88ac9112b8SAlp Dener gnorm2 = gnorm * gnorm; 89ac9112b8SAlp Dener 90c8bcdf1eSAdam Denchfield /* Initialize counters */ 91e031d6f5SAlp Dener tao->niter = 0; 9250b47da0SAdam Denchfield cg->ls_fails = cg->descent_error = 0; 93c8bcdf1eSAdam Denchfield cg->resets = -1; 94484c7b14SAdam Denchfield cg->skipped_updates = cg->pure_gd_steps = 0; 95c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 96c8bcdf1eSAdam Denchfield 97c8bcdf1eSAdam Denchfield /* Convergence test at the starting point. */ 98ac9112b8SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 99484c7b14SAdam Denchfield 1009566063dSJacob Faibussowitsch PetscCall(VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W)); 1019566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->W, NORM_2, &resnorm)); 1023c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 1039566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its)); 1049566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step)); 105dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 1063ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 107484c7b14SAdam Denchfield /* Calculate initial direction. */ 108414d97d3SAlp Dener if (!tao->recycle) { 109484c7b14SAdam Denchfield /* We are not recycling a solution/history from a past TaoSolve */ 1109566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2)); 111ac9112b8SAlp Dener } 112c8bcdf1eSAdam Denchfield /* Initial gradient descent step. Scaling by 1.0 also does a decent job for some problems. */ 113c8bcdf1eSAdam Denchfield while (1) { 114e1e80dc8SAlp Dener /* Call general purpose update function */ 115e1e80dc8SAlp Dener if (tao->ops->update) { 116dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, update, tao->niter, tao->user_update); 1177494f0b1SStefano Zampini PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient)); 118e1e80dc8SAlp Dener } 1199566063dSJacob Faibussowitsch PetscCall(TaoBNCGConductIteration(tao, gnorm)); 1203ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 121ac9112b8SAlp Dener } 122ac9112b8SAlp Dener } 123ac9112b8SAlp Dener 124d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_BNCG(Tao tao) 125d71ae5a4SJacob Faibussowitsch { 126ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 127ac9112b8SAlp Dener 128ac9112b8SAlp Dener PetscFunctionBegin; 12948a46eb9SPierre Jolivet if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 13048a46eb9SPierre Jolivet if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 13148a46eb9SPierre Jolivet if (!cg->W) PetscCall(VecDuplicate(tao->solution, &cg->W)); 13248a46eb9SPierre Jolivet if (!cg->work) PetscCall(VecDuplicate(tao->solution, &cg->work)); 13348a46eb9SPierre Jolivet if (!cg->sk) PetscCall(VecDuplicate(tao->solution, &cg->sk)); 13448a46eb9SPierre Jolivet if (!cg->yk) PetscCall(VecDuplicate(tao->gradient, &cg->yk)); 13548a46eb9SPierre Jolivet if (!cg->X_old) PetscCall(VecDuplicate(tao->solution, &cg->X_old)); 13648a46eb9SPierre Jolivet if (!cg->G_old) PetscCall(VecDuplicate(tao->gradient, &cg->G_old)); 137c8bcdf1eSAdam Denchfield if (cg->diag_scaling) { 1389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &cg->d_work)); 1399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &cg->y_work)); 1409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &cg->g_work)); 141c4b75bccSAlp Dener } 14248a46eb9SPierre Jolivet if (!cg->unprojected_gradient) PetscCall(VecDuplicate(tao->gradient, &cg->unprojected_gradient)); 14348a46eb9SPierre Jolivet if (!cg->unprojected_gradient_old) PetscCall(VecDuplicate(tao->gradient, &cg->unprojected_gradient_old)); 1449566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(cg->B, cg->sk, cg->yk)); 1451baa6e33SBarry Smith if (cg->pc) PetscCall(MatLMVMSetJ0(cg->B, cg->pc)); 1463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 147ac9112b8SAlp Dener } 148ac9112b8SAlp Dener 149d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_BNCG(Tao tao) 150d71ae5a4SJacob Faibussowitsch { 151ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 152ac9112b8SAlp Dener 153ac9112b8SAlp Dener PetscFunctionBegin; 154ac9112b8SAlp Dener if (tao->setupcalled) { 1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->W)); 1569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->work)); 1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->X_old)); 1589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->G_old)); 1599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->unprojected_gradient)); 1609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->unprojected_gradient_old)); 1619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->g_work)); 1629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->d_work)); 1639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->y_work)); 1649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->sk)); 1659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cg->yk)); 166ac9112b8SAlp Dener } 1679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->active_lower)); 1689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->active_upper)); 1699566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->active_fixed)); 1709566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->active_idx)); 1719566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->inactive_idx)); 1729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->inactive_old)); 1739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->new_inactives)); 1749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&cg->B)); 17548a46eb9SPierre Jolivet if (cg->pc) PetscCall(MatDestroy(&cg->pc)); 1769566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 178ac9112b8SAlp Dener } 179ac9112b8SAlp Dener 180d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetFromOptions_BNCG(Tao tao, PetscOptionItems *PetscOptionsObject) 181d71ae5a4SJacob Faibussowitsch { 182ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 183ac9112b8SAlp Dener 184ac9112b8SAlp Dener PetscFunctionBegin; 185d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Nonlinear Conjugate Gradient method for unconstrained optimization"); 186*d6e07cdcSHong Zhang PetscCall(PetscOptionsEnum("-tao_bncg_type", "CG update formula", "TaoBNCGTypes", TaoBNCGTypes, (PetscEnum)cg->cg_type, (PetscEnum *)&cg->cg_type, NULL)); 187*d6e07cdcSHong Zhang if (cg->cg_type != TAO_BNCG_SSML_BFGS) cg->alpha = -1.0; /* Setting defaults for non-BFGS methods. User can change it below. */ 188*d6e07cdcSHong Zhang if (TAO_BNCG_GD == cg->cg_type) { 189*d6e07cdcSHong Zhang cg->cg_type = TAO_BNCG_PCGD; 190484c7b14SAdam Denchfield /* Set scaling equal to none or, at best, scalar scaling. */ 191484c7b14SAdam Denchfield cg->unscaled_restart = PETSC_TRUE; 192484c7b14SAdam Denchfield cg->diag_scaling = PETSC_FALSE; 193484c7b14SAdam Denchfield } 1949566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_hz_eta", "(developer) cutoff tolerance for HZ", "", cg->hz_eta, &cg->hz_eta, NULL)); 1959566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_eps", "(developer) cutoff value for restarts", "", cg->epsilon, &cg->epsilon, NULL)); 1969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_dk_eta", "(developer) cutoff tolerance for DK", "", cg->dk_eta, &cg->dk_eta, NULL)); 1979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_xi", "(developer) Parameter in the KD method", "", cg->xi, &cg->xi, NULL)); 1989566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_theta", "(developer) update parameter for the Broyden method", "", cg->theta, &cg->theta, NULL)); 1999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_hz_theta", "(developer) parameter for the HZ (2006) method", "", cg->hz_theta, &cg->hz_theta, NULL)); 2009566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_alpha", "(developer) parameter for the scalar scaling", "", cg->alpha, &cg->alpha, NULL)); 2019566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_bfgs_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->bfgs_scale, &cg->bfgs_scale, NULL)); 2029566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_dfp_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->dfp_scale, &cg->dfp_scale, NULL)); 2039566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_bncg_diag_scaling", "Enable diagonal Broyden-like preconditioning", "", cg->diag_scaling, &cg->diag_scaling, NULL)); 2049566063dSJacob 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)); 2059566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_bncg_unscaled_restart", "(developer) use unscaled gradient restarts", "", cg->unscaled_restart, &cg->unscaled_restart, NULL)); 2069566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_zeta", "(developer) Free parameter for the Kou-Dai method", "", cg->zeta, &cg->zeta, NULL)); 2079566063dSJacob 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)); 2089566063dSJacob 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)); 2099566063dSJacob 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)); 2109566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_bncg_no_scaling", "Disable all scaling except in restarts", "", cg->no_scaling, &cg->no_scaling, NULL)); 211484c7b14SAdam Denchfield if (cg->no_scaling) { 212484c7b14SAdam Denchfield cg->diag_scaling = PETSC_FALSE; 213484c7b14SAdam Denchfield cg->alpha = -1.0; 214484c7b14SAdam Denchfield } 215*d6e07cdcSHong Zhang if (cg->alpha == -1.0 && cg->cg_type == TAO_BNCG_KD && !cg->diag_scaling) { /* Some more default options that appear to be good. */ 216484c7b14SAdam Denchfield cg->neg_xi = PETSC_TRUE; 217484c7b14SAdam Denchfield } 2189566063dSJacob 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)); 219*d6e07cdcSHong 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)); 2209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", cg->as_tol, &cg->as_tol, NULL)); 2219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_as_step", "(developer) step length used when estimating actively bounded variables", "", cg->as_step, &cg->as_step, NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_delta_min", "(developer) minimum scaling factor used for scaled gradient restarts", "", cg->delta_min, &cg->delta_min, NULL)); 2239566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bncg_delta_max", "(developer) maximum scaling factor used for scaled gradient restarts", "", cg->delta_max, &cg->delta_max, NULL)); 22450b47da0SAdam Denchfield 225d0609cedSBarry Smith PetscOptionsHeadEnd(); 2269566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(cg->B, ((PetscObject)tao)->prefix)); 2279566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(cg->B, "tao_bncg_")); 2289566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(cg->B)); 2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 230ac9112b8SAlp Dener } 231ac9112b8SAlp Dener 232d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_BNCG(Tao tao, PetscViewer viewer) 233d71ae5a4SJacob Faibussowitsch { 234ac9112b8SAlp Dener PetscBool isascii; 235ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 236ac9112b8SAlp Dener 237ac9112b8SAlp Dener PetscFunctionBegin; 2389566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 239ac9112b8SAlp Dener if (isascii) { 2409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 241*d6e07cdcSHong Zhang PetscCall(PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", TaoBNCGTypes[cg->cg_type])); 24263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Skipped Stepdirection Updates: %" PetscInt_FMT "\n", cg->skipped_updates)); 24363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", cg->resets)); 24463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Pure gradient steps: %" PetscInt_FMT "\n", cg->pure_gd_steps)); 24563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Not a descent direction: %" PetscInt_FMT "\n", cg->descent_error)); 24663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Line search fails: %" PetscInt_FMT "\n", cg->ls_fails)); 247484c7b14SAdam Denchfield if (cg->diag_scaling) { 2489566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 249484c7b14SAdam Denchfield if (isascii) { 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 2519566063dSJacob Faibussowitsch PetscCall(MatView(cg->B, viewer)); 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 253484c7b14SAdam Denchfield } 254484c7b14SAdam Denchfield } 2559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 256ac9112b8SAlp Dener } 2573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 258ac9112b8SAlp Dener } 259ac9112b8SAlp Dener 260d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal yty, PetscReal yts, PetscReal sts, PetscReal *scale, PetscReal alpha) 261d71ae5a4SJacob Faibussowitsch { 262c8bcdf1eSAdam Denchfield PetscReal a, b, c, sig1, sig2; 263c8bcdf1eSAdam Denchfield 264c8bcdf1eSAdam Denchfield PetscFunctionBegin; 265c8bcdf1eSAdam Denchfield *scale = 0.0; 2668ebe3e4eSStefano Zampini if (1.0 == alpha) *scale = yts / yty; 2678ebe3e4eSStefano Zampini else if (0.0 == alpha) *scale = sts / yts; 26850b47da0SAdam Denchfield else if (-1.0 == alpha) *scale = 1.0; 269c8bcdf1eSAdam Denchfield else { 270c8bcdf1eSAdam Denchfield a = yty; 271c8bcdf1eSAdam Denchfield b = yts; 272c8bcdf1eSAdam Denchfield c = sts; 273c8bcdf1eSAdam Denchfield a *= alpha; 274c8bcdf1eSAdam Denchfield b *= -(2.0 * alpha - 1.0); 275c8bcdf1eSAdam Denchfield c *= alpha - 1.0; 276c8bcdf1eSAdam Denchfield sig1 = (-b + PetscSqrtReal(b * b - 4.0 * a * c)) / (2.0 * a); 277c8bcdf1eSAdam Denchfield sig2 = (-b - PetscSqrtReal(b * b - 4.0 * a * c)) / (2.0 * a); 278c8bcdf1eSAdam Denchfield /* accept the positive root as the scalar */ 2798ebe3e4eSStefano Zampini if (sig1 > 0.0) *scale = sig1; 2808ebe3e4eSStefano Zampini else if (sig2 > 0.0) *scale = sig2; 2818ebe3e4eSStefano Zampini else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Cannot find positive scalar"); 282c8bcdf1eSAdam Denchfield } 2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 284c8bcdf1eSAdam Denchfield } 285c8bcdf1eSAdam Denchfield 286ac9112b8SAlp Dener /*MC 287ac9112b8SAlp Dener TAOBNCG - Bound-constrained Nonlinear Conjugate Gradient method. 288ac9112b8SAlp Dener 289ac9112b8SAlp Dener Options Database Keys: 29050b47da0SAdam Denchfield + -tao_bncg_recycle - enable recycling the latest calculated gradient vector in subsequent TaoSolve() calls (currently disabled) 291c4b75bccSAlp Dener . -tao_bncg_eta <r> - restart tolerance 29261be54a6SAlp Dener . -tao_bncg_type <taocg_type> - cg formula 293c4b75bccSAlp Dener . -tao_bncg_as_type <none,bertsekas> - active set estimation method 294c4b75bccSAlp Dener . -tao_bncg_as_tol <r> - tolerance used in Bertsekas active-set estimation 295c4b75bccSAlp Dener . -tao_bncg_as_step <r> - trial step length used in Bertsekas active-set estimation 29650b47da0SAdam 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. 29750b47da0SAdam Denchfield . -tao_bncg_theta <r> - convex combination parameter for the Broyden method 29850b47da0SAdam Denchfield . -tao_bncg_hz_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods 29950b47da0SAdam Denchfield . -tao_bncg_dk_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods 30050b47da0SAdam Denchfield . -tao_bncg_xi <r> - Multiplicative constant of the gamma term in the KD method 30150b47da0SAdam Denchfield . -tao_bncg_hz_theta <r> - Multiplicative constant of the theta term for the HZ method 30250b47da0SAdam Denchfield . -tao_bncg_bfgs_scale <r> - Scaling parameter of the bfgs contribution to the scalar Broyden method 30350b47da0SAdam Denchfield . -tao_bncg_dfp_scale <r> - Scaling parameter of the dfp contribution to the scalar Broyden method 30450b47da0SAdam Denchfield . -tao_bncg_diag_scaling <b> - Whether or not to use diagonal initialization/preconditioning for the CG methods. Default True. 30550b47da0SAdam Denchfield . -tao_bncg_dynamic_restart <b> - use dynamic restart strategy in the HZ, DK, KD methods 30650b47da0SAdam Denchfield . -tao_bncg_unscaled_restart <b> - whether or not to scale the gradient when doing gradient descent restarts 30750b47da0SAdam Denchfield . -tao_bncg_zeta <r> - Scaling parameter in the KD method 308484c7b14SAdam Denchfield . -tao_bncg_delta_min <r> - Minimum bound for rescaling during restarted gradient descent steps 309484c7b14SAdam Denchfield . -tao_bncg_delta_max <r> - Maximum bound for rescaling during restarted gradient descent steps 31050b47da0SAdam Denchfield . -tao_bncg_min_quad <i> - Number of quadratic-like steps in a row necessary to do a dynamic restart 31150b47da0SAdam 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 31250b47da0SAdam Denchfield . -tao_bncg_spaced_restart <b> - whether or not to do gradient descent steps every x*n iterations 313484c7b14SAdam Denchfield . -tao_bncg_no_scaling <b> - If true, eliminates all scaling, including defaults. 3143850be85SAlp Dener - -tao_bncg_neg_xi <b> - Whether or not to use negative xi in the KD method under certain conditions 315ac9112b8SAlp Dener 316ac9112b8SAlp Dener Notes: 317ac9112b8SAlp Dener CG formulas are: 3183850be85SAlp Dener + "gd" - Gradient Descent 3193850be85SAlp Dener . "fr" - Fletcher-Reeves 3203850be85SAlp Dener . "pr" - Polak-Ribiere-Polyak 3213850be85SAlp Dener . "prp" - Polak-Ribiere-Plus 3223850be85SAlp Dener . "hs" - Hestenes-Steifel 3233850be85SAlp Dener . "dy" - Dai-Yuan 3243850be85SAlp Dener . "ssml_bfgs" - Self-Scaling Memoryless BFGS 3253850be85SAlp Dener . "ssml_dfp" - Self-Scaling Memoryless DFP 3263850be85SAlp Dener . "ssml_brdn" - Self-Scaling Memoryless Broyden 3273850be85SAlp Dener . "hz" - Hager-Zhang (CG_DESCENT 5.3) 3283850be85SAlp Dener . "dk" - Dai-Kou (2013) 3293850be85SAlp Dener - "kd" - Kou-Dai (2015) 3309abc5736SPatrick Sanan 331ac9112b8SAlp Dener Level: beginner 332ac9112b8SAlp Dener M*/ 333ac9112b8SAlp Dener 334d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BNCG(Tao tao) 335d71ae5a4SJacob Faibussowitsch { 336ac9112b8SAlp Dener TAO_BNCG *cg; 337ac9112b8SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 338ac9112b8SAlp Dener 339ac9112b8SAlp Dener PetscFunctionBegin; 340ac9112b8SAlp Dener tao->ops->setup = TaoSetUp_BNCG; 341ac9112b8SAlp Dener tao->ops->solve = TaoSolve_BNCG; 342ac9112b8SAlp Dener tao->ops->view = TaoView_BNCG; 343ac9112b8SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNCG; 344ac9112b8SAlp Dener tao->ops->destroy = TaoDestroy_BNCG; 345ac9112b8SAlp Dener 346ac9112b8SAlp Dener /* Override default settings (unless already changed) */ 347ac9112b8SAlp Dener if (!tao->max_it_changed) tao->max_it = 2000; 348ac9112b8SAlp Dener if (!tao->max_funcs_changed) tao->max_funcs = 4000; 349ac9112b8SAlp Dener 350ac9112b8SAlp Dener /* Note: nondefault values should be used for nonlinear conjugate gradient */ 351ac9112b8SAlp Dener /* method. In particular, gtol should be less that 0.5; the value used in */ 352ac9112b8SAlp Dener /* Nocedal and Wright is 0.10. We use the default values for the */ 353ac9112b8SAlp Dener /* linesearch because it seems to work better. */ 3549566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 3559566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 3569566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type)); 3579566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao)); 358ac9112b8SAlp Dener 3594dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&cg)); 360ac9112b8SAlp Dener tao->data = (void *)cg; 3619566063dSJacob Faibussowitsch PetscCall(KSPInitializePackage()); 3629566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &cg->B)); 3639566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)cg->B, (PetscObject)tao, 1)); 3649566063dSJacob Faibussowitsch PetscCall(MatSetType(cg->B, MATLMVMDIAGBROYDEN)); 36550b47da0SAdam Denchfield 366484c7b14SAdam Denchfield cg->pc = NULL; 367484c7b14SAdam Denchfield 36850b47da0SAdam Denchfield cg->dk_eta = 0.5; 36950b47da0SAdam Denchfield cg->hz_eta = 0.4; 370c8bcdf1eSAdam Denchfield cg->dynamic_restart = PETSC_FALSE; 371c8bcdf1eSAdam Denchfield cg->unscaled_restart = PETSC_FALSE; 372484c7b14SAdam Denchfield cg->no_scaling = PETSC_FALSE; 373484c7b14SAdam Denchfield cg->delta_min = 1e-7; 374484c7b14SAdam Denchfield cg->delta_max = 100; 375c8bcdf1eSAdam Denchfield cg->theta = 1.0; 376c8bcdf1eSAdam Denchfield cg->hz_theta = 1.0; 377c8bcdf1eSAdam Denchfield cg->dfp_scale = 1.0; 378c8bcdf1eSAdam Denchfield cg->bfgs_scale = 1.0; 37950b47da0SAdam Denchfield cg->zeta = 0.1; 38050b47da0SAdam Denchfield cg->min_quad = 6; 381c8bcdf1eSAdam Denchfield cg->min_restart_num = 6; /* As in CG_DESCENT and KD2015*/ 382c8bcdf1eSAdam Denchfield cg->xi = 1.0; 38350b47da0SAdam Denchfield cg->neg_xi = PETSC_TRUE; 384c8bcdf1eSAdam Denchfield cg->spaced_restart = PETSC_FALSE; 385c8bcdf1eSAdam Denchfield cg->tol_quad = 1e-8; 38661be54a6SAlp Dener cg->as_step = 0.001; 38761be54a6SAlp Dener cg->as_tol = 0.001; 38850b47da0SAdam Denchfield cg->eps_23 = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0); /* Just a little tighter*/ 38961be54a6SAlp Dener cg->as_type = CG_AS_BERTSEKAS; 390*d6e07cdcSHong Zhang cg->cg_type = TAO_BNCG_SSML_BFGS; 391c8bcdf1eSAdam Denchfield cg->alpha = 1.0; 392c8bcdf1eSAdam Denchfield cg->diag_scaling = PETSC_TRUE; 3933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 394c8bcdf1eSAdam Denchfield } 395c8bcdf1eSAdam Denchfield 396d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNCGResetUpdate(Tao tao, PetscReal gnormsq) 397d71ae5a4SJacob Faibussowitsch { 398c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG *)tao->data; 399c8bcdf1eSAdam Denchfield PetscReal scaling; 400c8bcdf1eSAdam Denchfield 401c8bcdf1eSAdam Denchfield PetscFunctionBegin; 402c8bcdf1eSAdam Denchfield ++cg->resets; 403c8bcdf1eSAdam Denchfield scaling = 2.0 * PetscMax(1.0, PetscAbsScalar(cg->f)) / PetscMax(gnormsq, cg->eps_23); 404484c7b14SAdam Denchfield scaling = PetscMin(cg->delta_max, PetscMax(cg->delta_min, scaling)); 405484c7b14SAdam Denchfield if (cg->unscaled_restart) { 406484c7b14SAdam Denchfield scaling = 1.0; 407484c7b14SAdam Denchfield ++cg->pure_gd_steps; 408484c7b14SAdam Denchfield } 4099566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -scaling, 0.0, tao->gradient)); 410c8bcdf1eSAdam Denchfield /* Also want to reset our diagonal scaling with each restart */ 4111baa6e33SBarry Smith if (cg->diag_scaling) PetscCall(MatLMVMReset(cg->B, PETSC_FALSE)); 4123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 413c8bcdf1eSAdam Denchfield } 414c8bcdf1eSAdam Denchfield 415d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNCGCheckDynamicRestart(Tao tao, PetscReal stepsize, PetscReal gd, PetscReal gd_old, PetscBool *dynrestart, PetscReal fold) 416d71ae5a4SJacob Faibussowitsch { 417c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG *)tao->data; 418c8bcdf1eSAdam Denchfield PetscReal quadinterp; 419c8bcdf1eSAdam Denchfield 420c8bcdf1eSAdam Denchfield PetscFunctionBegin; 42150b47da0SAdam Denchfield if (cg->f < cg->min_quad / 10) { 42250b47da0SAdam Denchfield *dynrestart = PETSC_FALSE; 4233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42450b47da0SAdam Denchfield } /* just skip this since this strategy doesn't work well for functions near zero */ 425484c7b14SAdam Denchfield quadinterp = 2.0 * (cg->f - fold) / (stepsize * (gd + gd_old)); 42650b47da0SAdam Denchfield if (PetscAbs(quadinterp - 1.0) < cg->tol_quad) ++cg->iter_quad; 427c8bcdf1eSAdam Denchfield else { 428c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 429c8bcdf1eSAdam Denchfield *dynrestart = PETSC_FALSE; 430c8bcdf1eSAdam Denchfield } 431c8bcdf1eSAdam Denchfield if (cg->iter_quad >= cg->min_quad) { 432c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 433c8bcdf1eSAdam Denchfield *dynrestart = PETSC_TRUE; 434c8bcdf1eSAdam Denchfield } 4353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 436c8bcdf1eSAdam Denchfield } 437c8bcdf1eSAdam Denchfield 438d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao tao, PetscReal gnorm2, PetscReal step, PetscReal fold, PetscReal gnorm2_old, PetscReal dnorm, PetscBool pcgd_fallback) 439d71ae5a4SJacob Faibussowitsch { 440c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG *)tao->data; 44150b47da0SAdam Denchfield PetscReal gamma = 1.0, tau_k, beta; 442484c7b14SAdam Denchfield PetscReal tmp = 1.0, ynorm, ynorm2 = 1.0, snorm = 1.0, dk_yk = 1.0, gd; 44350b47da0SAdam Denchfield PetscReal gkp1_yk, gd_old, tau_bfgs, tau_dfp, gkp1D_yk, gtDg; 444c8bcdf1eSAdam Denchfield PetscInt dim; 445484c7b14SAdam Denchfield PetscBool cg_restart = PETSC_FALSE; 446c8bcdf1eSAdam Denchfield PetscFunctionBegin; 447c8bcdf1eSAdam Denchfield 44850b47da0SAdam Denchfield /* Local curvature check to see if we need to restart */ 449414d97d3SAlp Dener if (tao->niter >= 1 || tao->recycle) { 4509566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient)); 4519566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->yk, NORM_2, &ynorm)); 452c8bcdf1eSAdam Denchfield ynorm2 = ynorm * ynorm; 4539566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->stepdirection, &dk_yk)); 454484c7b14SAdam Denchfield if (step * dnorm < PETSC_MACHINE_EPSILON || step * dk_yk < PETSC_MACHINE_EPSILON) { 455e2570530SAlp Dener cg_restart = PETSC_TRUE; 456484c7b14SAdam Denchfield ++cg->skipped_updates; 457484c7b14SAdam Denchfield } 45850b47da0SAdam Denchfield if (cg->spaced_restart) { 4599566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->gradient, &dim)); 460e2570530SAlp Dener if (tao->niter % (dim * cg->min_restart_num)) cg_restart = PETSC_TRUE; 46150b47da0SAdam Denchfield } 46250b47da0SAdam Denchfield } 46350b47da0SAdam Denchfield /* If the user wants regular restarts, do it every 6n iterations, where n=dimension */ 46450b47da0SAdam Denchfield if (cg->spaced_restart) { 4659566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->gradient, &dim)); 466e2570530SAlp Dener if (0 == tao->niter % (6 * dim)) cg_restart = PETSC_TRUE; 46750b47da0SAdam Denchfield } 46850b47da0SAdam Denchfield /* Compute the diagonal scaling vector if applicable */ 4691baa6e33SBarry Smith if (cg->diag_scaling) PetscCall(MatLMVMUpdate(cg->B, tao->solution, tao->gradient)); 47050b47da0SAdam Denchfield 471484c7b14SAdam Denchfield /* A note on diagonal scaling (to be added to paper): 472484c7b14SAdam Denchfield For the FR, PR, PRP, and DY methods, the diagonally scaled versions 473484c7b14SAdam Denchfield must be derived as a preconditioned CG method rather than as 474484c7b14SAdam Denchfield a Hessian initialization like in the Broyden methods. */ 47550b47da0SAdam Denchfield 476484c7b14SAdam Denchfield /* In that case, one writes the objective function as 477484c7b14SAdam Denchfield f(x) \equiv f(Ay). Gradient evaluations yield g(x_k) = A g(Ay_k) = A g(x_k). 478484c7b14SAdam Denchfield Furthermore, the direction d_k \equiv (x_k - x_{k-1})/step according to 479484c7b14SAdam Denchfield HZ (2006) becomes A^{-1} d_k, such that d_k^T g_k remains the 480484c7b14SAdam Denchfield same under preconditioning. Note that A is diagonal, such that A^T = A. */ 48150b47da0SAdam Denchfield 482484c7b14SAdam Denchfield /* This yields questions like what the dot product d_k^T y_k 483484c7b14SAdam Denchfield should look like. HZ mistakenly treats that as the same under 484484c7b14SAdam Denchfield preconditioning, but that is not necessarily true. */ 48550b47da0SAdam Denchfield 486484c7b14SAdam Denchfield /* Observe y_k \equiv g_k - g_{k-1}, and under the P.C. transformation, 487484c7b14SAdam 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}), 488484c7b14SAdam 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 489484c7b14SAdam Denchfield NOT the same if our preconditioning matrix is updated between iterations. 490484c7b14SAdam Denchfield This same issue is found when considering dot products of the form g_{k+1}^T y_k. */ 49150b47da0SAdam Denchfield 49250b47da0SAdam Denchfield /* Compute CG step direction */ 49350b47da0SAdam Denchfield if (cg_restart) { 4949566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2)); 495484c7b14SAdam Denchfield } else if (pcgd_fallback) { 496484c7b14SAdam Denchfield /* Just like preconditioned CG */ 4979566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 4989566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work)); 49950b47da0SAdam Denchfield } else if (ynorm2 > PETSC_MACHINE_EPSILON) { 50050b47da0SAdam Denchfield switch (cg->cg_type) { 501*d6e07cdcSHong Zhang case TAO_BNCG_PCGD: 50250b47da0SAdam Denchfield if (!cg->diag_scaling) { 503484c7b14SAdam Denchfield if (!cg->no_scaling) { 50450b47da0SAdam Denchfield cg->sts = step * step * dnorm * dnorm; 5059566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, cg->sts, &tau_k, cg->alpha)); 506484c7b14SAdam Denchfield } else { 507484c7b14SAdam Denchfield tau_k = 1.0; 508484c7b14SAdam Denchfield ++cg->pure_gd_steps; 509484c7b14SAdam Denchfield } 5109566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, 0.0, tao->gradient)); 51150b47da0SAdam Denchfield } else { 5129566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 5139566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work)); 51450b47da0SAdam Denchfield } 51550b47da0SAdam Denchfield break; 516484c7b14SAdam Denchfield 517*d6e07cdcSHong Zhang case TAO_BNCG_HS: 51850b47da0SAdam Denchfield /* Classic Hestenes-Stiefel method, modified with scalar and diagonal preconditioning. */ 51950b47da0SAdam Denchfield if (!cg->diag_scaling) { 52050b47da0SAdam Denchfield cg->sts = step * step * dnorm * dnorm; 5219566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk)); 5229566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, cg->sts, &tau_k, cg->alpha)); 52350b47da0SAdam Denchfield beta = tau_k * gkp1_yk / dk_yk; 5249566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient)); 52550b47da0SAdam Denchfield } else { 5269566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 5279566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->g_work, &gkp1_yk)); 52850b47da0SAdam Denchfield beta = gkp1_yk / dk_yk; 5299566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work)); 53050b47da0SAdam Denchfield } 531c8bcdf1eSAdam Denchfield break; 532484c7b14SAdam Denchfield 533*d6e07cdcSHong Zhang case TAO_BNCG_FR: 5349566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old)); 5359566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient)); 5369566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->yk, NORM_2, &ynorm)); 53750b47da0SAdam Denchfield ynorm2 = ynorm * ynorm; 5389566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->stepdirection, &dk_yk)); 53950b47da0SAdam Denchfield if (!cg->diag_scaling) { 5409566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, step * step * dnorm * dnorm, &tau_k, cg->alpha)); 54150b47da0SAdam Denchfield beta = tau_k * gnorm2 / gnorm2_old; 5429566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient)); 54350b47da0SAdam Denchfield } else { 5449566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old)); /* Before it's updated */ 5459566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 5469566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, cg->g_work, &tmp)); 54750b47da0SAdam Denchfield beta = tmp / gnorm2_old; 5489566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work)); 54950b47da0SAdam Denchfield } 550c8bcdf1eSAdam Denchfield break; 551484c7b14SAdam Denchfield 552*d6e07cdcSHong Zhang case TAO_BNCG_PRP: 55350b47da0SAdam Denchfield snorm = step * dnorm; 55450b47da0SAdam Denchfield if (!cg->diag_scaling) { 5559566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old)); 5569566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk)); 5579566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, snorm * snorm, &tau_k, cg->alpha)); 55850b47da0SAdam Denchfield beta = tau_k * gkp1_yk / gnorm2_old; 5599566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient)); 56050b47da0SAdam Denchfield } else { 5619566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old)); 5629566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 5639566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk)); 56450b47da0SAdam Denchfield beta = gkp1_yk / gnorm2_old; 5659566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work)); 56650b47da0SAdam Denchfield } 567c8bcdf1eSAdam Denchfield break; 568484c7b14SAdam Denchfield 569*d6e07cdcSHong Zhang case TAO_BNCG_PRP_PLUS: 5709566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient)); 5719566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->yk, NORM_2, &ynorm)); 57250b47da0SAdam Denchfield ynorm2 = ynorm * ynorm; 57350b47da0SAdam Denchfield if (!cg->diag_scaling) { 5749566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->G_old, &gnorm2_old)); 5759566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk)); 5769566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, snorm * snorm, &tau_k, cg->alpha)); 57750b47da0SAdam Denchfield beta = tau_k * gkp1_yk / gnorm2_old; 57850b47da0SAdam Denchfield beta = PetscMax(beta, 0.0); 5799566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient)); 58050b47da0SAdam Denchfield } else { 5819566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->g_work, &gnorm2_old)); /* Old gtDg */ 5829566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 5839566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk)); 58450b47da0SAdam Denchfield beta = gkp1_yk / gnorm2_old; 58550b47da0SAdam Denchfield beta = PetscMax(beta, 0.0); 5869566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work)); 58750b47da0SAdam Denchfield } 588c8bcdf1eSAdam Denchfield break; 589484c7b14SAdam Denchfield 590*d6e07cdcSHong Zhang case TAO_BNCG_DY: 591484c7b14SAdam Denchfield /* Dai, Yu-Hong, and Yaxiang Yuan. "A nonlinear conjugate gradient method with a strong global convergence property." 592484c7b14SAdam Denchfield SIAM Journal on optimization 10, no. 1 (1999): 177-182. */ 59350b47da0SAdam Denchfield if (!cg->diag_scaling) { 5949566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, tao->gradient, &gd)); 5959566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old)); 5969566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, step * dk_yk, cg->yts, &tau_k, cg->alpha)); 59750b47da0SAdam Denchfield beta = tau_k * gnorm2 / (gd - gd_old); 5989566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient)); 59950b47da0SAdam Denchfield } else { 6009566063dSJacob Faibussowitsch PetscCall(MatMult(cg->B, tao->stepdirection, cg->d_work)); 6019566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 6029566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, tao->gradient, >Dg)); 6039566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, cg->G_old, &gd_old)); 6049566063dSJacob Faibussowitsch PetscCall(VecDot(cg->d_work, cg->g_work, &dk_yk)); 60550b47da0SAdam Denchfield dk_yk = dk_yk - gd_old; 60650b47da0SAdam Denchfield beta = gtDg / dk_yk; 6079566063dSJacob Faibussowitsch PetscCall(VecScale(cg->d_work, beta)); 6089566063dSJacob Faibussowitsch PetscCall(VecWAXPY(tao->stepdirection, -1.0, cg->g_work, cg->d_work)); 60950b47da0SAdam Denchfield } 610c8bcdf1eSAdam Denchfield break; 611484c7b14SAdam Denchfield 612*d6e07cdcSHong Zhang case TAO_BNCG_HZ: 613484c7b14SAdam Denchfield /* Hager, William W., and Hongchao Zhang. "Algorithm 851: CG_DESCENT, a conjugate gradient method with guaranteed descent." 614484c7b14SAdam Denchfield ACM Transactions on Mathematical Software (TOMS) 32, no. 1 (2006): 113-137. */ 6159566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd)); 6169566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old)); 6179566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution)); 61850b47da0SAdam Denchfield snorm = dnorm * step; 61950b47da0SAdam Denchfield cg->yts = step * dk_yk; 62048a46eb9SPierre Jolivet if (cg->use_dynamic_restart) PetscCall(TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold)); 62150b47da0SAdam Denchfield if (cg->dynamic_restart) { 6229566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2)); 623c8bcdf1eSAdam Denchfield } else { 624c8bcdf1eSAdam Denchfield if (!cg->diag_scaling) { 6259566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk)); 6269566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm * snorm, &tau_k, cg->alpha)); 627c8bcdf1eSAdam Denchfield /* Supplying cg->alpha = -1.0 will give the CG_DESCENT 5.3 special case of tau_k = 1.0 */ 628c8bcdf1eSAdam Denchfield tmp = gd / dk_yk; 629c8bcdf1eSAdam Denchfield beta = tau_k * (gkp1_yk / dk_yk - ynorm2 * gd / (dk_yk * dk_yk)); 630c8bcdf1eSAdam Denchfield /* Bound beta as in CG_DESCENT 5.3, as implemented, with the third comparison from DK 2013 */ 63150b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta * tau_k * gd_old / (dnorm * dnorm)), cg->dk_eta * tau_k * gd / (dnorm * dnorm)); 632c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d */ 6339566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient)); 634c8bcdf1eSAdam Denchfield } else { 635c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 636c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 637c8bcdf1eSAdam Denchfield cg->sts = snorm * snorm; 63850b47da0SAdam Denchfield /* Apply the diagonal scaling to all my vectors */ 6399566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 6409566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, cg->yk, cg->y_work)); 6419566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->stepdirection, cg->d_work)); 642c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 6439566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->g_work, &gkp1_yk)); 644c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */ 645c8bcdf1eSAdam Denchfield tmp = gd / dk_yk; 6469566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->y_work, &tau_k)); 647c8bcdf1eSAdam Denchfield tau_k = -tau_k * gd / (dk_yk * dk_yk); 648c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */ 649484c7b14SAdam Denchfield beta = gkp1_yk / dk_yk + cg->hz_theta * tau_k; /* HZ; (1.15) from DK 2013 */ 650c8bcdf1eSAdam Denchfield /* From HZ2013, modified to account for diagonal scaling*/ 6519566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->d_work, &gd_old)); 6529566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, cg->g_work, &gd)); 65350b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta * gd_old / (dnorm * dnorm)), cg->dk_eta * gd / (dnorm * dnorm)); 654c8bcdf1eSAdam Denchfield /* Do the update */ 6559566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work)); 65650b47da0SAdam Denchfield } 65750b47da0SAdam Denchfield } 658c8bcdf1eSAdam Denchfield break; 659484c7b14SAdam Denchfield 660*d6e07cdcSHong Zhang case TAO_BNCG_DK: 661484c7b14SAdam Denchfield /* Dai, Yu-Hong, and Cai-Xia Kou. "A nonlinear conjugate gradient algorithm with an optimal property and an improved Wolfe line search." 662484c7b14SAdam Denchfield SIAM Journal on Optimization 23, no. 1 (2013): 296-320. */ 6639566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd)); 6649566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, tao->stepdirection, &gd_old)); 6659566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution)); 66650b47da0SAdam Denchfield snorm = step * dnorm; 66750b47da0SAdam Denchfield cg->yts = dk_yk * step; 668c8bcdf1eSAdam Denchfield if (!cg->diag_scaling) { 6699566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk)); 6709566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm * snorm, &tau_k, cg->alpha)); 671c8bcdf1eSAdam Denchfield /* Use cg->alpha = -1.0 to get tau_k = 1.0 as in CG_DESCENT 5.3 */ 672c8bcdf1eSAdam Denchfield tmp = gd / dk_yk; 67350b47da0SAdam Denchfield beta = tau_k * (gkp1_yk / dk_yk - ynorm2 * gd / (dk_yk * dk_yk) + gd / (dnorm * dnorm)) - step * gd / dk_yk; 67450b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta * tau_k * gd_old / (dnorm * dnorm)), cg->dk_eta * tau_k * gd / (dnorm * dnorm)); 675c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d */ 6769566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, 0.0, beta, tao->gradient, cg->yk)); 677c8bcdf1eSAdam Denchfield } else { 678c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 679c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 680c8bcdf1eSAdam Denchfield cg->sts = snorm * snorm; 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)); 6869566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->y_work, &tau_k)); 687c8bcdf1eSAdam Denchfield tau_k = tau_k * gd / (dk_yk * dk_yk); 688c8bcdf1eSAdam Denchfield tmp = gd / dk_yk; 689c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */ 690484c7b14SAdam Denchfield beta = gkp1_yk / dk_yk - step * tmp - tau_k; 691c8bcdf1eSAdam Denchfield /* Update this for the last term in beta */ 6929566063dSJacob Faibussowitsch PetscCall(VecDot(cg->y_work, tao->stepdirection, &dk_yk)); 693c8bcdf1eSAdam Denchfield beta += tmp * dk_yk / (dnorm * dnorm); /* projection of y_work onto dk */ 6949566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, cg->g_work, &gd)); 6959566063dSJacob Faibussowitsch PetscCall(VecDot(cg->G_old, cg->d_work, &gd_old)); 69650b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta * gd_old / (dnorm * dnorm)), cg->dk_eta * gd / (dnorm * dnorm)); 697c8bcdf1eSAdam Denchfield /* Do the update */ 6989566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work)); 69950b47da0SAdam Denchfield } 700c8bcdf1eSAdam Denchfield break; 701484c7b14SAdam Denchfield 702*d6e07cdcSHong Zhang case TAO_BNCG_KD: 703110fc3b0SBarry Smith /* Kou, Cai-Xia, and Yu-Hong Dai. "A modified self-scaling memoryless Broyden-Fletcher-Goldfarb-Shanno method for unconstrained optimization." 704484c7b14SAdam Denchfield Journal of Optimization Theory and Applications 165, no. 1 (2015): 209-224. */ 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; 71048a46eb9SPierre Jolivet if (cg->use_dynamic_restart) PetscCall(TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold)); 71150b47da0SAdam Denchfield if (cg->dynamic_restart) { 7129566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2)); 713c8bcdf1eSAdam Denchfield } else { 714c8bcdf1eSAdam Denchfield if (!cg->diag_scaling) { 7159566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk)); 7169566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm * snorm, &tau_k, cg->alpha)); 717c8bcdf1eSAdam Denchfield beta = tau_k * (gkp1_yk / dk_yk - ynorm2 * gd / (dk_yk * dk_yk)) - step * gd / dk_yk; 718c8bcdf1eSAdam Denchfield if (beta < cg->zeta * tau_k * gd / (dnorm * dnorm)) /* 0.1 is KD's zeta parameter */ 719c8bcdf1eSAdam Denchfield { 720c8bcdf1eSAdam Denchfield beta = cg->zeta * tau_k * gd / (dnorm * dnorm); 721c8bcdf1eSAdam Denchfield gamma = 0.0; 722c8bcdf1eSAdam Denchfield } else { 723c8bcdf1eSAdam Denchfield if (gkp1_yk < 0 && cg->neg_xi) gamma = -1.0 * gd / dk_yk; 724484c7b14SAdam Denchfield /* This seems to be very effective when there's no tau_k scaling. 725484c7b14SAdam Denchfield This guarantees a large descent step every iteration, going through DK 2015 Lemma 3.1's proof but allowing for negative xi */ 726ad540459SPierre Jolivet else gamma = cg->xi * gd / dk_yk; 727c8bcdf1eSAdam Denchfield } 728c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d + t*tmp*yk */ 7299566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, gamma * tau_k, beta, tao->gradient, cg->yk)); 730c8bcdf1eSAdam Denchfield } else { 731c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 732c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 733c8bcdf1eSAdam Denchfield cg->sts = snorm * snorm; 7349566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 7359566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, cg->yk, cg->y_work)); 736c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 7379566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->g_work, &gkp1D_yk)); 738c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */ 739c8bcdf1eSAdam Denchfield gamma = gd / dk_yk; 740c8bcdf1eSAdam Denchfield /* tau_k = -ytDy/(ytd)^2 * gd */ 7419566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, cg->y_work, &tau_k)); 742c8bcdf1eSAdam Denchfield tau_k = tau_k * gd / (dk_yk * dk_yk); 743c8bcdf1eSAdam Denchfield /* beta is the constant which adds the d_k contribution */ 744c8bcdf1eSAdam Denchfield beta = gkp1D_yk / dk_yk - step * gamma - tau_k; 745c8bcdf1eSAdam Denchfield /* Here is the requisite check */ 7469566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, cg->g_work, &tmp)); 747c8bcdf1eSAdam Denchfield if (cg->neg_xi) { 748c8bcdf1eSAdam Denchfield /* modified KD implementation */ 749c8bcdf1eSAdam Denchfield if (gkp1D_yk / dk_yk < 0) gamma = -1.0 * gd / dk_yk; 750ad540459SPierre Jolivet else gamma = cg->xi * gd / dk_yk; 751c8bcdf1eSAdam Denchfield if (beta < cg->zeta * tmp / (dnorm * dnorm)) { 752c8bcdf1eSAdam Denchfield beta = cg->zeta * tmp / (dnorm * dnorm); 753c8bcdf1eSAdam Denchfield gamma = 0.0; 754c8bcdf1eSAdam Denchfield } 755c8bcdf1eSAdam Denchfield } else { /* original KD 2015 implementation */ 756c8bcdf1eSAdam Denchfield if (beta < cg->zeta * tmp / (dnorm * dnorm)) { 757c8bcdf1eSAdam Denchfield beta = cg->zeta * tmp / (dnorm * dnorm); 758c8bcdf1eSAdam Denchfield gamma = 0.0; 759ad540459SPierre Jolivet } else gamma = cg->xi * gd / dk_yk; 760c8bcdf1eSAdam Denchfield } 761c8bcdf1eSAdam Denchfield /* Do the update in two steps */ 7629566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work)); 7639566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->stepdirection, gamma, cg->y_work)); 76450b47da0SAdam Denchfield } 76550b47da0SAdam Denchfield } 766c8bcdf1eSAdam Denchfield break; 767484c7b14SAdam Denchfield 768*d6e07cdcSHong Zhang case TAO_BNCG_SSML_BFGS: 769484c7b14SAdam Denchfield /* Perry, J. M. "A class of conjugate gradient algorithms with a two-step variable-metric memory." 770484c7b14SAdam Denchfield Discussion Papers 269 (1977). */ 7719566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd)); 7729566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution)); 773484c7b14SAdam Denchfield snorm = step * dnorm; 774484c7b14SAdam Denchfield cg->yts = dk_yk * step; 775484c7b14SAdam Denchfield cg->yty = ynorm2; 776484c7b14SAdam Denchfield cg->sts = snorm * snorm; 777484c7b14SAdam Denchfield if (!cg->diag_scaling) { 7789566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk)); 7799566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha)); 780484c7b14SAdam Denchfield tmp = gd / dk_yk; 781484c7b14SAdam Denchfield beta = tau_k * (gkp1_yk / dk_yk - cg->yty * gd / (dk_yk * dk_yk)) - step * tmp; 782484c7b14SAdam Denchfield /* d <- -t*g + beta*t*d + t*tmp*yk */ 7839566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp * tau_k, beta, tao->gradient, cg->yk)); 784484c7b14SAdam Denchfield } else { 785484c7b14SAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 7869566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 7879566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, cg->yk, cg->y_work)); 788484c7b14SAdam Denchfield /* compute scalar gamma */ 7899566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk)); 7909566063dSJacob Faibussowitsch PetscCall(VecDot(cg->y_work, cg->yk, &tmp)); 791484c7b14SAdam Denchfield gamma = gd / dk_yk; 792484c7b14SAdam Denchfield /* Compute scalar beta */ 793484c7b14SAdam Denchfield beta = (gkp1_yk / dk_yk - gd * tmp / (dk_yk * dk_yk)) - step * gd / dk_yk; 794484c7b14SAdam Denchfield /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 7959566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work)); 796484c7b14SAdam Denchfield } 797484c7b14SAdam Denchfield break; 798484c7b14SAdam Denchfield 799*d6e07cdcSHong Zhang case TAO_BNCG_SSML_DFP: 8009566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd)); 8019566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution)); 802484c7b14SAdam Denchfield snorm = step * dnorm; 803484c7b14SAdam Denchfield cg->yts = dk_yk * step; 804484c7b14SAdam Denchfield cg->yty = ynorm2; 805484c7b14SAdam Denchfield cg->sts = snorm * snorm; 806484c7b14SAdam Denchfield if (!cg->diag_scaling) { 807484c7b14SAdam Denchfield /* Instead of a regular convex combination, we will solve a quadratic formula. */ 8089566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha)); 8099566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk)); 810484c7b14SAdam Denchfield tau_k = cg->dfp_scale * tau_k; 811484c7b14SAdam Denchfield tmp = tau_k * gkp1_yk / cg->yty; 812484c7b14SAdam Denchfield beta = -step * gd / dk_yk; 813484c7b14SAdam Denchfield /* d <- -t*g + beta*d + tmp*yk */ 8149566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk)); 815484c7b14SAdam Denchfield } else { 816484c7b14SAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless DFP step */ 8179566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 8189566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, cg->yk, cg->y_work)); 819484c7b14SAdam Denchfield /* compute scalar gamma */ 8209566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk)); 8219566063dSJacob Faibussowitsch PetscCall(VecDot(cg->y_work, cg->yk, &tmp)); 822484c7b14SAdam Denchfield gamma = (gkp1_yk / tmp); 823484c7b14SAdam Denchfield /* Compute scalar beta */ 824484c7b14SAdam Denchfield beta = -step * gd / dk_yk; 825484c7b14SAdam Denchfield /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 8269566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work)); 827484c7b14SAdam Denchfield } 828484c7b14SAdam Denchfield break; 829484c7b14SAdam Denchfield 830*d6e07cdcSHong Zhang case TAO_BNCG_SSML_BRDN: 8319566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd)); 8329566063dSJacob Faibussowitsch PetscCall(VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution)); 833484c7b14SAdam Denchfield snorm = step * dnorm; 834484c7b14SAdam Denchfield cg->yts = step * dk_yk; 835484c7b14SAdam Denchfield cg->yty = ynorm2; 836484c7b14SAdam Denchfield cg->sts = snorm * snorm; 837484c7b14SAdam Denchfield if (!cg->diag_scaling) { 838484c7b14SAdam Denchfield /* Instead of a regular convex combination, we will solve a quadratic formula. */ 8399566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(cg->yty, step * dk_yk, snorm * snorm, &tau_bfgs, cg->bfgs_scale)); 8409566063dSJacob Faibussowitsch PetscCall(TaoBNCGComputeScalarScaling(cg->yty, step * dk_yk, snorm * snorm, &tau_dfp, cg->dfp_scale)); 8419566063dSJacob Faibussowitsch PetscCall(VecDot(cg->yk, tao->gradient, &gkp1_yk)); 842484c7b14SAdam Denchfield tau_k = cg->theta * tau_bfgs + (1.0 - cg->theta) * tau_dfp; 843484c7b14SAdam Denchfield /* If bfgs_scale = 1.0, it should reproduce the bfgs tau_bfgs. If bfgs_scale = 0.0, 844484c7b14SAdam Denchfield it should reproduce the tau_dfp scaling. Same with dfp_scale. */ 845484c7b14SAdam Denchfield tmp = cg->theta * tau_bfgs * gd / dk_yk + (1 - cg->theta) * tau_dfp * gkp1_yk / cg->yty; 846484c7b14SAdam Denchfield beta = cg->theta * tau_bfgs * (gkp1_yk / dk_yk - cg->yty * gd / (dk_yk * dk_yk)) - step * gd / dk_yk; 847484c7b14SAdam Denchfield /* d <- -t*g + beta*d + tmp*yk */ 8489566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk)); 849484c7b14SAdam Denchfield } else { 850484c7b14SAdam Denchfield /* We have diagonal scaling enabled */ 8519566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, tao->gradient, cg->g_work)); 8529566063dSJacob Faibussowitsch PetscCall(MatSolve(cg->B, cg->yk, cg->y_work)); 853484c7b14SAdam Denchfield /* compute scalar gamma */ 8549566063dSJacob Faibussowitsch PetscCall(VecDot(cg->g_work, cg->yk, &gkp1_yk)); 8559566063dSJacob Faibussowitsch PetscCall(VecDot(cg->y_work, cg->yk, &tmp)); 856484c7b14SAdam Denchfield gamma = cg->theta * gd / dk_yk + (1 - cg->theta) * (gkp1_yk / tmp); 857484c7b14SAdam Denchfield /* Compute scalar beta */ 858484c7b14SAdam Denchfield beta = cg->theta * (gkp1_yk / dk_yk - gd * tmp / (dk_yk * dk_yk)) - step * gd / dk_yk; 859484c7b14SAdam Denchfield /* Compute stepdirection dkp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 8609566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work)); 861484c7b14SAdam Denchfield } 862484c7b14SAdam Denchfield break; 863484c7b14SAdam Denchfield 864d71ae5a4SJacob Faibussowitsch default: 865d71ae5a4SJacob Faibussowitsch break; 866c8bcdf1eSAdam Denchfield } 867c8bcdf1eSAdam Denchfield } 8683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 869c8bcdf1eSAdam Denchfield } 870c8bcdf1eSAdam Denchfield 871d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao tao, PetscReal gnorm) 872d71ae5a4SJacob Faibussowitsch { 873c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG *)tao->data; 874c8bcdf1eSAdam Denchfield TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 8758ca2df50S PetscReal step = 1.0, gnorm2, gd, dnorm = 0.0; 876c8bcdf1eSAdam Denchfield PetscReal gnorm2_old, f_old, resnorm, gnorm_old; 877c624ebd3SAlp Dener PetscBool pcgd_fallback = PETSC_FALSE; 878c8bcdf1eSAdam Denchfield 879c8bcdf1eSAdam Denchfield PetscFunctionBegin; 880c8bcdf1eSAdam Denchfield /* We are now going to perform a line search along the direction. */ 881c8bcdf1eSAdam Denchfield /* Store solution and gradient info before it changes */ 8829566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, cg->X_old)); 8839566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, cg->G_old)); 8849566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->unprojected_gradient, cg->unprojected_gradient_old)); 885c8bcdf1eSAdam Denchfield 886c8bcdf1eSAdam Denchfield gnorm_old = gnorm; 887c8bcdf1eSAdam Denchfield gnorm2_old = gnorm_old * gnorm_old; 888c8bcdf1eSAdam Denchfield f_old = cg->f; 889484c7b14SAdam Denchfield /* Perform bounded line search. If we are recycling a solution from a previous */ 890484c7b14SAdam Denchfield /* TaoSolve, then we want to immediately skip to calculating a new direction rather than performing a linesearch */ 891414d97d3SAlp Dener if (!(tao->recycle && 0 == tao->niter)) { 892484c7b14SAdam Denchfield /* Above logic: the below code happens every iteration, except for the first iteration of a recycled TaoSolve */ 8939566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 8949566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status)); 8959566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 896c8bcdf1eSAdam Denchfield 897c8bcdf1eSAdam Denchfield /* Check linesearch failure */ 898c8bcdf1eSAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 899c8bcdf1eSAdam Denchfield ++cg->ls_fails; 900*d6e07cdcSHong Zhang if (cg->cg_type == TAO_BNCG_GD) { 901c8bcdf1eSAdam Denchfield /* Nothing left to do but fail out of the optimization */ 902c8bcdf1eSAdam Denchfield step = 0.0; 903c8bcdf1eSAdam Denchfield tao->reason = TAO_DIVERGED_LS_FAILURE; 904c8bcdf1eSAdam Denchfield } else { 905484c7b14SAdam Denchfield /* Restore previous point, perform preconditioned GD and regular GD steps at the last good point */ 9069566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->X_old, tao->solution)); 9079566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->G_old, tao->gradient)); 9089566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->unprojected_gradient_old, cg->unprojected_gradient)); 909c8bcdf1eSAdam Denchfield gnorm = gnorm_old; 910c8bcdf1eSAdam Denchfield gnorm2 = gnorm2_old; 911c8bcdf1eSAdam Denchfield cg->f = f_old; 912c8bcdf1eSAdam Denchfield 913484c7b14SAdam Denchfield /* Fall back on preconditioned CG (so long as you're not already using it) */ 914*d6e07cdcSHong Zhang if (cg->cg_type != TAO_BNCG_PCGD && cg->diag_scaling) { 915e2570530SAlp Dener pcgd_fallback = PETSC_TRUE; 9169566063dSJacob Faibussowitsch PetscCall(TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback)); 917484c7b14SAdam Denchfield 9189566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2)); 9199566063dSJacob Faibussowitsch PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection)); 920c8bcdf1eSAdam Denchfield 9219566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 9229566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status)); 9239566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 924c8bcdf1eSAdam Denchfield 925484c7b14SAdam Denchfield pcgd_fallback = PETSC_FALSE; 926484c7b14SAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 927484c7b14SAdam Denchfield /* Going to perform a regular gradient descent step. */ 928484c7b14SAdam Denchfield ++cg->ls_fails; 929484c7b14SAdam Denchfield step = 0.0; 930484c7b14SAdam Denchfield } 931484c7b14SAdam Denchfield } 932484c7b14SAdam Denchfield /* Fall back on the scaled gradient step */ 933484c7b14SAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 934484c7b14SAdam Denchfield ++cg->ls_fails; 9359566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2)); 9369566063dSJacob Faibussowitsch PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection)); 9379566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 9389566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status)); 9399566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 940484c7b14SAdam Denchfield } 941484c7b14SAdam Denchfield 942c8bcdf1eSAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 943c8bcdf1eSAdam Denchfield /* Nothing left to do but fail out of the optimization */ 94450b47da0SAdam Denchfield ++cg->ls_fails; 945c8bcdf1eSAdam Denchfield step = 0.0; 946c8bcdf1eSAdam Denchfield tao->reason = TAO_DIVERGED_LS_FAILURE; 947484c7b14SAdam Denchfield } else { 948484c7b14SAdam Denchfield /* One of the fallbacks worked. Set them both back equal to false. */ 949484c7b14SAdam Denchfield pcgd_fallback = PETSC_FALSE; 950c8bcdf1eSAdam Denchfield } 951c8bcdf1eSAdam Denchfield } 952c8bcdf1eSAdam Denchfield } 953c8bcdf1eSAdam Denchfield /* Convergence test for line search failure */ 9543ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 955c8bcdf1eSAdam Denchfield 956c8bcdf1eSAdam Denchfield /* Standard convergence test */ 9579566063dSJacob Faibussowitsch PetscCall(VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W)); 9589566063dSJacob Faibussowitsch PetscCall(VecNorm(cg->W, NORM_2, &resnorm)); 9593c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 9609566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its)); 9619566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step)); 962dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 9633ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 964484c7b14SAdam Denchfield } 965c8bcdf1eSAdam Denchfield /* Assert we have an updated step and we need at least one more iteration. */ 966c8bcdf1eSAdam Denchfield /* Calculate the next direction */ 967c8bcdf1eSAdam Denchfield /* Estimate the active set at the new solution */ 9689566063dSJacob Faibussowitsch PetscCall(TaoBNCGEstimateActiveSet(tao, cg->as_type)); 969c8bcdf1eSAdam Denchfield /* Compute the projected gradient and its norm */ 9709566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->unprojected_gradient, tao->gradient)); 9719566063dSJacob Faibussowitsch PetscCall(VecISSet(tao->gradient, cg->active_idx, 0.0)); 9729566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm)); 973c8bcdf1eSAdam Denchfield gnorm2 = gnorm * gnorm; 974c8bcdf1eSAdam Denchfield 975484c7b14SAdam Denchfield /* Calculate some quantities used in the StepDirectionUpdate. */ 9769566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection, NORM_2, &dnorm)); 977484c7b14SAdam Denchfield /* Update the step direction. */ 9789566063dSJacob Faibussowitsch PetscCall(TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback)); 979484c7b14SAdam Denchfield ++tao->niter; 9809566063dSJacob Faibussowitsch PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection)); 981c8bcdf1eSAdam Denchfield 982*d6e07cdcSHong Zhang if (cg->cg_type != TAO_BNCG_GD) { 983c8bcdf1eSAdam Denchfield /* Figure out which previously active variables became inactive this iteration */ 9849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cg->new_inactives)); 98548a46eb9SPierre Jolivet if (cg->inactive_idx && cg->inactive_old) PetscCall(ISDifference(cg->inactive_idx, cg->inactive_old, &cg->new_inactives)); 986c8bcdf1eSAdam Denchfield /* Selectively reset the CG step those freshly inactive variables */ 987c8bcdf1eSAdam Denchfield if (cg->new_inactives) { 9889566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step)); 9899566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad)); 9909566063dSJacob Faibussowitsch PetscCall(VecCopy(cg->inactive_grad, cg->inactive_step)); 9919566063dSJacob Faibussowitsch PetscCall(VecScale(cg->inactive_step, -1.0)); 9929566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step)); 9939566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad)); 994c8bcdf1eSAdam Denchfield } 995c8bcdf1eSAdam Denchfield /* Verify that this is a descent direction */ 9969566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd)); 9979566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection, NORM_2, &dnorm)); 99850b47da0SAdam Denchfield if (PetscIsInfOrNanReal(gd) || (gd / (dnorm * dnorm) <= -1e10 || gd / (dnorm * dnorm) >= -1e-10)) { 999c8bcdf1eSAdam Denchfield /* Not a descent direction, so we reset back to projected gradient descent */ 10009566063dSJacob Faibussowitsch PetscCall(TaoBNCGResetUpdate(tao, gnorm2)); 10019566063dSJacob Faibussowitsch PetscCall(TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection)); 1002c8bcdf1eSAdam Denchfield ++cg->descent_error; 1003c8bcdf1eSAdam Denchfield } else { 1004c8bcdf1eSAdam Denchfield } 1005c8bcdf1eSAdam Denchfield } 10063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1007ac9112b8SAlp Dener } 1008484c7b14SAdam Denchfield 1009*d6e07cdcSHong Zhang PETSC_INTERN PetscErrorCode TaoBNCGSetH0(Tao tao, Mat H0) 1010d71ae5a4SJacob Faibussowitsch { 1011484c7b14SAdam Denchfield TAO_BNCG *cg = (TAO_BNCG *)tao->data; 1012*d6e07cdcSHong Zhang PetscBool same; 1013484c7b14SAdam Denchfield 1014484c7b14SAdam Denchfield PetscFunctionBegin; 1015*d6e07cdcSHong Zhang PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBNCG, &same)); 1016*d6e07cdcSHong Zhang if (same) { 10179566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H0)); 1018484c7b14SAdam Denchfield cg->pc = H0; 1019*d6e07cdcSHong Zhang } 1020*d6e07cdcSHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 1021*d6e07cdcSHong Zhang } 1022*d6e07cdcSHong Zhang 1023*d6e07cdcSHong Zhang /*@ 1024*d6e07cdcSHong Zhang TaoBNCGGetType - Return the type for the `TAOBNCG` solver 1025*d6e07cdcSHong Zhang 1026*d6e07cdcSHong Zhang Input Parameter: 1027*d6e07cdcSHong Zhang . tao - the `Tao` solver context 1028*d6e07cdcSHong Zhang 1029*d6e07cdcSHong Zhang Output Parameter: 1030*d6e07cdcSHong Zhang . type - `TAOBNCG` type 1031*d6e07cdcSHong Zhang 1032*d6e07cdcSHong Zhang Level: advanced 1033*d6e07cdcSHong Zhang 1034*d6e07cdcSHong Zhang .seealso: `Tao`, `TAOBNCG`, `TaoBNCGSetType()`, `TaoBNCGType` 1035*d6e07cdcSHong Zhang @*/ 1036*d6e07cdcSHong Zhang PetscErrorCode TaoBNCGGetType(Tao tao, TaoBNCGType *type) 1037*d6e07cdcSHong Zhang { 1038*d6e07cdcSHong Zhang TAO_BNCG *cg = (TAO_BNCG *)tao->data; 1039*d6e07cdcSHong Zhang PetscBool same; 1040*d6e07cdcSHong Zhang 1041*d6e07cdcSHong Zhang PetscFunctionBegin; 1042*d6e07cdcSHong Zhang PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBNCG, &same)); 1043*d6e07cdcSHong Zhang PetscCheck(same, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "TAO solver is not BNCG type"); 1044*d6e07cdcSHong Zhang *type = cg->cg_type; 1045*d6e07cdcSHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 1046*d6e07cdcSHong Zhang } 1047*d6e07cdcSHong Zhang 1048*d6e07cdcSHong Zhang /*@ 1049*d6e07cdcSHong Zhang TaoBNCGSetType - Set the type for the `TAOBNCG` solver 1050*d6e07cdcSHong Zhang 1051*d6e07cdcSHong Zhang Input Parameters: 1052*d6e07cdcSHong Zhang + tao - the `Tao` solver context 1053*d6e07cdcSHong Zhang - type - `TAOBNCG` type 1054*d6e07cdcSHong Zhang 1055*d6e07cdcSHong Zhang Level: advanced 1056*d6e07cdcSHong Zhang 1057*d6e07cdcSHong Zhang .seealso: `Tao`, `TAOBNCG`, `TaoBNCGGetType()`, `TaoBNCGType` 1058*d6e07cdcSHong Zhang @*/ 1059*d6e07cdcSHong Zhang PetscErrorCode TaoBNCGSetType(Tao tao, TaoBNCGType type) 1060*d6e07cdcSHong Zhang { 1061*d6e07cdcSHong Zhang TAO_BNCG *cg = (TAO_BNCG *)tao->data; 1062*d6e07cdcSHong Zhang PetscBool same; 1063*d6e07cdcSHong Zhang 1064*d6e07cdcSHong Zhang PetscFunctionBegin; 1065*d6e07cdcSHong Zhang PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBNCG, &same)); 1066*d6e07cdcSHong Zhang if (same) cg->cg_type = type; 10673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1068484c7b14SAdam Denchfield } 1069