1ac9112b8SAlp Dener #include <petsctaolinesearch.h> 2ac9112b8SAlp Dener #include <../src/tao/bound/impls/bncg/bncg.h> 350b47da0SAdam Denchfield #include <petscksp.h> 4ac9112b8SAlp Dener 5c8bcdf1eSAdam Denchfield #define CG_GradientDescent 0 6c8bcdf1eSAdam Denchfield #define CG_HestenesStiefel 1 7c8bcdf1eSAdam Denchfield #define CG_FletcherReeves 2 850b47da0SAdam Denchfield #define CG_PolakRibierePolyak 3 9c8bcdf1eSAdam Denchfield #define CG_PolakRibierePlus 4 10c8bcdf1eSAdam Denchfield #define CG_DaiYuan 5 11c8bcdf1eSAdam Denchfield #define CG_HagerZhang 6 12c8bcdf1eSAdam Denchfield #define CG_DaiKou 7 13c8bcdf1eSAdam Denchfield #define CG_KouDai 8 14c8bcdf1eSAdam Denchfield #define CG_SSML_BFGS 9 15c8bcdf1eSAdam Denchfield #define CG_SSML_DFP 10 16c8bcdf1eSAdam Denchfield #define CG_SSML_BROYDEN 11 17484c7b14SAdam Denchfield #define CG_PCGradientDescent 12 18484c7b14SAdam Denchfield #define CGTypes 13 19ac9112b8SAlp Dener 20484c7b14SAdam Denchfield static const char *CG_Table[64] = {"gd", "hs", "fr", "pr", "prp", "dy", "hz", "dk", "kd", "ssml_bfgs", "ssml_dfp", "ssml_brdn", "pcgd"}; 21ac9112b8SAlp Dener 2261be54a6SAlp Dener #define CG_AS_NONE 0 2361be54a6SAlp Dener #define CG_AS_BERTSEKAS 1 2461be54a6SAlp Dener #define CG_AS_SIZE 2 25ac9112b8SAlp Dener 2661be54a6SAlp Dener static const char *CG_AS_TYPE[64] = {"none", "bertsekas"}; 27ac9112b8SAlp Dener 28c0f10754SAlp Dener PetscErrorCode TaoBNCGSetRecycleFlag(Tao tao, PetscBool recycle) 29c0f10754SAlp Dener { 30c0f10754SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 31c0f10754SAlp Dener 32c0f10754SAlp Dener PetscFunctionBegin; 33c0f10754SAlp Dener cg->recycle = recycle; 34c0f10754SAlp Dener PetscFunctionReturn(0); 35c0f10754SAlp Dener } 36c0f10754SAlp Dener 3761be54a6SAlp Dener PetscErrorCode TaoBNCGEstimateActiveSet(Tao tao, PetscInt asType) 3861be54a6SAlp Dener { 3961be54a6SAlp Dener PetscErrorCode ierr; 4061be54a6SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 4161be54a6SAlp Dener 4261be54a6SAlp Dener PetscFunctionBegin; 4361be54a6SAlp Dener ierr = ISDestroy(&cg->inactive_old);CHKERRQ(ierr); 4461be54a6SAlp Dener if (cg->inactive_idx) { 4561be54a6SAlp Dener ierr = ISDuplicate(cg->inactive_idx, &cg->inactive_old);CHKERRQ(ierr); 4661be54a6SAlp Dener ierr = ISCopy(cg->inactive_idx, cg->inactive_old);CHKERRQ(ierr); 4761be54a6SAlp Dener } 4861be54a6SAlp Dener switch (asType) { 4961be54a6SAlp Dener case CG_AS_NONE: 5061be54a6SAlp Dener ierr = ISDestroy(&cg->inactive_idx);CHKERRQ(ierr); 5161be54a6SAlp Dener ierr = VecWhichInactive(tao->XL, tao->solution, cg->unprojected_gradient, tao->XU, PETSC_TRUE, &cg->inactive_idx);CHKERRQ(ierr); 5261be54a6SAlp Dener ierr = ISDestroy(&cg->active_idx);CHKERRQ(ierr); 5361be54a6SAlp Dener ierr = ISComplementVec(cg->inactive_idx, tao->solution, &cg->active_idx);CHKERRQ(ierr); 5461be54a6SAlp Dener break; 5561be54a6SAlp Dener 5661be54a6SAlp Dener case CG_AS_BERTSEKAS: 5761be54a6SAlp Dener /* Use gradient descent to estimate the active set */ 5861be54a6SAlp Dener ierr = VecCopy(cg->unprojected_gradient, cg->W);CHKERRQ(ierr); 5961be54a6SAlp Dener ierr = VecScale(cg->W, -1.0);CHKERRQ(ierr); 6089da521bSAlp Dener ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, cg->unprojected_gradient, cg->W, cg->work, cg->as_step, &cg->as_tol, 6189da521bSAlp Dener &cg->active_lower, &cg->active_upper, &cg->active_fixed, &cg->active_idx, &cg->inactive_idx);CHKERRQ(ierr); 62c4b75bccSAlp Dener break; 6361be54a6SAlp Dener 6461be54a6SAlp Dener default: 6561be54a6SAlp Dener break; 6661be54a6SAlp Dener } 6761be54a6SAlp Dener PetscFunctionReturn(0); 6861be54a6SAlp Dener } 6961be54a6SAlp Dener 70a1318120SAlp Dener PetscErrorCode TaoBNCGBoundStep(Tao tao, PetscInt asType, Vec step) 7161be54a6SAlp Dener { 7261be54a6SAlp Dener PetscErrorCode ierr; 7361be54a6SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 7461be54a6SAlp Dener 7561be54a6SAlp Dener PetscFunctionBegin; 76a1318120SAlp Dener switch (asType) { 7761be54a6SAlp Dener case CG_AS_NONE: 78c4b75bccSAlp Dener ierr = VecISSet(step, cg->active_idx, 0.0);CHKERRQ(ierr); 7961be54a6SAlp Dener break; 8061be54a6SAlp Dener 8161be54a6SAlp Dener case CG_AS_BERTSEKAS: 82c4b75bccSAlp Dener ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, cg->active_lower, cg->active_upper, cg->active_fixed, 1.0, step);CHKERRQ(ierr); 8361be54a6SAlp Dener break; 8461be54a6SAlp Dener 8561be54a6SAlp Dener default: 8661be54a6SAlp Dener break; 8761be54a6SAlp Dener } 8861be54a6SAlp Dener PetscFunctionReturn(0); 8961be54a6SAlp Dener } 9061be54a6SAlp Dener 91ac9112b8SAlp Dener static PetscErrorCode TaoSolve_BNCG(Tao tao) 92ac9112b8SAlp Dener { 93ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 94ac9112b8SAlp Dener PetscErrorCode ierr; 95484c7b14SAdam Denchfield PetscReal step=1.0,gnorm,gnorm2, resnorm; 96c4b75bccSAlp Dener PetscInt nDiff; 97ac9112b8SAlp Dener 98ac9112b8SAlp Dener PetscFunctionBegin; 99ac9112b8SAlp Dener /* Project the current point onto the feasible set */ 100ac9112b8SAlp Dener ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 101cd929ea3SAlp Dener ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 102ac9112b8SAlp Dener 103c8bcdf1eSAdam Denchfield /* Project the initial point onto the feasible region */ 104c8bcdf1eSAdam Denchfield ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr); 105484c7b14SAdam Denchfield 106484c7b14SAdam Denchfield if (nDiff > 0 || !cg->recycle){ 107c0f10754SAlp Dener ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient);CHKERRQ(ierr); 108484c7b14SAdam Denchfield } 109ac9112b8SAlp Dener ierr = VecNorm(cg->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr); 110c0f10754SAlp Dener if (PetscIsInfOrNanReal(cg->f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 111ac9112b8SAlp Dener 11261be54a6SAlp Dener /* Estimate the active set and compute the projected gradient */ 11361be54a6SAlp Dener ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr); 11461be54a6SAlp Dener 115ac9112b8SAlp Dener /* Project the gradient and calculate the norm */ 11661be54a6SAlp Dener ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 11761be54a6SAlp Dener ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr); 118ac9112b8SAlp Dener ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 119ac9112b8SAlp Dener gnorm2 = gnorm*gnorm; 120ac9112b8SAlp Dener 121c8bcdf1eSAdam Denchfield /* Initialize counters */ 122e031d6f5SAlp Dener tao->niter = 0; 12350b47da0SAdam Denchfield cg->ls_fails = cg->descent_error = 0; 124c8bcdf1eSAdam Denchfield cg->resets = -1; 125484c7b14SAdam Denchfield cg->skipped_updates = cg->pure_gd_steps = 0; 126c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 127c8bcdf1eSAdam Denchfield 128c8bcdf1eSAdam Denchfield /* Convergence test at the starting point. */ 129ac9112b8SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 130484c7b14SAdam Denchfield 131484c7b14SAdam Denchfield ierr = VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W);CHKERRQ(ierr); 132484c7b14SAdam Denchfield ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr); 133484c7b14SAdam Denchfield if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 134484c7b14SAdam Denchfield ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 135484c7b14SAdam Denchfield ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr); 136ac9112b8SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 137ac9112b8SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 138484c7b14SAdam Denchfield /* Calculate initial direction. */ 139484c7b14SAdam Denchfield if (!cg->recycle) { 140484c7b14SAdam Denchfield /* We are not recycling a solution/history from a past TaoSolve */ 141484c7b14SAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 142ac9112b8SAlp Dener } 143c8bcdf1eSAdam Denchfield /* Initial gradient descent step. Scaling by 1.0 also does a decent job for some problems. */ 144c8bcdf1eSAdam Denchfield while(1) { 145c8bcdf1eSAdam Denchfield ierr = TaoBNCGConductIteration(tao, gnorm);CHKERRQ(ierr); 146c8bcdf1eSAdam Denchfield if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 147ac9112b8SAlp Dener } 148ac9112b8SAlp Dener PetscFunctionReturn(0); 149ac9112b8SAlp Dener } 150ac9112b8SAlp Dener 151ac9112b8SAlp Dener static PetscErrorCode TaoSetUp_BNCG(Tao tao) 152ac9112b8SAlp Dener { 153ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 154ac9112b8SAlp Dener PetscErrorCode ierr; 155ac9112b8SAlp Dener 156ac9112b8SAlp Dener PetscFunctionBegin; 157c4b75bccSAlp Dener if (!tao->gradient) { 158c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 159c4b75bccSAlp Dener } 160c4b75bccSAlp Dener if (!tao->stepdirection) { 161c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); 162c4b75bccSAlp Dener } 163c4b75bccSAlp Dener if (!cg->W) { 164c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&cg->W);CHKERRQ(ierr); 165c4b75bccSAlp Dener } 166c4b75bccSAlp Dener if (!cg->work) { 167c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&cg->work);CHKERRQ(ierr); 168c4b75bccSAlp Dener } 169c8bcdf1eSAdam Denchfield if (!cg->sk) { 170c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->sk);CHKERRQ(ierr); 171c8bcdf1eSAdam Denchfield } 172c8bcdf1eSAdam Denchfield if (!cg->yk) { 173c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->gradient,&cg->yk);CHKERRQ(ierr); 174c8bcdf1eSAdam Denchfield } 175c4b75bccSAlp Dener if (!cg->X_old) { 176c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&cg->X_old);CHKERRQ(ierr); 177c4b75bccSAlp Dener } 178c4b75bccSAlp Dener if (!cg->G_old) { 179c4b75bccSAlp Dener ierr = VecDuplicate(tao->gradient,&cg->G_old);CHKERRQ(ierr); 180c8bcdf1eSAdam Denchfield } 181c8bcdf1eSAdam Denchfield if (cg->diag_scaling){ 182c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->d_work);CHKERRQ(ierr); 183c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->y_work);CHKERRQ(ierr); 18450b47da0SAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->g_work);CHKERRQ(ierr); 185c4b75bccSAlp Dener } 186c4b75bccSAlp Dener if (!cg->unprojected_gradient) { 187c4b75bccSAlp Dener ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient);CHKERRQ(ierr); 188c4b75bccSAlp Dener } 189c4b75bccSAlp Dener if (!cg->unprojected_gradient_old) { 190c4b75bccSAlp Dener ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient_old);CHKERRQ(ierr); 191c4b75bccSAlp Dener } 19250b47da0SAdam Denchfield ierr = MatLMVMAllocate(cg->B, cg->sk, cg->yk);CHKERRQ(ierr); 193484c7b14SAdam Denchfield if (cg->pc){ 194484c7b14SAdam Denchfield ierr = MatLMVMSetJ0(cg->B, cg->pc);CHKERRQ(ierr); 195484c7b14SAdam Denchfield } 196ac9112b8SAlp Dener PetscFunctionReturn(0); 197ac9112b8SAlp Dener } 198ac9112b8SAlp Dener 199ac9112b8SAlp Dener static PetscErrorCode TaoDestroy_BNCG(Tao tao) 200ac9112b8SAlp Dener { 201ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*) tao->data; 202ac9112b8SAlp Dener PetscErrorCode ierr; 203ac9112b8SAlp Dener 204ac9112b8SAlp Dener PetscFunctionBegin; 205ac9112b8SAlp Dener if (tao->setupcalled) { 20661be54a6SAlp Dener ierr = VecDestroy(&cg->W);CHKERRQ(ierr); 207c4b75bccSAlp Dener ierr = VecDestroy(&cg->work);CHKERRQ(ierr); 208ac9112b8SAlp Dener ierr = VecDestroy(&cg->X_old);CHKERRQ(ierr); 209ac9112b8SAlp Dener ierr = VecDestroy(&cg->G_old);CHKERRQ(ierr); 210ac9112b8SAlp Dener ierr = VecDestroy(&cg->unprojected_gradient);CHKERRQ(ierr); 211ac9112b8SAlp Dener ierr = VecDestroy(&cg->unprojected_gradient_old);CHKERRQ(ierr); 21250b47da0SAdam Denchfield ierr = VecDestroy(&cg->g_work);CHKERRQ(ierr); 21350b47da0SAdam Denchfield ierr = VecDestroy(&cg->d_work);CHKERRQ(ierr); 21450b47da0SAdam Denchfield ierr = VecDestroy(&cg->y_work);CHKERRQ(ierr); 21550b47da0SAdam Denchfield ierr = VecDestroy(&cg->sk);CHKERRQ(ierr); 21650b47da0SAdam Denchfield ierr = VecDestroy(&cg->yk);CHKERRQ(ierr); 217ac9112b8SAlp Dener } 218ca964c49SAlp Dener ierr = ISDestroy(&cg->active_lower);CHKERRQ(ierr); 219ca964c49SAlp Dener ierr = ISDestroy(&cg->active_upper);CHKERRQ(ierr); 220ca964c49SAlp Dener ierr = ISDestroy(&cg->active_fixed);CHKERRQ(ierr); 221ca964c49SAlp Dener ierr = ISDestroy(&cg->active_idx);CHKERRQ(ierr); 222ca964c49SAlp Dener ierr = ISDestroy(&cg->inactive_idx);CHKERRQ(ierr); 223ca964c49SAlp Dener ierr = ISDestroy(&cg->inactive_old);CHKERRQ(ierr); 224ca964c49SAlp Dener ierr = ISDestroy(&cg->new_inactives);CHKERRQ(ierr); 225*01e1e359SAlp Dener ierr = MatDestroy(&cg->B);CHKERRQ(ierr); 226484c7b14SAdam Denchfield if (cg->pc) { 227*01e1e359SAlp Dener ierr = MatDestroy(&cg->pc);CHKERRQ(ierr); 228484c7b14SAdam Denchfield } 229ac9112b8SAlp Dener ierr = PetscFree(tao->data);CHKERRQ(ierr); 230ac9112b8SAlp Dener PetscFunctionReturn(0); 231ac9112b8SAlp Dener } 232ac9112b8SAlp Dener 233ac9112b8SAlp Dener static PetscErrorCode TaoSetFromOptions_BNCG(PetscOptionItems *PetscOptionsObject,Tao tao) 234ac9112b8SAlp Dener { 235ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 236ac9112b8SAlp Dener PetscErrorCode ierr; 237ac9112b8SAlp Dener 238ac9112b8SAlp Dener PetscFunctionBegin; 239ac9112b8SAlp Dener ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 240ac9112b8SAlp Dener ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Conjugate Gradient method for unconstrained optimization");CHKERRQ(ierr); 241484c7b14SAdam Denchfield ierr = PetscOptionsEList("-tao_bncg_type","cg formula", "", CG_Table, CGTypes, CG_Table[cg->cg_type], &cg->cg_type,NULL);CHKERRQ(ierr); 242484c7b14SAdam Denchfield if (cg->cg_type != CG_SSML_BFGS){ 243484c7b14SAdam Denchfield cg->alpha = -1.0; /* Setting defaults for non-BFGS methods. User can change it below. */ 244484c7b14SAdam Denchfield } 245484c7b14SAdam Denchfield if (CG_GradientDescent == cg->cg_type){ 246484c7b14SAdam Denchfield cg->cg_type = CG_PCGradientDescent; 247484c7b14SAdam Denchfield /* Set scaling equal to none or, at best, scalar scaling. */ 248484c7b14SAdam Denchfield cg->unscaled_restart = PETSC_TRUE; 249484c7b14SAdam Denchfield cg->diag_scaling = PETSC_FALSE; 250484c7b14SAdam Denchfield } 25150b47da0SAdam Denchfield ierr = PetscOptionsEList("-tao_bncg_as_type","active set estimation method", "", CG_AS_TYPE, CG_AS_SIZE, CG_AS_TYPE[cg->cg_type], &cg->cg_type,NULL);CHKERRQ(ierr); 25250b47da0SAdam Denchfield 25350b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_hz_eta","(developer) cutoff tolerance for HZ", "", cg->hz_eta,&cg->hz_eta,NULL);CHKERRQ(ierr); 25450b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_eps","(developer) cutoff value for restarts", "", cg->epsilon,&cg->epsilon,NULL);CHKERRQ(ierr); 25550b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_dk_eta","(developer) cutoff tolerance for DK", "", cg->dk_eta,&cg->dk_eta,NULL);CHKERRQ(ierr); 25650b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_xi","(developer) Parameter in the KD method", "", cg->xi,&cg->xi,NULL);CHKERRQ(ierr); 25750b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_theta", "(developer) update parameter for the Broyden method", "", cg->theta, &cg->theta, NULL);CHKERRQ(ierr); 25850b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_hz_theta", "(developer) parameter for the HZ (2006) method", "", cg->hz_theta, &cg->hz_theta, NULL);CHKERRQ(ierr); 25950b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_alpha","(developer) parameter for the scalar scaling","",cg->alpha,&cg->alpha,NULL);CHKERRQ(ierr); 260c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_bfgs_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->bfgs_scale, &cg->bfgs_scale, NULL);CHKERRQ(ierr); 261c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_dfp_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->dfp_scale, &cg->dfp_scale, NULL);CHKERRQ(ierr); 262c8bcdf1eSAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_diag_scaling","Enable diagonal Broyden-like preconditioning","",cg->diag_scaling,&cg->diag_scaling,NULL);CHKERRQ(ierr); 26350b47da0SAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_dynamic_restart","(developer) use dynamic restarts as in HZ, DK, KD","",cg->use_dynamic_restart,&cg->use_dynamic_restart,NULL);CHKERRQ(ierr); 26450b47da0SAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_unscaled_restart","(developer) use unscaled gradient restarts","",cg->unscaled_restart,&cg->unscaled_restart,NULL);CHKERRQ(ierr); 26550b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_zeta", "(developer) Free parameter for the Kou-Dai method", "", cg->zeta, &cg->zeta, NULL);CHKERRQ(ierr); 266c8bcdf1eSAdam Denchfield ierr = PetscOptionsInt("-tao_bncg_min_quad", "(developer) Number of iterations with approximate quadratic behavior needed for restart", "", cg->min_quad, &cg->min_quad, NULL);CHKERRQ(ierr); 26750b47da0SAdam Denchfield ierr = PetscOptionsInt("-tao_bncg_min_restart_num", "(developer) Number of iterations between restarts (times dimension)", "", cg->min_restart_num, &cg->min_restart_num, NULL);CHKERRQ(ierr); 268484c7b14SAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_recycle","enable recycling the existing solution, gradient, and diagonal scaling vector at the start of a new TaoSolve() call","",cg->recycle,&cg->recycle,NULL);CHKERRQ(ierr); 26950b47da0SAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_spaced_restart","(developer) Enable regular steepest descent restarting every fixed number of iterations","",cg->spaced_restart,&cg->spaced_restart,NULL);CHKERRQ(ierr); 270484c7b14SAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_no_scaling","Disable all scaling except in restarts","",cg->no_scaling,&cg->no_scaling,NULL);CHKERRQ(ierr); 271484c7b14SAdam Denchfield if (cg->no_scaling){ 272484c7b14SAdam Denchfield cg->diag_scaling = PETSC_FALSE; 273484c7b14SAdam Denchfield cg->alpha = -1.0; 274484c7b14SAdam Denchfield } 275b474139fSKarl Rupp if (cg->alpha == -1.0 && cg->cg_type == CG_KouDai && !cg->diag_scaling){ /* Some more default options that appear to be good. */ 276484c7b14SAdam Denchfield cg->neg_xi = PETSC_TRUE; 277484c7b14SAdam Denchfield } 27850b47da0SAdam Denchfield ierr = 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);CHKERRQ(ierr); 27950b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_as_tol", "(developer) initial tolerance used when estimating actively bounded variables","",cg->as_tol,&cg->as_tol,NULL);CHKERRQ(ierr); 28050b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_as_step", "(developer) step length used when estimating actively bounded variables","",cg->as_step,&cg->as_step,NULL);CHKERRQ(ierr); 281484c7b14SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_delta_min", "(developer) minimum scaling factor used for scaled gradient restarts","",cg->delta_min,&cg->delta_min,NULL);CHKERRQ(ierr); 282484c7b14SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_delta_max", "(developer) maximum scaling factor used for scaled gradient restarts","",cg->delta_max,&cg->delta_max,NULL);CHKERRQ(ierr); 28350b47da0SAdam Denchfield 284ac9112b8SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 28550b47da0SAdam Denchfield ierr = MatSetFromOptions(cg->B);CHKERRQ(ierr); 286ac9112b8SAlp Dener PetscFunctionReturn(0); 287ac9112b8SAlp Dener } 288ac9112b8SAlp Dener 289ac9112b8SAlp Dener static PetscErrorCode TaoView_BNCG(Tao tao, PetscViewer viewer) 290ac9112b8SAlp Dener { 291ac9112b8SAlp Dener PetscBool isascii; 292ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 293ac9112b8SAlp Dener PetscErrorCode ierr; 294ac9112b8SAlp Dener 295ac9112b8SAlp Dener PetscFunctionBegin; 296ac9112b8SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 297ac9112b8SAlp Dener if (isascii) { 298ac9112b8SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 299ac9112b8SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cg->cg_type]);CHKERRQ(ierr); 300484c7b14SAdam Denchfield ierr = PetscViewerASCIIPrintf(viewer, "Skipped Stepdirection Updates: %i\n", cg->skipped_updates);CHKERRQ(ierr); 301484c7b14SAdam Denchfield ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %i\n", cg->resets);CHKERRQ(ierr); 302484c7b14SAdam Denchfield ierr = PetscViewerASCIIPrintf(viewer, "Pure gradient steps: %i\n", cg->pure_gd_steps);CHKERRQ(ierr); 303484c7b14SAdam Denchfield ierr = PetscViewerASCIIPrintf(viewer, "Not a descent direction: %i\n", cg->descent_error);CHKERRQ(ierr); 304ac9112b8SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Line search fails: %i\n", cg->ls_fails);CHKERRQ(ierr); 305484c7b14SAdam Denchfield if (cg->diag_scaling){ 306484c7b14SAdam Denchfield ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 307484c7b14SAdam Denchfield if (isascii) { 308484c7b14SAdam Denchfield ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 309484c7b14SAdam Denchfield ierr = MatView(cg->B, viewer);CHKERRQ(ierr); 310484c7b14SAdam Denchfield ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 311484c7b14SAdam Denchfield } 312484c7b14SAdam Denchfield } 313ac9112b8SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 314ac9112b8SAlp Dener } 315ac9112b8SAlp Dener PetscFunctionReturn(0); 316ac9112b8SAlp Dener } 317ac9112b8SAlp Dener 318c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal yty, PetscReal yts, PetscReal sts, PetscReal *scale, PetscReal alpha) 319c8bcdf1eSAdam Denchfield { 320c8bcdf1eSAdam Denchfield PetscReal a, b, c, sig1, sig2; 321c8bcdf1eSAdam Denchfield 322c8bcdf1eSAdam Denchfield PetscFunctionBegin; 323c8bcdf1eSAdam Denchfield *scale = 0.0; 324c8bcdf1eSAdam Denchfield 32550b47da0SAdam Denchfield if (1.0 == alpha){ 326c8bcdf1eSAdam Denchfield *scale = yts/yty; 32750b47da0SAdam Denchfield } else if (0.0 == alpha){ 328c8bcdf1eSAdam Denchfield *scale = sts/yts; 329c8bcdf1eSAdam Denchfield } 33050b47da0SAdam Denchfield else if (-1.0 == alpha) *scale = 1.0; 331c8bcdf1eSAdam Denchfield else { 332c8bcdf1eSAdam Denchfield a = yty; 333c8bcdf1eSAdam Denchfield b = yts; 334c8bcdf1eSAdam Denchfield c = sts; 335c8bcdf1eSAdam Denchfield a *= alpha; 336c8bcdf1eSAdam Denchfield b *= -(2.0*alpha - 1.0); 337c8bcdf1eSAdam Denchfield c *= alpha - 1.0; 338c8bcdf1eSAdam Denchfield sig1 = (-b + PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a); 339c8bcdf1eSAdam Denchfield sig2 = (-b - PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a); 340c8bcdf1eSAdam Denchfield /* accept the positive root as the scalar */ 341c8bcdf1eSAdam Denchfield if (sig1 > 0.0) { 342c8bcdf1eSAdam Denchfield *scale = sig1; 343c8bcdf1eSAdam Denchfield } else if (sig2 > 0.0) { 344c8bcdf1eSAdam Denchfield *scale = sig2; 345c8bcdf1eSAdam Denchfield } else { 346c8bcdf1eSAdam Denchfield SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Cannot find positive scalar"); 347c8bcdf1eSAdam Denchfield } 348c8bcdf1eSAdam Denchfield } 349c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 350c8bcdf1eSAdam Denchfield } 351c8bcdf1eSAdam Denchfield 352ac9112b8SAlp Dener /*MC 353ac9112b8SAlp Dener TAOBNCG - Bound-constrained Nonlinear Conjugate Gradient method. 354ac9112b8SAlp Dener 355ac9112b8SAlp Dener Options Database Keys: 35650b47da0SAdam Denchfield + -tao_bncg_recycle - enable recycling the latest calculated gradient vector in subsequent TaoSolve() calls (currently disabled) 357c4b75bccSAlp Dener . -tao_bncg_eta <r> - restart tolerance 35861be54a6SAlp Dener . -tao_bncg_type <taocg_type> - cg formula 359c4b75bccSAlp Dener . -tao_bncg_as_type <none,bertsekas> - active set estimation method 360c4b75bccSAlp Dener . -tao_bncg_as_tol <r> - tolerance used in Bertsekas active-set estimation 361c4b75bccSAlp Dener . -tao_bncg_as_step <r> - trial step length used in Bertsekas active-set estimation 36250b47da0SAdam 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. 36350b47da0SAdam Denchfield . -tao_bncg_theta <r> - convex combination parameter for the Broyden method 36450b47da0SAdam Denchfield . -tao_bncg_hz_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods 36550b47da0SAdam Denchfield . -tao_bncg_dk_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods 36650b47da0SAdam Denchfield . -tao_bncg_xi <r> - Multiplicative constant of the gamma term in the KD method 36750b47da0SAdam Denchfield . -tao_bncg_hz_theta <r> - Multiplicative constant of the theta term for the HZ method 36850b47da0SAdam Denchfield . -tao_bncg_bfgs_scale <r> - Scaling parameter of the bfgs contribution to the scalar Broyden method 36950b47da0SAdam Denchfield . -tao_bncg_dfp_scale <r> - Scaling parameter of the dfp contribution to the scalar Broyden method 37050b47da0SAdam Denchfield . -tao_bncg_diag_scaling <b> - Whether or not to use diagonal initialization/preconditioning for the CG methods. Default True. 37150b47da0SAdam Denchfield . -tao_bncg_dynamic_restart <b> - use dynamic restart strategy in the HZ, DK, KD methods 37250b47da0SAdam Denchfield . -tao_bncg_unscaled_restart <b> - whether or not to scale the gradient when doing gradient descent restarts 37350b47da0SAdam Denchfield . -tao_bncg_zeta <r> - Scaling parameter in the KD method 374484c7b14SAdam Denchfield . -tao_bncg_delta_min <r> - Minimum bound for rescaling during restarted gradient descent steps 375484c7b14SAdam Denchfield . -tao_bncg_delta_max <r> - Maximum bound for rescaling during restarted gradient descent steps 37650b47da0SAdam Denchfield . -tao_bncg_min_quad <i> - Number of quadratic-like steps in a row necessary to do a dynamic restart 37750b47da0SAdam 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 37850b47da0SAdam Denchfield . -tao_bncg_spaced_restart <b> - whether or not to do gradient descent steps every x*n iterations 379484c7b14SAdam Denchfield . -tao_bncg_no_scaling <b> - If true, eliminates all scaling, including defaults. 38050b47da0SAdam Denchfield . -tao_bncg_neg_xi <b> - Whether or not to use negative xi in the KD method under certain conditions 381ac9112b8SAlp Dener 382ac9112b8SAlp Dener Notes: 383ac9112b8SAlp Dener CG formulas are: 38450b47da0SAdam Denchfield "gd" - Gradient Descent 385ac9112b8SAlp Dener "fr" - Fletcher-Reeves 38650b47da0SAdam Denchfield "pr" - Polak-Ribiere-Polyak 387ac9112b8SAlp Dener "prp" - Polak-Ribiere-Plus 388ac9112b8SAlp Dener "hs" - Hestenes-Steifel 389ac9112b8SAlp Dener "dy" - Dai-Yuan 39050b47da0SAdam Denchfield "ssml_bfgs" - Self-Scaling Memoryless BFGS 39150b47da0SAdam Denchfield "ssml_dfp" - Self-Scaling Memoryless DFP 39250b47da0SAdam Denchfield "ssml_brdn" - Self-Scaling Memoryless Broyden 39350b47da0SAdam Denchfield "hz" - Hager-Zhang (CG_DESCENT 5.3) 39450b47da0SAdam Denchfield "dk" - Dai-Kou (2013) 39550b47da0SAdam Denchfield "kd" - Kou-Dai (2015) 396ac9112b8SAlp Dener Level: beginner 397ac9112b8SAlp Dener M*/ 398ac9112b8SAlp Dener 399ac9112b8SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNCG(Tao tao) 400ac9112b8SAlp Dener { 401ac9112b8SAlp Dener TAO_BNCG *cg; 402ac9112b8SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 403ac9112b8SAlp Dener PetscErrorCode ierr; 404ac9112b8SAlp Dener 405ac9112b8SAlp Dener PetscFunctionBegin; 406ac9112b8SAlp Dener tao->ops->setup = TaoSetUp_BNCG; 407ac9112b8SAlp Dener tao->ops->solve = TaoSolve_BNCG; 408ac9112b8SAlp Dener tao->ops->view = TaoView_BNCG; 409ac9112b8SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNCG; 410ac9112b8SAlp Dener tao->ops->destroy = TaoDestroy_BNCG; 411ac9112b8SAlp Dener 412ac9112b8SAlp Dener /* Override default settings (unless already changed) */ 413ac9112b8SAlp Dener if (!tao->max_it_changed) tao->max_it = 2000; 414ac9112b8SAlp Dener if (!tao->max_funcs_changed) tao->max_funcs = 4000; 415ac9112b8SAlp Dener 416ac9112b8SAlp Dener /* Note: nondefault values should be used for nonlinear conjugate gradient */ 417ac9112b8SAlp Dener /* method. In particular, gtol should be less that 0.5; the value used in */ 418ac9112b8SAlp Dener /* Nocedal and Wright is 0.10. We use the default values for the */ 419ac9112b8SAlp Dener /* linesearch because it seems to work better. */ 420ac9112b8SAlp Dener ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 421ac9112b8SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 422ac9112b8SAlp Dener ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 423ac9112b8SAlp Dener ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr); 424ac9112b8SAlp Dener ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 425ac9112b8SAlp Dener 426ac9112b8SAlp Dener ierr = PetscNewLog(tao,&cg);CHKERRQ(ierr); 427ac9112b8SAlp Dener tao->data = (void*)cg; 428484c7b14SAdam Denchfield ierr = KSPInitializePackage();CHKERRQ(ierr); 42950b47da0SAdam Denchfield ierr = MatCreate(PetscObjectComm((PetscObject)tao), &cg->B);CHKERRQ(ierr); 43050b47da0SAdam Denchfield ierr = PetscObjectIncrementTabLevel((PetscObject)cg->B, (PetscObject)tao, 1);CHKERRQ(ierr); 43150b47da0SAdam Denchfield ierr = MatSetOptionsPrefix(cg->B, "tao_bncg_");CHKERRQ(ierr); 43250b47da0SAdam Denchfield ierr = MatSetType(cg->B, MATLMVMDIAGBRDN);CHKERRQ(ierr); 43350b47da0SAdam Denchfield 434484c7b14SAdam Denchfield cg->pc = NULL; 435484c7b14SAdam Denchfield 43650b47da0SAdam Denchfield cg->dk_eta = 0.5; 43750b47da0SAdam Denchfield cg->hz_eta = 0.4; 438c8bcdf1eSAdam Denchfield cg->dynamic_restart = PETSC_FALSE; 439c8bcdf1eSAdam Denchfield cg->unscaled_restart = PETSC_FALSE; 440484c7b14SAdam Denchfield cg->no_scaling = PETSC_FALSE; 441484c7b14SAdam Denchfield cg->delta_min = 1e-7; 442484c7b14SAdam Denchfield cg->delta_max = 100; 443c8bcdf1eSAdam Denchfield cg->theta = 1.0; 444c8bcdf1eSAdam Denchfield cg->hz_theta = 1.0; 445c8bcdf1eSAdam Denchfield cg->dfp_scale = 1.0; 446c8bcdf1eSAdam Denchfield cg->bfgs_scale = 1.0; 44750b47da0SAdam Denchfield cg->zeta = 0.1; 44850b47da0SAdam Denchfield cg->min_quad = 6; 449c8bcdf1eSAdam Denchfield cg->min_restart_num = 6; /* As in CG_DESCENT and KD2015*/ 450c8bcdf1eSAdam Denchfield cg->xi = 1.0; 45150b47da0SAdam Denchfield cg->neg_xi = PETSC_TRUE; 452c8bcdf1eSAdam Denchfield cg->spaced_restart = PETSC_FALSE; 453c8bcdf1eSAdam Denchfield cg->tol_quad = 1e-8; 45461be54a6SAlp Dener cg->as_step = 0.001; 45561be54a6SAlp Dener cg->as_tol = 0.001; 45650b47da0SAdam Denchfield cg->eps_23 = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); /* Just a little tighter*/ 45761be54a6SAlp Dener cg->as_type = CG_AS_BERTSEKAS; 458c8bcdf1eSAdam Denchfield cg->cg_type = CG_SSML_BFGS; 459c0f10754SAlp Dener cg->recycle = PETSC_FALSE; 460c8bcdf1eSAdam Denchfield cg->alpha = 1.0; 461c8bcdf1eSAdam Denchfield cg->diag_scaling = PETSC_TRUE; 462c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 463c8bcdf1eSAdam Denchfield } 464c8bcdf1eSAdam Denchfield 465c8bcdf1eSAdam Denchfield 466c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGResetUpdate(Tao tao, PetscReal gnormsq) 467c8bcdf1eSAdam Denchfield { 468c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 469c8bcdf1eSAdam Denchfield PetscErrorCode ierr; 470c8bcdf1eSAdam Denchfield PetscReal scaling; 471c8bcdf1eSAdam Denchfield 472c8bcdf1eSAdam Denchfield PetscFunctionBegin; 473c8bcdf1eSAdam Denchfield ++cg->resets; 474c8bcdf1eSAdam Denchfield scaling = 2.0 * PetscMax(1.0, PetscAbsScalar(cg->f)) / PetscMax(gnormsq, cg->eps_23); 475484c7b14SAdam Denchfield scaling = PetscMin(cg->delta_max, PetscMax(cg->delta_min, scaling)); 476484c7b14SAdam Denchfield if (cg->unscaled_restart) { 477484c7b14SAdam Denchfield scaling = 1.0; 478484c7b14SAdam Denchfield ++cg->pure_gd_steps; 479484c7b14SAdam Denchfield } 480c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -scaling, 0.0, tao->gradient);CHKERRQ(ierr); 481c8bcdf1eSAdam Denchfield /* Also want to reset our diagonal scaling with each restart */ 482c8bcdf1eSAdam Denchfield if (cg->diag_scaling) { 48350b47da0SAdam Denchfield ierr = MatLMVMReset(cg->B, PETSC_FALSE);CHKERRQ(ierr); 484c8bcdf1eSAdam Denchfield } 485c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 486c8bcdf1eSAdam Denchfield } 487c8bcdf1eSAdam Denchfield 488c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGCheckDynamicRestart(Tao tao, PetscReal stepsize, PetscReal gd, PetscReal gd_old, PetscBool *dynrestart, PetscReal fold) 489c8bcdf1eSAdam Denchfield { 490c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 491c8bcdf1eSAdam Denchfield PetscReal quadinterp; 492c8bcdf1eSAdam Denchfield 493c8bcdf1eSAdam Denchfield PetscFunctionBegin; 49450b47da0SAdam Denchfield if (cg->f < cg->min_quad/10) { 49550b47da0SAdam Denchfield *dynrestart = PETSC_FALSE; 49650b47da0SAdam Denchfield PetscFunctionReturn(0); 49750b47da0SAdam Denchfield } /* just skip this since this strategy doesn't work well for functions near zero */ 498484c7b14SAdam Denchfield quadinterp = 2.0*(cg->f - fold)/(stepsize*(gd + gd_old)); 49950b47da0SAdam Denchfield if (PetscAbs(quadinterp - 1.0) < cg->tol_quad) ++cg->iter_quad; 500c8bcdf1eSAdam Denchfield else { 501c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 502c8bcdf1eSAdam Denchfield *dynrestart = PETSC_FALSE; 503c8bcdf1eSAdam Denchfield } 504c8bcdf1eSAdam Denchfield if (cg->iter_quad >= cg->min_quad) { 505c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 506c8bcdf1eSAdam Denchfield *dynrestart = PETSC_TRUE; 507c8bcdf1eSAdam Denchfield } 508c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 509c8bcdf1eSAdam Denchfield } 510c8bcdf1eSAdam Denchfield 511484c7b14SAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao tao, PetscReal gnorm2, PetscReal step, PetscReal fold, PetscReal gnorm2_old, PetscReal dnorm, PetscReal ginner, PetscBool pcgd_fallback) 51250b47da0SAdam Denchfield { 513c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 514c8bcdf1eSAdam Denchfield PetscErrorCode ierr; 51550b47da0SAdam Denchfield PetscReal gamma = 1.0, tau_k, beta; 516484c7b14SAdam Denchfield PetscReal tmp = 1.0, ynorm, ynorm2 = 1.0, snorm = 1.0, dk_yk=1.0, gd; 51750b47da0SAdam Denchfield PetscReal gkp1_yk, gd_old, tau_bfgs, tau_dfp, gkp1D_yk, gtDg; 518c8bcdf1eSAdam Denchfield PetscInt dim; 519484c7b14SAdam Denchfield PetscBool cg_restart = PETSC_FALSE; 520c8bcdf1eSAdam Denchfield PetscFunctionBegin; 521c8bcdf1eSAdam Denchfield 52250b47da0SAdam Denchfield /* Local curvature check to see if we need to restart */ 523484c7b14SAdam Denchfield if (tao->niter >= 1 || cg->recycle){ 524c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr); 525c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr); 526c8bcdf1eSAdam Denchfield ynorm2 = ynorm*ynorm; 527c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk);CHKERRQ(ierr); 528484c7b14SAdam Denchfield if (step*dnorm < PETSC_MACHINE_EPSILON || step*dk_yk < PETSC_MACHINE_EPSILON){ 529e2570530SAlp Dener cg_restart = PETSC_TRUE; 530484c7b14SAdam Denchfield ++cg->skipped_updates; 531484c7b14SAdam Denchfield } 53250b47da0SAdam Denchfield if (cg->spaced_restart) { 53350b47da0SAdam Denchfield ierr = VecGetSize(tao->gradient, &dim);CHKERRQ(ierr); 534e2570530SAlp Dener if (tao->niter % (dim*cg->min_restart_num)) cg_restart = PETSC_TRUE; 53550b47da0SAdam Denchfield } 53650b47da0SAdam Denchfield } 53750b47da0SAdam Denchfield /* If the user wants regular restarts, do it every 6n iterations, where n=dimension */ 53850b47da0SAdam Denchfield if (cg->spaced_restart){ 53950b47da0SAdam Denchfield ierr = VecGetSize(tao->gradient, &dim);CHKERRQ(ierr); 540e2570530SAlp Dener if (0 == tao->niter % (6*dim)) cg_restart = PETSC_TRUE; 54150b47da0SAdam Denchfield } 54250b47da0SAdam Denchfield /* Compute the diagonal scaling vector if applicable */ 54350b47da0SAdam Denchfield if (cg->diag_scaling) { 54450b47da0SAdam Denchfield ierr = MatLMVMUpdate(cg->B, tao->solution, tao->gradient);CHKERRQ(ierr); 54550b47da0SAdam Denchfield } 54650b47da0SAdam Denchfield 547484c7b14SAdam Denchfield /* A note on diagonal scaling (to be added to paper): 548484c7b14SAdam Denchfield For the FR, PR, PRP, and DY methods, the diagonally scaled versions 549484c7b14SAdam Denchfield must be derived as a preconditioned CG method rather than as 550484c7b14SAdam Denchfield a Hessian initialization like in the Broyden methods. */ 55150b47da0SAdam Denchfield 552484c7b14SAdam Denchfield /* In that case, one writes the objective function as 553484c7b14SAdam Denchfield f(x) \equiv f(Ay). Gradient evaluations yield g(x_k) = A g(Ay_k) = A g(x_k). 554484c7b14SAdam Denchfield Furthermore, the direction d_k \equiv (x_k - x_{k-1})/step according to 555484c7b14SAdam Denchfield HZ (2006) becomes A^{-1} d_k, such that d_k^T g_k remains the 556484c7b14SAdam Denchfield same under preconditioning. Note that A is diagonal, such that A^T = A. */ 55750b47da0SAdam Denchfield 558484c7b14SAdam Denchfield /* This yields questions like what the dot product d_k^T y_k 559484c7b14SAdam Denchfield should look like. HZ mistakenly treats that as the same under 560484c7b14SAdam Denchfield preconditioning, but that is not necessarily true. */ 56150b47da0SAdam Denchfield 562484c7b14SAdam Denchfield /* Observe y_k \equiv g_k - g_{k-1}, and under the P.C. transformation, 563484c7b14SAdam 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}), 564484c7b14SAdam 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 565484c7b14SAdam Denchfield NOT the same if our preconditioning matrix is updated between iterations. 566484c7b14SAdam Denchfield This same issue is found when considering dot products of the form g_{k+1}^T y_k. */ 56750b47da0SAdam Denchfield 56850b47da0SAdam Denchfield /* Compute CG step direction */ 56950b47da0SAdam Denchfield if (cg_restart) { 57050b47da0SAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 571484c7b14SAdam Denchfield } else if (pcgd_fallback) { 572484c7b14SAdam Denchfield /* Just like preconditioned CG */ 573484c7b14SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 574484c7b14SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work);CHKERRQ(ierr); 57550b47da0SAdam Denchfield } else if (ynorm2 > PETSC_MACHINE_EPSILON) { 57650b47da0SAdam Denchfield switch (cg->cg_type) { 577484c7b14SAdam Denchfield case CG_PCGradientDescent: 57850b47da0SAdam Denchfield if (!cg->diag_scaling){ 579484c7b14SAdam Denchfield if (!cg->no_scaling){ 58050b47da0SAdam Denchfield cg->sts = step*step*dnorm*dnorm; 58150b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 582484c7b14SAdam Denchfield } else { 583484c7b14SAdam Denchfield tau_k = 1.0; 584484c7b14SAdam Denchfield ++cg->pure_gd_steps; 585484c7b14SAdam Denchfield } 58650b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, 0.0, tao->gradient);CHKERRQ(ierr); 58750b47da0SAdam Denchfield } else { 58850b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 58950b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work);CHKERRQ(ierr); 59050b47da0SAdam Denchfield } 59150b47da0SAdam Denchfield break; 592484c7b14SAdam Denchfield 59350b47da0SAdam Denchfield case CG_HestenesStiefel: 59450b47da0SAdam Denchfield /* Classic Hestenes-Stiefel method, modified with scalar and diagonal preconditioning. */ 59550b47da0SAdam Denchfield if (!cg->diag_scaling){ 59650b47da0SAdam Denchfield cg->sts = step*step*dnorm*dnorm; 59750b47da0SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 59850b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 59950b47da0SAdam Denchfield beta = tau_k*gkp1_yk/dk_yk; 60050b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 60150b47da0SAdam Denchfield } else { 60250b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 60350b47da0SAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr); 60450b47da0SAdam Denchfield beta = gkp1_yk/dk_yk; 60550b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 60650b47da0SAdam Denchfield } 607c8bcdf1eSAdam Denchfield break; 608484c7b14SAdam Denchfield 609c8bcdf1eSAdam Denchfield case CG_FletcherReeves: 61050b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr); 61150b47da0SAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr); 61250b47da0SAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr); 61350b47da0SAdam Denchfield ynorm2 = ynorm*ynorm; 61450b47da0SAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk);CHKERRQ(ierr); 61550b47da0SAdam Denchfield if (!cg->diag_scaling){ 61650b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, step*step*dnorm*dnorm, &tau_k, cg->alpha);CHKERRQ(ierr); 61750b47da0SAdam Denchfield beta = tau_k*gnorm2/gnorm2_old; 61850b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 61950b47da0SAdam Denchfield } else { 62050b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->g_work, &gnorm2_old);CHKERRQ(ierr); /* Before it's updated */ 62150b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 62250b47da0SAdam Denchfield ierr = VecDot(tao->gradient, cg->g_work, &tmp);CHKERRQ(ierr); 62350b47da0SAdam Denchfield beta = tmp/gnorm2_old; 62450b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 62550b47da0SAdam Denchfield } 626c8bcdf1eSAdam Denchfield break; 627484c7b14SAdam Denchfield 62850b47da0SAdam Denchfield case CG_PolakRibierePolyak: 62950b47da0SAdam Denchfield snorm = step*dnorm; 63050b47da0SAdam Denchfield if (!cg->diag_scaling){ 63150b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr); 63250b47da0SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 63350b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 63450b47da0SAdam Denchfield beta = tau_k*gkp1_yk/gnorm2_old; 63550b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 63650b47da0SAdam Denchfield } else { 63750b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->g_work, &gnorm2_old);CHKERRQ(ierr); 63850b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 63950b47da0SAdam Denchfield ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 64050b47da0SAdam Denchfield beta = gkp1_yk/gnorm2_old; 64150b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 64250b47da0SAdam Denchfield } 643c8bcdf1eSAdam Denchfield break; 644484c7b14SAdam Denchfield 645c8bcdf1eSAdam Denchfield case CG_PolakRibierePlus: 64650b47da0SAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr); 64750b47da0SAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr); 64850b47da0SAdam Denchfield ynorm2 = ynorm*ynorm; 64950b47da0SAdam Denchfield if (!cg->diag_scaling){ 65050b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr); 65150b47da0SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 65250b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 65350b47da0SAdam Denchfield beta = tau_k*gkp1_yk/gnorm2_old; 65450b47da0SAdam Denchfield beta = PetscMax(beta, 0.0); 65550b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 65650b47da0SAdam Denchfield } else { 65750b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->g_work, &gnorm2_old);CHKERRQ(ierr); /* Old gtDg */ 65850b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 65950b47da0SAdam Denchfield ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 66050b47da0SAdam Denchfield beta = gkp1_yk/gnorm2_old; 66150b47da0SAdam Denchfield beta = PetscMax(beta, 0.0); 66250b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 66350b47da0SAdam Denchfield } 664c8bcdf1eSAdam Denchfield break; 665484c7b14SAdam Denchfield 666484c7b14SAdam Denchfield case CG_DaiYuan: 667484c7b14SAdam Denchfield /* Dai, Yu-Hong, and Yaxiang Yuan. "A nonlinear conjugate gradient method with a strong global convergence property." 668484c7b14SAdam Denchfield SIAM Journal on optimization 10, no. 1 (1999): 177-182. */ 66950b47da0SAdam Denchfield if (!cg->diag_scaling){ 67050b47da0SAdam Denchfield ierr = VecDot(tao->stepdirection, tao->gradient, &gd);CHKERRQ(ierr); 671c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 67250b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->yts, &tau_k, cg->alpha);CHKERRQ(ierr); 67350b47da0SAdam Denchfield beta = tau_k*gnorm2/(gd - gd_old); 67450b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 67550b47da0SAdam Denchfield } else { 676484c7b14SAdam Denchfield ierr = MatMult(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr); 677484c7b14SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 67850b47da0SAdam Denchfield ierr = VecDot(cg->g_work, tao->gradient, >Dg);CHKERRQ(ierr); 67950b47da0SAdam Denchfield ierr = VecDot(tao->stepdirection, cg->G_old, &gd_old);CHKERRQ(ierr); 68050b47da0SAdam Denchfield ierr = VecDot(cg->d_work, cg->g_work, &dk_yk);CHKERRQ(ierr); 68150b47da0SAdam Denchfield dk_yk = dk_yk - gd_old; 68250b47da0SAdam Denchfield beta = gtDg/dk_yk; 68350b47da0SAdam Denchfield ierr = VecScale(cg->d_work, beta); 68450b47da0SAdam Denchfield ierr = VecWAXPY(tao->stepdirection, -1.0, cg->g_work, cg->d_work);CHKERRQ(ierr); 68550b47da0SAdam Denchfield } 686c8bcdf1eSAdam Denchfield break; 687484c7b14SAdam Denchfield 688c8bcdf1eSAdam Denchfield case CG_HagerZhang: 689484c7b14SAdam Denchfield /* Hager, William W., and Hongchao Zhang. "Algorithm 851: CG_DESCENT, a conjugate gradient method with guaranteed descent." 690484c7b14SAdam Denchfield ACM Transactions on Mathematical Software (TOMS) 32, no. 1 (2006): 113-137. */ 691c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 692c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 693c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 69450b47da0SAdam Denchfield snorm = dnorm*step; 69550b47da0SAdam Denchfield cg->yts = step*dk_yk; 696c8bcdf1eSAdam Denchfield if (cg->use_dynamic_restart){ 697c8bcdf1eSAdam Denchfield ierr = TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold);CHKERRQ(ierr); 698c8bcdf1eSAdam Denchfield } 69950b47da0SAdam Denchfield if (cg->dynamic_restart){ 700c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 701c8bcdf1eSAdam Denchfield } else { 702c8bcdf1eSAdam Denchfield if (!cg->diag_scaling){ 703c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 704c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 705c8bcdf1eSAdam Denchfield /* Supplying cg->alpha = -1.0 will give the CG_DESCENT 5.3 special case of tau_k = 1.0 */ 706c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 707c8bcdf1eSAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk)); 708c8bcdf1eSAdam Denchfield /* Bound beta as in CG_DESCENT 5.3, as implemented, with the third comparison from DK 2013 */ 70950b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta*tau_k*gd_old/(dnorm*dnorm)), cg->dk_eta*tau_k*gd/(dnorm*dnorm)); 710c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d */ 71150b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 712c8bcdf1eSAdam Denchfield } else { 713c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 714c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 715c8bcdf1eSAdam Denchfield cg->sts = snorm*snorm; 71650b47da0SAdam Denchfield /* Apply the diagonal scaling to all my vectors */ 71750b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 71850b47da0SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 71950b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr); 720c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 721c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr); 722c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */ 723c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 72450b47da0SAdam Denchfield ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr); 725c8bcdf1eSAdam Denchfield tau_k = -tau_k*gd/(dk_yk*dk_yk); 726c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */ 727484c7b14SAdam Denchfield beta = gkp1_yk/dk_yk + cg->hz_theta*tau_k; /* HZ; (1.15) from DK 2013 */ 728c8bcdf1eSAdam Denchfield /* From HZ2013, modified to account for diagonal scaling*/ 72950b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->d_work, &gd_old);CHKERRQ(ierr); 73050b47da0SAdam Denchfield ierr = VecDot(tao->stepdirection, cg->g_work, &gd);CHKERRQ(ierr); 73150b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta*gd_old/(dnorm*dnorm)), cg->dk_eta*gd/(dnorm*dnorm)); 732c8bcdf1eSAdam Denchfield /* Do the update */ 733484c7b14SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 73450b47da0SAdam Denchfield } 73550b47da0SAdam Denchfield } 736c8bcdf1eSAdam Denchfield break; 737484c7b14SAdam Denchfield 738c8bcdf1eSAdam Denchfield case CG_DaiKou: 739484c7b14SAdam Denchfield /* Dai, Yu-Hong, and Cai-Xia Kou. "A nonlinear conjugate gradient algorithm with an optimal property and an improved Wolfe line search." 740484c7b14SAdam Denchfield SIAM Journal on Optimization 23, no. 1 (2013): 296-320. */ 741c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 742c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 743c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 74450b47da0SAdam Denchfield snorm = step*dnorm; 74550b47da0SAdam Denchfield cg->yts = dk_yk*step; 746c8bcdf1eSAdam Denchfield if (!cg->diag_scaling){ 747c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 74850b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 749c8bcdf1eSAdam Denchfield /* Use cg->alpha = -1.0 to get tau_k = 1.0 as in CG_DESCENT 5.3 */ 750c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 75150b47da0SAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk) + gd/(dnorm*dnorm)) - step*gd/dk_yk; 75250b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta*tau_k*gd_old/(dnorm*dnorm)), cg->dk_eta*tau_k*gd/(dnorm*dnorm)); 753c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d */ 754c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, 0.0, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 755c8bcdf1eSAdam Denchfield } else { 756c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 757c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 758c8bcdf1eSAdam Denchfield cg->sts = snorm*snorm; 75950b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 76050b47da0SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 76150b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr); 762c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 763c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr); 76450b47da0SAdam Denchfield ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr); 765c8bcdf1eSAdam Denchfield tau_k = tau_k*gd/(dk_yk*dk_yk); 766c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 767c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */ 768484c7b14SAdam Denchfield beta = gkp1_yk/dk_yk - step*tmp - tau_k; 769c8bcdf1eSAdam Denchfield /* Update this for the last term in beta */ 770c8bcdf1eSAdam Denchfield ierr = VecDot(cg->y_work, tao->stepdirection, &dk_yk);CHKERRQ(ierr); 771c8bcdf1eSAdam Denchfield beta += tmp*dk_yk/(dnorm*dnorm); /* projection of y_work onto dk */ 77250b47da0SAdam Denchfield ierr = VecDot(tao->stepdirection, cg->g_work, &gd);CHKERRQ(ierr); 77350b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->d_work, &gd_old);CHKERRQ(ierr); 77450b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta*gd_old/(dnorm*dnorm)), cg->dk_eta*gd/(dnorm*dnorm)); 775c8bcdf1eSAdam Denchfield /* Do the update */ 776484c7b14SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 77750b47da0SAdam Denchfield } 778c8bcdf1eSAdam Denchfield break; 779484c7b14SAdam Denchfield 780c8bcdf1eSAdam Denchfield case CG_KouDai: 781484c7b14SAdam Denchfield /* Kou, Cai-Xia, and Yu-Hong Dai. "A modified self-scaling memoryless Broyden–Fletcher–Goldfarb–Shanno method for unconstrained optimization." 782484c7b14SAdam Denchfield Journal of Optimization Theory and Applications 165, no. 1 (2015): 209-224. */ 783c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 784c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 785c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 78650b47da0SAdam Denchfield snorm = step*dnorm; 78750b47da0SAdam Denchfield cg->yts = dk_yk*step; 788c8bcdf1eSAdam Denchfield if (cg->use_dynamic_restart){ 789c8bcdf1eSAdam Denchfield ierr = TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold);CHKERRQ(ierr); 790c8bcdf1eSAdam Denchfield } 79150b47da0SAdam Denchfield if (cg->dynamic_restart){ 792c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 793c8bcdf1eSAdam Denchfield } else { 794c8bcdf1eSAdam Denchfield if (!cg->diag_scaling){ 795c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 796c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 797c8bcdf1eSAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk)) - step*gd/dk_yk; 798c8bcdf1eSAdam Denchfield if (beta < cg->zeta*tau_k*gd/(dnorm*dnorm)) /* 0.1 is KD's zeta parameter */ 799c8bcdf1eSAdam Denchfield { 800c8bcdf1eSAdam Denchfield beta = cg->zeta*tau_k*gd/(dnorm*dnorm); 801c8bcdf1eSAdam Denchfield gamma = 0.0; 802c8bcdf1eSAdam Denchfield } else { 803c8bcdf1eSAdam Denchfield if (gkp1_yk < 0 && cg->neg_xi) gamma = -1.0*gd/dk_yk; 804484c7b14SAdam Denchfield /* This seems to be very effective when there's no tau_k scaling. 805484c7b14SAdam Denchfield This guarantees a large descent step every iteration, going through DK 2015 Lemma 3.1's proof but allowing for negative xi */ 80650b47da0SAdam Denchfield else { 80750b47da0SAdam Denchfield gamma = cg->xi*gd/dk_yk; 80850b47da0SAdam Denchfield } 809c8bcdf1eSAdam Denchfield } 810c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d + t*tmp*yk */ 811c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, gamma*tau_k, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 812c8bcdf1eSAdam Denchfield } else { 813c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 814c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 815c8bcdf1eSAdam Denchfield cg->sts = snorm*snorm; 81650b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 81750b47da0SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 818c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 819c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1D_yk);CHKERRQ(ierr); 820c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */ 821c8bcdf1eSAdam Denchfield gamma = gd/dk_yk; 822c8bcdf1eSAdam Denchfield /* tau_k = -ytDy/(ytd)^2 * gd */ 82350b47da0SAdam Denchfield ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr); 824c8bcdf1eSAdam Denchfield tau_k = tau_k*gd/(dk_yk*dk_yk); 825c8bcdf1eSAdam Denchfield /* beta is the constant which adds the d_k contribution */ 826c8bcdf1eSAdam Denchfield beta = gkp1D_yk/dk_yk - step*gamma - tau_k; 827c8bcdf1eSAdam Denchfield /* Here is the requisite check */ 82850b47da0SAdam Denchfield ierr = VecDot(tao->stepdirection, cg->g_work, &tmp);CHKERRQ(ierr); 829c8bcdf1eSAdam Denchfield if (cg->neg_xi){ 830c8bcdf1eSAdam Denchfield /* modified KD implementation */ 831c8bcdf1eSAdam Denchfield if (gkp1D_yk/dk_yk < 0) gamma = -1.0*gd/dk_yk; 83250b47da0SAdam Denchfield else { 83350b47da0SAdam Denchfield gamma = cg->xi*gd/dk_yk; 83450b47da0SAdam Denchfield } 835c8bcdf1eSAdam Denchfield if (beta < cg->zeta*tmp/(dnorm*dnorm)){ 836c8bcdf1eSAdam Denchfield beta = cg->zeta*tmp/(dnorm*dnorm); 837c8bcdf1eSAdam Denchfield gamma = 0.0; 838c8bcdf1eSAdam Denchfield } 839c8bcdf1eSAdam Denchfield } else { /* original KD 2015 implementation */ 840c8bcdf1eSAdam Denchfield if (beta < cg->zeta*tmp/(dnorm*dnorm)) { 841c8bcdf1eSAdam Denchfield beta = cg->zeta*tmp/(dnorm*dnorm); 842c8bcdf1eSAdam Denchfield gamma = 0.0; 843c8bcdf1eSAdam Denchfield } else { 844c8bcdf1eSAdam Denchfield gamma = cg->xi*gd/dk_yk; 845c8bcdf1eSAdam Denchfield } 846c8bcdf1eSAdam Denchfield } 847c8bcdf1eSAdam Denchfield /* Do the update in two steps */ 848c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 849c8bcdf1eSAdam Denchfield ierr = VecAXPY(tao->stepdirection, gamma, cg->y_work);CHKERRQ(ierr); 85050b47da0SAdam Denchfield } 85150b47da0SAdam Denchfield } 852c8bcdf1eSAdam Denchfield break; 853484c7b14SAdam Denchfield 854484c7b14SAdam Denchfield case CG_SSML_BFGS: 855484c7b14SAdam Denchfield /* Perry, J. M. "A class of conjugate gradient algorithms with a two-step variable-metric memory." 856484c7b14SAdam Denchfield Discussion Papers 269 (1977). */ 857484c7b14SAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 858484c7b14SAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 859484c7b14SAdam Denchfield snorm = step*dnorm; 860484c7b14SAdam Denchfield cg->yts = dk_yk*step; 861484c7b14SAdam Denchfield cg->yty = ynorm2; 862484c7b14SAdam Denchfield cg->sts = snorm*snorm; 863484c7b14SAdam Denchfield if (!cg->diag_scaling){ 864484c7b14SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 865484c7b14SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 866484c7b14SAdam Denchfield tmp = gd/dk_yk; 867484c7b14SAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - cg->yty*gd/(dk_yk*dk_yk)) - step*tmp; 868484c7b14SAdam Denchfield /* d <- -t*g + beta*t*d + t*tmp*yk */ 869484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp*tau_k, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 870484c7b14SAdam Denchfield } else { 871484c7b14SAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 872484c7b14SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 873484c7b14SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 874484c7b14SAdam Denchfield /* compute scalar gamma */ 875484c7b14SAdam Denchfield ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 876484c7b14SAdam Denchfield ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr); 877484c7b14SAdam Denchfield gamma = gd/dk_yk; 878484c7b14SAdam Denchfield /* Compute scalar beta */ 879484c7b14SAdam Denchfield beta = (gkp1_yk/dk_yk - gd*tmp/(dk_yk*dk_yk)) - step*gd/dk_yk; 880484c7b14SAdam Denchfield /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 881484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr); 882484c7b14SAdam Denchfield } 883484c7b14SAdam Denchfield break; 884484c7b14SAdam Denchfield 885484c7b14SAdam Denchfield case CG_SSML_DFP: 886484c7b14SAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 887484c7b14SAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 888484c7b14SAdam Denchfield snorm = step*dnorm; 889484c7b14SAdam Denchfield cg->yts = dk_yk*step; 890484c7b14SAdam Denchfield cg->yty = ynorm2; 891484c7b14SAdam Denchfield cg->sts = snorm*snorm; 892484c7b14SAdam Denchfield if (!cg->diag_scaling){ 893484c7b14SAdam Denchfield /* Instead of a regular convex combination, we will solve a quadratic formula. */ 894484c7b14SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 895484c7b14SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 896484c7b14SAdam Denchfield tau_k = cg->dfp_scale*tau_k; 897484c7b14SAdam Denchfield tmp = tau_k*gkp1_yk/cg->yty; 898484c7b14SAdam Denchfield beta = -step*gd/dk_yk; 899484c7b14SAdam Denchfield /* d <- -t*g + beta*d + tmp*yk */ 900484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 901484c7b14SAdam Denchfield } else { 902484c7b14SAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless DFP step */ 903484c7b14SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 904484c7b14SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 905484c7b14SAdam Denchfield /* compute scalar gamma */ 906484c7b14SAdam Denchfield ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 907484c7b14SAdam Denchfield ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr); 908484c7b14SAdam Denchfield gamma = (gkp1_yk/tmp); 909484c7b14SAdam Denchfield /* Compute scalar beta */ 910484c7b14SAdam Denchfield beta = -step*gd/dk_yk; 911484c7b14SAdam Denchfield /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 912484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr); 913484c7b14SAdam Denchfield } 914484c7b14SAdam Denchfield break; 915484c7b14SAdam Denchfield 916484c7b14SAdam Denchfield case CG_SSML_BROYDEN: 917484c7b14SAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 918484c7b14SAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 919484c7b14SAdam Denchfield snorm = step*dnorm; 920484c7b14SAdam Denchfield cg->yts = step*dk_yk; 921484c7b14SAdam Denchfield cg->yty = ynorm2; 922484c7b14SAdam Denchfield cg->sts = snorm*snorm; 923484c7b14SAdam Denchfield if (!cg->diag_scaling){ 924484c7b14SAdam Denchfield /* Instead of a regular convex combination, we will solve a quadratic formula. */ 925484c7b14SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(cg->yty, step*dk_yk, snorm*snorm, &tau_bfgs, cg->bfgs_scale);CHKERRQ(ierr); 926484c7b14SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(cg->yty, step*dk_yk, snorm*snorm, &tau_dfp, cg->dfp_scale);CHKERRQ(ierr); 927484c7b14SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 928484c7b14SAdam Denchfield tau_k = cg->theta*tau_bfgs + (1.0-cg->theta)*tau_dfp; 929484c7b14SAdam Denchfield /* If bfgs_scale = 1.0, it should reproduce the bfgs tau_bfgs. If bfgs_scale = 0.0, 930484c7b14SAdam Denchfield it should reproduce the tau_dfp scaling. Same with dfp_scale. */ 931484c7b14SAdam Denchfield tmp = cg->theta*tau_bfgs*gd/dk_yk + (1-cg->theta)*tau_dfp*gkp1_yk/cg->yty; 932484c7b14SAdam Denchfield beta = cg->theta*tau_bfgs*(gkp1_yk/dk_yk - cg->yty*gd/(dk_yk*dk_yk)) - step*gd/dk_yk; 933484c7b14SAdam Denchfield /* d <- -t*g + beta*d + tmp*yk */ 934484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 935484c7b14SAdam Denchfield } else { 936484c7b14SAdam Denchfield /* We have diagonal scaling enabled */ 937484c7b14SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 938484c7b14SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 939484c7b14SAdam Denchfield /* compute scalar gamma */ 940484c7b14SAdam Denchfield ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 941484c7b14SAdam Denchfield ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr); 942484c7b14SAdam Denchfield gamma = cg->theta*gd/dk_yk + (1-cg->theta)*(gkp1_yk/tmp); 943484c7b14SAdam Denchfield /* Compute scalar beta */ 944484c7b14SAdam Denchfield beta = cg->theta*(gkp1_yk/dk_yk - gd*tmp/(dk_yk*dk_yk)) - step*gd/dk_yk; 945484c7b14SAdam Denchfield /* Compute stepdirection dkp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 946484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr); 947484c7b14SAdam Denchfield } 948484c7b14SAdam Denchfield break; 949484c7b14SAdam Denchfield 950c8bcdf1eSAdam Denchfield default: 951c8bcdf1eSAdam Denchfield beta = 0.0; 952c8bcdf1eSAdam Denchfield break; 953484c7b14SAdam Denchfield 954c8bcdf1eSAdam Denchfield } 955c8bcdf1eSAdam Denchfield } 956c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 957c8bcdf1eSAdam Denchfield } 958c8bcdf1eSAdam Denchfield 959c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao tao, PetscReal gnorm) 960c8bcdf1eSAdam Denchfield { 961c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 962c8bcdf1eSAdam Denchfield PetscErrorCode ierr; 963c8bcdf1eSAdam Denchfield TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 964e2570530SAlp Dener PetscReal step=1.0,gnorm2,gd,ginner=0.0,dnorm=0.0; 965c8bcdf1eSAdam Denchfield PetscReal gnorm2_old,f_old,resnorm, gnorm_old; 966484c7b14SAdam Denchfield PetscBool gd_fallback = PETSC_FALSE, pcgd_fallback = PETSC_FALSE; 967c8bcdf1eSAdam Denchfield 968c8bcdf1eSAdam Denchfield PetscFunctionBegin; 969c8bcdf1eSAdam Denchfield /* We are now going to perform a line search along the direction. */ 970c8bcdf1eSAdam Denchfield /* Store solution and gradient info before it changes */ 971c8bcdf1eSAdam Denchfield ierr = VecCopy(tao->solution, cg->X_old);CHKERRQ(ierr); 972c8bcdf1eSAdam Denchfield ierr = VecCopy(tao->gradient, cg->G_old);CHKERRQ(ierr); 973c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->unprojected_gradient, cg->unprojected_gradient_old);CHKERRQ(ierr); 974c8bcdf1eSAdam Denchfield 975c8bcdf1eSAdam Denchfield gnorm_old = gnorm; 976c8bcdf1eSAdam Denchfield gnorm2_old = gnorm_old*gnorm_old; 977c8bcdf1eSAdam Denchfield f_old = cg->f; 978484c7b14SAdam Denchfield /* Perform bounded line search. If we are recycling a solution from a previous */ 979484c7b14SAdam Denchfield /* TaoSolve, then we want to immediately skip to calculating a new direction rather than performing a linesearch */ 980484c7b14SAdam Denchfield if (!(cg->recycle && 0 == tao->niter)){ 981484c7b14SAdam Denchfield /* Above logic: the below code happens every iteration, except for the first iteration of a recycled TaoSolve */ 982c8bcdf1eSAdam Denchfield ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 983c8bcdf1eSAdam Denchfield ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 984c8bcdf1eSAdam Denchfield ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 985c8bcdf1eSAdam Denchfield 986c8bcdf1eSAdam Denchfield /* Check linesearch failure */ 987c8bcdf1eSAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 988c8bcdf1eSAdam Denchfield ++cg->ls_fails; 989c8bcdf1eSAdam Denchfield if (cg->cg_type == CG_GradientDescent || gd_fallback){ 990c8bcdf1eSAdam Denchfield /* Nothing left to do but fail out of the optimization */ 991c8bcdf1eSAdam Denchfield step = 0.0; 992c8bcdf1eSAdam Denchfield tao->reason = TAO_DIVERGED_LS_FAILURE; 993c8bcdf1eSAdam Denchfield } else { 994484c7b14SAdam Denchfield /* Restore previous point, perform preconditioned GD and regular GD steps at the last good point */ 995c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->X_old, tao->solution);CHKERRQ(ierr); 996c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->G_old, tao->gradient);CHKERRQ(ierr); 997c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->unprojected_gradient_old, cg->unprojected_gradient);CHKERRQ(ierr); 998c8bcdf1eSAdam Denchfield gnorm = gnorm_old; 999c8bcdf1eSAdam Denchfield gnorm2 = gnorm2_old; 1000c8bcdf1eSAdam Denchfield cg->f = f_old; 1001c8bcdf1eSAdam Denchfield 1002484c7b14SAdam Denchfield /* Fall back on preconditioned CG (so long as you're not already using it) */ 1003484c7b14SAdam Denchfield if (cg->cg_type != CG_PCGradientDescent && cg->diag_scaling){ 1004e2570530SAlp Dener pcgd_fallback = PETSC_TRUE; 1005484c7b14SAdam Denchfield ierr = TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, ginner, pcgd_fallback);CHKERRQ(ierr); 1006484c7b14SAdam Denchfield 100750b47da0SAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 1008c8bcdf1eSAdam Denchfield ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1009c8bcdf1eSAdam Denchfield 1010c8bcdf1eSAdam Denchfield ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 1011c8bcdf1eSAdam Denchfield ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 1012c8bcdf1eSAdam Denchfield ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 1013c8bcdf1eSAdam Denchfield 1014484c7b14SAdam Denchfield pcgd_fallback = PETSC_FALSE; 1015484c7b14SAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER){ 1016484c7b14SAdam Denchfield /* Going to perform a regular gradient descent step. */ 1017484c7b14SAdam Denchfield ++cg->ls_fails; 1018484c7b14SAdam Denchfield step = 0.0; 1019e2570530SAlp Dener gd_fallback = PETSC_TRUE; 1020484c7b14SAdam Denchfield } 1021484c7b14SAdam Denchfield } 1022484c7b14SAdam Denchfield /* Fall back on the scaled gradient step */ 1023484c7b14SAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER){ 1024484c7b14SAdam Denchfield ++cg->ls_fails; 1025484c7b14SAdam Denchfield gd_fallback = PETSC_TRUE; 1026484c7b14SAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 1027484c7b14SAdam Denchfield ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1028484c7b14SAdam Denchfield ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 1029484c7b14SAdam Denchfield ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 1030484c7b14SAdam Denchfield ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 1031484c7b14SAdam Denchfield } 1032484c7b14SAdam Denchfield 1033c8bcdf1eSAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER){ 1034c8bcdf1eSAdam Denchfield /* Nothing left to do but fail out of the optimization */ 103550b47da0SAdam Denchfield ++cg->ls_fails; 1036c8bcdf1eSAdam Denchfield step = 0.0; 1037c8bcdf1eSAdam Denchfield tao->reason = TAO_DIVERGED_LS_FAILURE; 1038484c7b14SAdam Denchfield } else { 1039484c7b14SAdam Denchfield /* One of the fallbacks worked. Set them both back equal to false. */ 1040484c7b14SAdam Denchfield gd_fallback = PETSC_FALSE; 1041484c7b14SAdam Denchfield pcgd_fallback = PETSC_FALSE; 1042c8bcdf1eSAdam Denchfield } 1043c8bcdf1eSAdam Denchfield } 1044c8bcdf1eSAdam Denchfield } 1045c8bcdf1eSAdam Denchfield /* Convergence test for line search failure */ 1046c8bcdf1eSAdam Denchfield if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 1047c8bcdf1eSAdam Denchfield 1048c8bcdf1eSAdam Denchfield /* Standard convergence test */ 1049c8bcdf1eSAdam Denchfield ierr = VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W);CHKERRQ(ierr); 1050c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr); 1051c8bcdf1eSAdam Denchfield if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 1052c8bcdf1eSAdam Denchfield ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 1053c8bcdf1eSAdam Denchfield ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr); 1054c8bcdf1eSAdam Denchfield ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 1055c8bcdf1eSAdam Denchfield if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 1056484c7b14SAdam Denchfield } 1057c8bcdf1eSAdam Denchfield /* Assert we have an updated step and we need at least one more iteration. */ 1058c8bcdf1eSAdam Denchfield /* Calculate the next direction */ 1059c8bcdf1eSAdam Denchfield /* Estimate the active set at the new solution */ 1060c8bcdf1eSAdam Denchfield ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr); 1061c8bcdf1eSAdam Denchfield /* Compute the projected gradient and its norm */ 1062c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 1063c8bcdf1eSAdam Denchfield ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr); 1064c8bcdf1eSAdam Denchfield ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 1065c8bcdf1eSAdam Denchfield gnorm2 = gnorm*gnorm; 1066c8bcdf1eSAdam Denchfield 1067484c7b14SAdam Denchfield /* Calculate some quantities used in the StepDirectionUpdate. */ 1068c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, cg->G_old, &ginner);CHKERRQ(ierr); 1069c8bcdf1eSAdam Denchfield ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm);CHKERRQ(ierr); 1070484c7b14SAdam Denchfield /* Update the step direction. */ 1071484c7b14SAdam Denchfield ierr = TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, ginner, pcgd_fallback);CHKERRQ(ierr); 1072484c7b14SAdam Denchfield ++tao->niter; 1073c8bcdf1eSAdam Denchfield ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1074c8bcdf1eSAdam Denchfield 1075c8bcdf1eSAdam Denchfield if (cg->cg_type != CG_GradientDescent) { 1076c8bcdf1eSAdam Denchfield /* Figure out which previously active variables became inactive this iteration */ 1077c8bcdf1eSAdam Denchfield ierr = ISDestroy(&cg->new_inactives);CHKERRQ(ierr); 1078c8bcdf1eSAdam Denchfield if (cg->inactive_idx && cg->inactive_old) { 1079c8bcdf1eSAdam Denchfield ierr = ISDifference(cg->inactive_idx, cg->inactive_old, &cg->new_inactives);CHKERRQ(ierr); 1080c8bcdf1eSAdam Denchfield } 1081c8bcdf1eSAdam Denchfield /* Selectively reset the CG step those freshly inactive variables */ 1082c8bcdf1eSAdam Denchfield if (cg->new_inactives) { 1083c8bcdf1eSAdam Denchfield ierr = VecGetSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr); 1084c8bcdf1eSAdam Denchfield ierr = VecGetSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr); 1085c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->inactive_grad, cg->inactive_step);CHKERRQ(ierr); 1086c8bcdf1eSAdam Denchfield ierr = VecScale(cg->inactive_step, -1.0);CHKERRQ(ierr); 1087c8bcdf1eSAdam Denchfield ierr = VecRestoreSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr); 1088c8bcdf1eSAdam Denchfield ierr = VecRestoreSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr); 1089c8bcdf1eSAdam Denchfield } 1090c8bcdf1eSAdam Denchfield /* Verify that this is a descent direction */ 1091c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 109250b47da0SAdam Denchfield ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm);CHKERRQ(ierr); 109350b47da0SAdam Denchfield if (PetscIsInfOrNanReal(gd) || (gd/(dnorm*dnorm) <= -1e10 || gd/(dnorm*dnorm) >= -1e-10)) { 1094c8bcdf1eSAdam Denchfield /* Not a descent direction, so we reset back to projected gradient descent */ 109550b47da0SAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 1096c8bcdf1eSAdam Denchfield ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1097c8bcdf1eSAdam Denchfield ++cg->descent_error; 1098c8bcdf1eSAdam Denchfield gd_fallback = PETSC_TRUE; 1099c8bcdf1eSAdam Denchfield } else { 1100c8bcdf1eSAdam Denchfield gd_fallback = PETSC_FALSE; 1101c8bcdf1eSAdam Denchfield } 1102c8bcdf1eSAdam Denchfield } 1103ac9112b8SAlp Dener PetscFunctionReturn(0); 1104ac9112b8SAlp Dener } 1105484c7b14SAdam Denchfield 1106484c7b14SAdam Denchfield PetscErrorCode TaoBNCGSetH0(Tao tao, Mat H0) 1107484c7b14SAdam Denchfield { 1108484c7b14SAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 1109484c7b14SAdam Denchfield PetscErrorCode ierr; 1110484c7b14SAdam Denchfield 1111484c7b14SAdam Denchfield PetscFunctionBegin; 1112484c7b14SAdam Denchfield ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 1113484c7b14SAdam Denchfield cg->pc = H0; 1114484c7b14SAdam Denchfield PetscFunctionReturn(0); 1115484c7b14SAdam Denchfield } 1116