1ac9112b8SAlp Dener #include <petsctaolinesearch.h> 2414d97d3SAlp Dener #include <../src/tao/bound/impls/bncg/bncg.h> /*I "petsctao.h" I*/ 350b47da0SAdam Denchfield #include <petscksp.h> 4ac9112b8SAlp Dener 5c8bcdf1eSAdam Denchfield #define CG_GradientDescent 0 6c8bcdf1eSAdam Denchfield #define CG_HestenesStiefel 1 7c8bcdf1eSAdam Denchfield #define CG_FletcherReeves 2 850b47da0SAdam Denchfield #define CG_PolakRibierePolyak 3 9c8bcdf1eSAdam Denchfield #define CG_PolakRibierePlus 4 10c8bcdf1eSAdam Denchfield #define CG_DaiYuan 5 11c8bcdf1eSAdam Denchfield #define CG_HagerZhang 6 12c8bcdf1eSAdam Denchfield #define CG_DaiKou 7 13c8bcdf1eSAdam Denchfield #define CG_KouDai 8 14c8bcdf1eSAdam Denchfield #define CG_SSML_BFGS 9 15c8bcdf1eSAdam Denchfield #define CG_SSML_DFP 10 16c8bcdf1eSAdam Denchfield #define CG_SSML_BROYDEN 11 17484c7b14SAdam Denchfield #define CG_PCGradientDescent 12 18484c7b14SAdam Denchfield #define CGTypes 13 19ac9112b8SAlp Dener 20484c7b14SAdam Denchfield static const char *CG_Table[64] = {"gd", "hs", "fr", "pr", "prp", "dy", "hz", "dk", "kd", "ssml_bfgs", "ssml_dfp", "ssml_brdn", "pcgd"}; 21ac9112b8SAlp Dener 2261be54a6SAlp Dener #define CG_AS_NONE 0 2361be54a6SAlp Dener #define CG_AS_BERTSEKAS 1 2461be54a6SAlp Dener #define CG_AS_SIZE 2 25ac9112b8SAlp Dener 2661be54a6SAlp Dener static const char *CG_AS_TYPE[64] = {"none", "bertsekas"}; 27ac9112b8SAlp Dener 2861be54a6SAlp Dener PetscErrorCode TaoBNCGEstimateActiveSet(Tao tao, PetscInt asType) 2961be54a6SAlp Dener { 3061be54a6SAlp Dener PetscErrorCode ierr; 3161be54a6SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 3261be54a6SAlp Dener 3361be54a6SAlp Dener PetscFunctionBegin; 3461be54a6SAlp Dener ierr = ISDestroy(&cg->inactive_old);CHKERRQ(ierr); 3561be54a6SAlp Dener if (cg->inactive_idx) { 3661be54a6SAlp Dener ierr = ISDuplicate(cg->inactive_idx, &cg->inactive_old);CHKERRQ(ierr); 3761be54a6SAlp Dener ierr = ISCopy(cg->inactive_idx, cg->inactive_old);CHKERRQ(ierr); 3861be54a6SAlp Dener } 3961be54a6SAlp Dener switch (asType) { 4061be54a6SAlp Dener case CG_AS_NONE: 4161be54a6SAlp Dener ierr = ISDestroy(&cg->inactive_idx);CHKERRQ(ierr); 4261be54a6SAlp Dener ierr = VecWhichInactive(tao->XL, tao->solution, cg->unprojected_gradient, tao->XU, PETSC_TRUE, &cg->inactive_idx);CHKERRQ(ierr); 4361be54a6SAlp Dener ierr = ISDestroy(&cg->active_idx);CHKERRQ(ierr); 4461be54a6SAlp Dener ierr = ISComplementVec(cg->inactive_idx, tao->solution, &cg->active_idx);CHKERRQ(ierr); 4561be54a6SAlp Dener break; 4661be54a6SAlp Dener 4761be54a6SAlp Dener case CG_AS_BERTSEKAS: 4861be54a6SAlp Dener /* Use gradient descent to estimate the active set */ 4961be54a6SAlp Dener ierr = VecCopy(cg->unprojected_gradient, cg->W);CHKERRQ(ierr); 5061be54a6SAlp Dener ierr = VecScale(cg->W, -1.0);CHKERRQ(ierr); 5189da521bSAlp Dener ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, cg->unprojected_gradient, cg->W, cg->work, cg->as_step, &cg->as_tol, 5289da521bSAlp Dener &cg->active_lower, &cg->active_upper, &cg->active_fixed, &cg->active_idx, &cg->inactive_idx);CHKERRQ(ierr); 53c4b75bccSAlp Dener break; 5461be54a6SAlp Dener 5561be54a6SAlp Dener default: 5661be54a6SAlp Dener break; 5761be54a6SAlp Dener } 5861be54a6SAlp Dener PetscFunctionReturn(0); 5961be54a6SAlp Dener } 6061be54a6SAlp Dener 61a1318120SAlp Dener PetscErrorCode TaoBNCGBoundStep(Tao tao, PetscInt asType, Vec step) 6261be54a6SAlp Dener { 6361be54a6SAlp Dener PetscErrorCode ierr; 6461be54a6SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 6561be54a6SAlp Dener 6661be54a6SAlp Dener PetscFunctionBegin; 67a1318120SAlp Dener switch (asType) { 6861be54a6SAlp Dener case CG_AS_NONE: 69c4b75bccSAlp Dener ierr = VecISSet(step, cg->active_idx, 0.0);CHKERRQ(ierr); 7061be54a6SAlp Dener break; 7161be54a6SAlp Dener 7261be54a6SAlp Dener case CG_AS_BERTSEKAS: 73c4b75bccSAlp Dener ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, cg->active_lower, cg->active_upper, cg->active_fixed, 1.0, step);CHKERRQ(ierr); 7461be54a6SAlp Dener break; 7561be54a6SAlp Dener 7661be54a6SAlp Dener default: 7761be54a6SAlp Dener break; 7861be54a6SAlp Dener } 7961be54a6SAlp Dener PetscFunctionReturn(0); 8061be54a6SAlp Dener } 8161be54a6SAlp Dener 82ac9112b8SAlp Dener static PetscErrorCode TaoSolve_BNCG(Tao tao) 83ac9112b8SAlp Dener { 84ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 85ac9112b8SAlp Dener PetscErrorCode ierr; 86484c7b14SAdam Denchfield PetscReal step=1.0,gnorm,gnorm2, resnorm; 87c4b75bccSAlp Dener PetscInt nDiff; 88ac9112b8SAlp Dener 89ac9112b8SAlp Dener PetscFunctionBegin; 90ac9112b8SAlp Dener /* Project the current point onto the feasible set */ 91ac9112b8SAlp Dener ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 92cd929ea3SAlp Dener ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 93ac9112b8SAlp Dener 94c8bcdf1eSAdam Denchfield /* Project the initial point onto the feasible region */ 95c8bcdf1eSAdam Denchfield ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr); 96484c7b14SAdam Denchfield 97414d97d3SAlp Dener if (nDiff > 0 || !tao->recycle) { 98c0f10754SAlp Dener ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient);CHKERRQ(ierr); 99484c7b14SAdam Denchfield } 100ac9112b8SAlp Dener ierr = VecNorm(cg->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr); 101*3c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(cg->f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 102ac9112b8SAlp Dener 10361be54a6SAlp Dener /* Estimate the active set and compute the projected gradient */ 10461be54a6SAlp Dener ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr); 10561be54a6SAlp Dener 106ac9112b8SAlp Dener /* Project the gradient and calculate the norm */ 10761be54a6SAlp Dener ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 10861be54a6SAlp Dener ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr); 109ac9112b8SAlp Dener ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 110ac9112b8SAlp Dener gnorm2 = gnorm*gnorm; 111ac9112b8SAlp Dener 112c8bcdf1eSAdam Denchfield /* Initialize counters */ 113e031d6f5SAlp Dener tao->niter = 0; 11450b47da0SAdam Denchfield cg->ls_fails = cg->descent_error = 0; 115c8bcdf1eSAdam Denchfield cg->resets = -1; 116484c7b14SAdam Denchfield cg->skipped_updates = cg->pure_gd_steps = 0; 117c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 118c8bcdf1eSAdam Denchfield 119c8bcdf1eSAdam Denchfield /* Convergence test at the starting point. */ 120ac9112b8SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 121484c7b14SAdam Denchfield 122484c7b14SAdam Denchfield ierr = VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W);CHKERRQ(ierr); 123484c7b14SAdam Denchfield ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr); 124*3c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(resnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 125484c7b14SAdam Denchfield ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 126484c7b14SAdam Denchfield ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr); 127ac9112b8SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 128ac9112b8SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 129484c7b14SAdam Denchfield /* Calculate initial direction. */ 130414d97d3SAlp Dener if (!tao->recycle) { 131484c7b14SAdam Denchfield /* We are not recycling a solution/history from a past TaoSolve */ 132484c7b14SAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 133ac9112b8SAlp Dener } 134c8bcdf1eSAdam Denchfield /* Initial gradient descent step. Scaling by 1.0 also does a decent job for some problems. */ 135c8bcdf1eSAdam Denchfield while (1) { 136e1e80dc8SAlp Dener /* Call general purpose update function */ 137e1e80dc8SAlp Dener if (tao->ops->update) { 1388fcddce6SStefano Zampini ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr); 139e1e80dc8SAlp Dener } 140c8bcdf1eSAdam Denchfield ierr = TaoBNCGConductIteration(tao, gnorm);CHKERRQ(ierr); 141c8bcdf1eSAdam Denchfield if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 142ac9112b8SAlp Dener } 143ac9112b8SAlp Dener } 144ac9112b8SAlp Dener 145ac9112b8SAlp Dener static PetscErrorCode TaoSetUp_BNCG(Tao tao) 146ac9112b8SAlp Dener { 147ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 148ac9112b8SAlp Dener PetscErrorCode ierr; 149ac9112b8SAlp Dener 150ac9112b8SAlp Dener PetscFunctionBegin; 151c4b75bccSAlp Dener if (!tao->gradient) { 152c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 153c4b75bccSAlp Dener } 154c4b75bccSAlp Dener if (!tao->stepdirection) { 155c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); 156c4b75bccSAlp Dener } 157c4b75bccSAlp Dener if (!cg->W) { 158c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&cg->W);CHKERRQ(ierr); 159c4b75bccSAlp Dener } 160c4b75bccSAlp Dener if (!cg->work) { 161c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&cg->work);CHKERRQ(ierr); 162c4b75bccSAlp Dener } 163c8bcdf1eSAdam Denchfield if (!cg->sk) { 164c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->sk);CHKERRQ(ierr); 165c8bcdf1eSAdam Denchfield } 166c8bcdf1eSAdam Denchfield if (!cg->yk) { 167c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->gradient,&cg->yk);CHKERRQ(ierr); 168c8bcdf1eSAdam Denchfield } 169c4b75bccSAlp Dener if (!cg->X_old) { 170c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&cg->X_old);CHKERRQ(ierr); 171c4b75bccSAlp Dener } 172c4b75bccSAlp Dener if (!cg->G_old) { 173c4b75bccSAlp Dener ierr = VecDuplicate(tao->gradient,&cg->G_old);CHKERRQ(ierr); 174c8bcdf1eSAdam Denchfield } 175c8bcdf1eSAdam Denchfield if (cg->diag_scaling) { 176c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->d_work);CHKERRQ(ierr); 177c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->y_work);CHKERRQ(ierr); 17850b47da0SAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->g_work);CHKERRQ(ierr); 179c4b75bccSAlp Dener } 180c4b75bccSAlp Dener if (!cg->unprojected_gradient) { 181c4b75bccSAlp Dener ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient);CHKERRQ(ierr); 182c4b75bccSAlp Dener } 183c4b75bccSAlp Dener if (!cg->unprojected_gradient_old) { 184c4b75bccSAlp Dener ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient_old);CHKERRQ(ierr); 185c4b75bccSAlp Dener } 18650b47da0SAdam Denchfield ierr = MatLMVMAllocate(cg->B, cg->sk, cg->yk);CHKERRQ(ierr); 187484c7b14SAdam Denchfield if (cg->pc) { 188484c7b14SAdam Denchfield ierr = MatLMVMSetJ0(cg->B, cg->pc);CHKERRQ(ierr); 189484c7b14SAdam Denchfield } 190ac9112b8SAlp Dener PetscFunctionReturn(0); 191ac9112b8SAlp Dener } 192ac9112b8SAlp Dener 193ac9112b8SAlp Dener static PetscErrorCode TaoDestroy_BNCG(Tao tao) 194ac9112b8SAlp Dener { 195ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*) tao->data; 196ac9112b8SAlp Dener PetscErrorCode ierr; 197ac9112b8SAlp Dener 198ac9112b8SAlp Dener PetscFunctionBegin; 199ac9112b8SAlp Dener if (tao->setupcalled) { 20061be54a6SAlp Dener ierr = VecDestroy(&cg->W);CHKERRQ(ierr); 201c4b75bccSAlp Dener ierr = VecDestroy(&cg->work);CHKERRQ(ierr); 202ac9112b8SAlp Dener ierr = VecDestroy(&cg->X_old);CHKERRQ(ierr); 203ac9112b8SAlp Dener ierr = VecDestroy(&cg->G_old);CHKERRQ(ierr); 204ac9112b8SAlp Dener ierr = VecDestroy(&cg->unprojected_gradient);CHKERRQ(ierr); 205ac9112b8SAlp Dener ierr = VecDestroy(&cg->unprojected_gradient_old);CHKERRQ(ierr); 20650b47da0SAdam Denchfield ierr = VecDestroy(&cg->g_work);CHKERRQ(ierr); 20750b47da0SAdam Denchfield ierr = VecDestroy(&cg->d_work);CHKERRQ(ierr); 20850b47da0SAdam Denchfield ierr = VecDestroy(&cg->y_work);CHKERRQ(ierr); 20950b47da0SAdam Denchfield ierr = VecDestroy(&cg->sk);CHKERRQ(ierr); 21050b47da0SAdam Denchfield ierr = VecDestroy(&cg->yk);CHKERRQ(ierr); 211ac9112b8SAlp Dener } 212ca964c49SAlp Dener ierr = ISDestroy(&cg->active_lower);CHKERRQ(ierr); 213ca964c49SAlp Dener ierr = ISDestroy(&cg->active_upper);CHKERRQ(ierr); 214ca964c49SAlp Dener ierr = ISDestroy(&cg->active_fixed);CHKERRQ(ierr); 215ca964c49SAlp Dener ierr = ISDestroy(&cg->active_idx);CHKERRQ(ierr); 216ca964c49SAlp Dener ierr = ISDestroy(&cg->inactive_idx);CHKERRQ(ierr); 217ca964c49SAlp Dener ierr = ISDestroy(&cg->inactive_old);CHKERRQ(ierr); 218ca964c49SAlp Dener ierr = ISDestroy(&cg->new_inactives);CHKERRQ(ierr); 21901e1e359SAlp Dener ierr = MatDestroy(&cg->B);CHKERRQ(ierr); 220484c7b14SAdam Denchfield if (cg->pc) { 22101e1e359SAlp Dener ierr = MatDestroy(&cg->pc);CHKERRQ(ierr); 222484c7b14SAdam Denchfield } 223ac9112b8SAlp Dener ierr = PetscFree(tao->data);CHKERRQ(ierr); 224ac9112b8SAlp Dener PetscFunctionReturn(0); 225ac9112b8SAlp Dener } 226ac9112b8SAlp Dener 227ac9112b8SAlp Dener static PetscErrorCode TaoSetFromOptions_BNCG(PetscOptionItems *PetscOptionsObject,Tao tao) 228ac9112b8SAlp Dener { 229ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 230ac9112b8SAlp Dener PetscErrorCode ierr; 231ac9112b8SAlp Dener 232ac9112b8SAlp Dener PetscFunctionBegin; 233ac9112b8SAlp Dener ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Conjugate Gradient method for unconstrained optimization");CHKERRQ(ierr); 234484c7b14SAdam Denchfield ierr = PetscOptionsEList("-tao_bncg_type","cg formula", "", CG_Table, CGTypes, CG_Table[cg->cg_type], &cg->cg_type,NULL);CHKERRQ(ierr); 2358ebe3e4eSStefano Zampini if (cg->cg_type != CG_SSML_BFGS) cg->alpha = -1.0; /* Setting defaults for non-BFGS methods. User can change it below. */ 236484c7b14SAdam Denchfield if (CG_GradientDescent == cg->cg_type) { 237484c7b14SAdam Denchfield cg->cg_type = CG_PCGradientDescent; 238484c7b14SAdam Denchfield /* Set scaling equal to none or, at best, scalar scaling. */ 239484c7b14SAdam Denchfield cg->unscaled_restart = PETSC_TRUE; 240484c7b14SAdam Denchfield cg->diag_scaling = PETSC_FALSE; 241484c7b14SAdam Denchfield } 24250b47da0SAdam 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); 24350b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_hz_eta","(developer) cutoff tolerance for HZ", "", cg->hz_eta,&cg->hz_eta,NULL);CHKERRQ(ierr); 24450b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_eps","(developer) cutoff value for restarts", "", cg->epsilon,&cg->epsilon,NULL);CHKERRQ(ierr); 24550b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_dk_eta","(developer) cutoff tolerance for DK", "", cg->dk_eta,&cg->dk_eta,NULL);CHKERRQ(ierr); 24650b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_xi","(developer) Parameter in the KD method", "", cg->xi,&cg->xi,NULL);CHKERRQ(ierr); 24750b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_theta", "(developer) update parameter for the Broyden method", "", cg->theta, &cg->theta, NULL);CHKERRQ(ierr); 24850b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_hz_theta", "(developer) parameter for the HZ (2006) method", "", cg->hz_theta, &cg->hz_theta, NULL);CHKERRQ(ierr); 24950b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_alpha","(developer) parameter for the scalar scaling","",cg->alpha,&cg->alpha,NULL);CHKERRQ(ierr); 250c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_bfgs_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->bfgs_scale, &cg->bfgs_scale, NULL);CHKERRQ(ierr); 251c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_dfp_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->dfp_scale, &cg->dfp_scale, NULL);CHKERRQ(ierr); 252c8bcdf1eSAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_diag_scaling","Enable diagonal Broyden-like preconditioning","",cg->diag_scaling,&cg->diag_scaling,NULL);CHKERRQ(ierr); 25350b47da0SAdam 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); 25450b47da0SAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_unscaled_restart","(developer) use unscaled gradient restarts","",cg->unscaled_restart,&cg->unscaled_restart,NULL);CHKERRQ(ierr); 25550b47da0SAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_zeta", "(developer) Free parameter for the Kou-Dai method", "", cg->zeta, &cg->zeta, NULL);CHKERRQ(ierr); 256c8bcdf1eSAdam 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); 25750b47da0SAdam 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); 25850b47da0SAdam 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); 259484c7b14SAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_no_scaling","Disable all scaling except in restarts","",cg->no_scaling,&cg->no_scaling,NULL);CHKERRQ(ierr); 260484c7b14SAdam Denchfield if (cg->no_scaling) { 261484c7b14SAdam Denchfield cg->diag_scaling = PETSC_FALSE; 262484c7b14SAdam Denchfield cg->alpha = -1.0; 263484c7b14SAdam Denchfield } 264b474139fSKarl Rupp if (cg->alpha == -1.0 && cg->cg_type == CG_KouDai && !cg->diag_scaling) { /* Some more default options that appear to be good. */ 265484c7b14SAdam Denchfield cg->neg_xi = PETSC_TRUE; 266484c7b14SAdam Denchfield } 26750b47da0SAdam 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); 26850b47da0SAdam 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); 26950b47da0SAdam 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); 270484c7b14SAdam 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); 271484c7b14SAdam 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); 27250b47da0SAdam Denchfield 273ac9112b8SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 2748ebe3e4eSStefano Zampini ierr = MatSetOptionsPrefix(cg->B, ((PetscObject)tao)->prefix);CHKERRQ(ierr); 2758ebe3e4eSStefano Zampini ierr = MatAppendOptionsPrefix(cg->B, "tao_bncg_");CHKERRQ(ierr); 27650b47da0SAdam Denchfield ierr = MatSetFromOptions(cg->B);CHKERRQ(ierr); 277ac9112b8SAlp Dener PetscFunctionReturn(0); 278ac9112b8SAlp Dener } 279ac9112b8SAlp Dener 280ac9112b8SAlp Dener static PetscErrorCode TaoView_BNCG(Tao tao, PetscViewer viewer) 281ac9112b8SAlp Dener { 282ac9112b8SAlp Dener PetscBool isascii; 283ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 284ac9112b8SAlp Dener PetscErrorCode ierr; 285ac9112b8SAlp Dener 286ac9112b8SAlp Dener PetscFunctionBegin; 287ac9112b8SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 288ac9112b8SAlp Dener if (isascii) { 289ac9112b8SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 290ac9112b8SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cg->cg_type]);CHKERRQ(ierr); 291484c7b14SAdam Denchfield ierr = PetscViewerASCIIPrintf(viewer, "Skipped Stepdirection Updates: %i\n", cg->skipped_updates);CHKERRQ(ierr); 292484c7b14SAdam Denchfield ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %i\n", cg->resets);CHKERRQ(ierr); 293484c7b14SAdam Denchfield ierr = PetscViewerASCIIPrintf(viewer, "Pure gradient steps: %i\n", cg->pure_gd_steps);CHKERRQ(ierr); 294484c7b14SAdam Denchfield ierr = PetscViewerASCIIPrintf(viewer, "Not a descent direction: %i\n", cg->descent_error);CHKERRQ(ierr); 295ac9112b8SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Line search fails: %i\n", cg->ls_fails);CHKERRQ(ierr); 296484c7b14SAdam Denchfield if (cg->diag_scaling) { 297484c7b14SAdam Denchfield ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 298484c7b14SAdam Denchfield if (isascii) { 299484c7b14SAdam Denchfield ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 300484c7b14SAdam Denchfield ierr = MatView(cg->B, viewer);CHKERRQ(ierr); 301484c7b14SAdam Denchfield ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 302484c7b14SAdam Denchfield } 303484c7b14SAdam Denchfield } 304ac9112b8SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 305ac9112b8SAlp Dener } 306ac9112b8SAlp Dener PetscFunctionReturn(0); 307ac9112b8SAlp Dener } 308ac9112b8SAlp Dener 309c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal yty, PetscReal yts, PetscReal sts, PetscReal *scale, PetscReal alpha) 310c8bcdf1eSAdam Denchfield { 311c8bcdf1eSAdam Denchfield PetscReal a, b, c, sig1, sig2; 312c8bcdf1eSAdam Denchfield 313c8bcdf1eSAdam Denchfield PetscFunctionBegin; 314c8bcdf1eSAdam Denchfield *scale = 0.0; 3158ebe3e4eSStefano Zampini if (1.0 == alpha) *scale = yts/yty; 3168ebe3e4eSStefano Zampini else if (0.0 == alpha) *scale = sts/yts; 31750b47da0SAdam Denchfield else if (-1.0 == alpha) *scale = 1.0; 318c8bcdf1eSAdam Denchfield else { 319c8bcdf1eSAdam Denchfield a = yty; 320c8bcdf1eSAdam Denchfield b = yts; 321c8bcdf1eSAdam Denchfield c = sts; 322c8bcdf1eSAdam Denchfield a *= alpha; 323c8bcdf1eSAdam Denchfield b *= -(2.0*alpha - 1.0); 324c8bcdf1eSAdam Denchfield c *= alpha - 1.0; 325c8bcdf1eSAdam Denchfield sig1 = (-b + PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a); 326c8bcdf1eSAdam Denchfield sig2 = (-b - PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a); 327c8bcdf1eSAdam Denchfield /* accept the positive root as the scalar */ 3288ebe3e4eSStefano Zampini if (sig1 > 0.0) *scale = sig1; 3298ebe3e4eSStefano Zampini else if (sig2 > 0.0) *scale = sig2; 3308ebe3e4eSStefano Zampini else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Cannot find positive scalar"); 331c8bcdf1eSAdam Denchfield } 332c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 333c8bcdf1eSAdam Denchfield } 334c8bcdf1eSAdam Denchfield 335ac9112b8SAlp Dener /*MC 336ac9112b8SAlp Dener TAOBNCG - Bound-constrained Nonlinear Conjugate Gradient method. 337ac9112b8SAlp Dener 338ac9112b8SAlp Dener Options Database Keys: 33950b47da0SAdam Denchfield + -tao_bncg_recycle - enable recycling the latest calculated gradient vector in subsequent TaoSolve() calls (currently disabled) 340c4b75bccSAlp Dener . -tao_bncg_eta <r> - restart tolerance 34161be54a6SAlp Dener . -tao_bncg_type <taocg_type> - cg formula 342c4b75bccSAlp Dener . -tao_bncg_as_type <none,bertsekas> - active set estimation method 343c4b75bccSAlp Dener . -tao_bncg_as_tol <r> - tolerance used in Bertsekas active-set estimation 344c4b75bccSAlp Dener . -tao_bncg_as_step <r> - trial step length used in Bertsekas active-set estimation 34550b47da0SAdam 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. 34650b47da0SAdam Denchfield . -tao_bncg_theta <r> - convex combination parameter for the Broyden method 34750b47da0SAdam Denchfield . -tao_bncg_hz_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods 34850b47da0SAdam Denchfield . -tao_bncg_dk_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods 34950b47da0SAdam Denchfield . -tao_bncg_xi <r> - Multiplicative constant of the gamma term in the KD method 35050b47da0SAdam Denchfield . -tao_bncg_hz_theta <r> - Multiplicative constant of the theta term for the HZ method 35150b47da0SAdam Denchfield . -tao_bncg_bfgs_scale <r> - Scaling parameter of the bfgs contribution to the scalar Broyden method 35250b47da0SAdam Denchfield . -tao_bncg_dfp_scale <r> - Scaling parameter of the dfp contribution to the scalar Broyden method 35350b47da0SAdam Denchfield . -tao_bncg_diag_scaling <b> - Whether or not to use diagonal initialization/preconditioning for the CG methods. Default True. 35450b47da0SAdam Denchfield . -tao_bncg_dynamic_restart <b> - use dynamic restart strategy in the HZ, DK, KD methods 35550b47da0SAdam Denchfield . -tao_bncg_unscaled_restart <b> - whether or not to scale the gradient when doing gradient descent restarts 35650b47da0SAdam Denchfield . -tao_bncg_zeta <r> - Scaling parameter in the KD method 357484c7b14SAdam Denchfield . -tao_bncg_delta_min <r> - Minimum bound for rescaling during restarted gradient descent steps 358484c7b14SAdam Denchfield . -tao_bncg_delta_max <r> - Maximum bound for rescaling during restarted gradient descent steps 35950b47da0SAdam Denchfield . -tao_bncg_min_quad <i> - Number of quadratic-like steps in a row necessary to do a dynamic restart 36050b47da0SAdam 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 36150b47da0SAdam Denchfield . -tao_bncg_spaced_restart <b> - whether or not to do gradient descent steps every x*n iterations 362484c7b14SAdam Denchfield . -tao_bncg_no_scaling <b> - If true, eliminates all scaling, including defaults. 3633850be85SAlp Dener - -tao_bncg_neg_xi <b> - Whether or not to use negative xi in the KD method under certain conditions 364ac9112b8SAlp Dener 365ac9112b8SAlp Dener Notes: 366ac9112b8SAlp Dener CG formulas are: 3673850be85SAlp Dener + "gd" - Gradient Descent 3683850be85SAlp Dener . "fr" - Fletcher-Reeves 3693850be85SAlp Dener . "pr" - Polak-Ribiere-Polyak 3703850be85SAlp Dener . "prp" - Polak-Ribiere-Plus 3713850be85SAlp Dener . "hs" - Hestenes-Steifel 3723850be85SAlp Dener . "dy" - Dai-Yuan 3733850be85SAlp Dener . "ssml_bfgs" - Self-Scaling Memoryless BFGS 3743850be85SAlp Dener . "ssml_dfp" - Self-Scaling Memoryless DFP 3753850be85SAlp Dener . "ssml_brdn" - Self-Scaling Memoryless Broyden 3763850be85SAlp Dener . "hz" - Hager-Zhang (CG_DESCENT 5.3) 3773850be85SAlp Dener . "dk" - Dai-Kou (2013) 3783850be85SAlp Dener - "kd" - Kou-Dai (2015) 3799abc5736SPatrick Sanan 380ac9112b8SAlp Dener Level: beginner 381ac9112b8SAlp Dener M*/ 382ac9112b8SAlp Dener 383ac9112b8SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNCG(Tao tao) 384ac9112b8SAlp Dener { 385ac9112b8SAlp Dener TAO_BNCG *cg; 386ac9112b8SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 387ac9112b8SAlp Dener PetscErrorCode ierr; 388ac9112b8SAlp Dener 389ac9112b8SAlp Dener PetscFunctionBegin; 390ac9112b8SAlp Dener tao->ops->setup = TaoSetUp_BNCG; 391ac9112b8SAlp Dener tao->ops->solve = TaoSolve_BNCG; 392ac9112b8SAlp Dener tao->ops->view = TaoView_BNCG; 393ac9112b8SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNCG; 394ac9112b8SAlp Dener tao->ops->destroy = TaoDestroy_BNCG; 395ac9112b8SAlp Dener 396ac9112b8SAlp Dener /* Override default settings (unless already changed) */ 397ac9112b8SAlp Dener if (!tao->max_it_changed) tao->max_it = 2000; 398ac9112b8SAlp Dener if (!tao->max_funcs_changed) tao->max_funcs = 4000; 399ac9112b8SAlp Dener 400ac9112b8SAlp Dener /* Note: nondefault values should be used for nonlinear conjugate gradient */ 401ac9112b8SAlp Dener /* method. In particular, gtol should be less that 0.5; the value used in */ 402ac9112b8SAlp Dener /* Nocedal and Wright is 0.10. We use the default values for the */ 403ac9112b8SAlp Dener /* linesearch because it seems to work better. */ 404ac9112b8SAlp Dener ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 405ac9112b8SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 406ac9112b8SAlp Dener ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 407ac9112b8SAlp Dener ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr); 408ac9112b8SAlp Dener 409ac9112b8SAlp Dener ierr = PetscNewLog(tao,&cg);CHKERRQ(ierr); 410ac9112b8SAlp Dener tao->data = (void*)cg; 411484c7b14SAdam Denchfield ierr = KSPInitializePackage();CHKERRQ(ierr); 41250b47da0SAdam Denchfield ierr = MatCreate(PetscObjectComm((PetscObject)tao), &cg->B);CHKERRQ(ierr); 41350b47da0SAdam Denchfield ierr = PetscObjectIncrementTabLevel((PetscObject)cg->B, (PetscObject)tao, 1);CHKERRQ(ierr); 414864588a7SAlp Dener ierr = MatSetType(cg->B, MATLMVMDIAGBROYDEN);CHKERRQ(ierr); 41550b47da0SAdam Denchfield 416484c7b14SAdam Denchfield cg->pc = NULL; 417484c7b14SAdam Denchfield 41850b47da0SAdam Denchfield cg->dk_eta = 0.5; 41950b47da0SAdam Denchfield cg->hz_eta = 0.4; 420c8bcdf1eSAdam Denchfield cg->dynamic_restart = PETSC_FALSE; 421c8bcdf1eSAdam Denchfield cg->unscaled_restart = PETSC_FALSE; 422484c7b14SAdam Denchfield cg->no_scaling = PETSC_FALSE; 423484c7b14SAdam Denchfield cg->delta_min = 1e-7; 424484c7b14SAdam Denchfield cg->delta_max = 100; 425c8bcdf1eSAdam Denchfield cg->theta = 1.0; 426c8bcdf1eSAdam Denchfield cg->hz_theta = 1.0; 427c8bcdf1eSAdam Denchfield cg->dfp_scale = 1.0; 428c8bcdf1eSAdam Denchfield cg->bfgs_scale = 1.0; 42950b47da0SAdam Denchfield cg->zeta = 0.1; 43050b47da0SAdam Denchfield cg->min_quad = 6; 431c8bcdf1eSAdam Denchfield cg->min_restart_num = 6; /* As in CG_DESCENT and KD2015*/ 432c8bcdf1eSAdam Denchfield cg->xi = 1.0; 43350b47da0SAdam Denchfield cg->neg_xi = PETSC_TRUE; 434c8bcdf1eSAdam Denchfield cg->spaced_restart = PETSC_FALSE; 435c8bcdf1eSAdam Denchfield cg->tol_quad = 1e-8; 43661be54a6SAlp Dener cg->as_step = 0.001; 43761be54a6SAlp Dener cg->as_tol = 0.001; 43850b47da0SAdam Denchfield cg->eps_23 = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); /* Just a little tighter*/ 43961be54a6SAlp Dener cg->as_type = CG_AS_BERTSEKAS; 440c8bcdf1eSAdam Denchfield cg->cg_type = CG_SSML_BFGS; 441c8bcdf1eSAdam Denchfield cg->alpha = 1.0; 442c8bcdf1eSAdam Denchfield cg->diag_scaling = PETSC_TRUE; 443c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 444c8bcdf1eSAdam Denchfield } 445c8bcdf1eSAdam Denchfield 446c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGResetUpdate(Tao tao, PetscReal gnormsq) 447c8bcdf1eSAdam Denchfield { 448c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 449c8bcdf1eSAdam Denchfield PetscErrorCode ierr; 450c8bcdf1eSAdam Denchfield PetscReal scaling; 451c8bcdf1eSAdam Denchfield 452c8bcdf1eSAdam Denchfield PetscFunctionBegin; 453c8bcdf1eSAdam Denchfield ++cg->resets; 454c8bcdf1eSAdam Denchfield scaling = 2.0 * PetscMax(1.0, PetscAbsScalar(cg->f)) / PetscMax(gnormsq, cg->eps_23); 455484c7b14SAdam Denchfield scaling = PetscMin(cg->delta_max, PetscMax(cg->delta_min, scaling)); 456484c7b14SAdam Denchfield if (cg->unscaled_restart) { 457484c7b14SAdam Denchfield scaling = 1.0; 458484c7b14SAdam Denchfield ++cg->pure_gd_steps; 459484c7b14SAdam Denchfield } 460c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -scaling, 0.0, tao->gradient);CHKERRQ(ierr); 461c8bcdf1eSAdam Denchfield /* Also want to reset our diagonal scaling with each restart */ 462c8bcdf1eSAdam Denchfield if (cg->diag_scaling) { 46350b47da0SAdam Denchfield ierr = MatLMVMReset(cg->B, PETSC_FALSE);CHKERRQ(ierr); 464c8bcdf1eSAdam Denchfield } 465c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 466c8bcdf1eSAdam Denchfield } 467c8bcdf1eSAdam Denchfield 468c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGCheckDynamicRestart(Tao tao, PetscReal stepsize, PetscReal gd, PetscReal gd_old, PetscBool *dynrestart, PetscReal fold) 469c8bcdf1eSAdam Denchfield { 470c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 471c8bcdf1eSAdam Denchfield PetscReal quadinterp; 472c8bcdf1eSAdam Denchfield 473c8bcdf1eSAdam Denchfield PetscFunctionBegin; 47450b47da0SAdam Denchfield if (cg->f < cg->min_quad/10) { 47550b47da0SAdam Denchfield *dynrestart = PETSC_FALSE; 47650b47da0SAdam Denchfield PetscFunctionReturn(0); 47750b47da0SAdam Denchfield } /* just skip this since this strategy doesn't work well for functions near zero */ 478484c7b14SAdam Denchfield quadinterp = 2.0*(cg->f - fold)/(stepsize*(gd + gd_old)); 47950b47da0SAdam Denchfield if (PetscAbs(quadinterp - 1.0) < cg->tol_quad) ++cg->iter_quad; 480c8bcdf1eSAdam Denchfield else { 481c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 482c8bcdf1eSAdam Denchfield *dynrestart = PETSC_FALSE; 483c8bcdf1eSAdam Denchfield } 484c8bcdf1eSAdam Denchfield if (cg->iter_quad >= cg->min_quad) { 485c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 486c8bcdf1eSAdam Denchfield *dynrestart = PETSC_TRUE; 487c8bcdf1eSAdam Denchfield } 488c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 489c8bcdf1eSAdam Denchfield } 490c8bcdf1eSAdam Denchfield 4918ca2df50S PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao tao, PetscReal gnorm2, PetscReal step, PetscReal fold, PetscReal gnorm2_old, PetscReal dnorm, PetscBool pcgd_fallback) 49250b47da0SAdam Denchfield { 493c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 494c8bcdf1eSAdam Denchfield PetscErrorCode ierr; 49550b47da0SAdam Denchfield PetscReal gamma = 1.0, tau_k, beta; 496484c7b14SAdam Denchfield PetscReal tmp = 1.0, ynorm, ynorm2 = 1.0, snorm = 1.0, dk_yk=1.0, gd; 49750b47da0SAdam Denchfield PetscReal gkp1_yk, gd_old, tau_bfgs, tau_dfp, gkp1D_yk, gtDg; 498c8bcdf1eSAdam Denchfield PetscInt dim; 499484c7b14SAdam Denchfield PetscBool cg_restart = PETSC_FALSE; 500c8bcdf1eSAdam Denchfield PetscFunctionBegin; 501c8bcdf1eSAdam Denchfield 50250b47da0SAdam Denchfield /* Local curvature check to see if we need to restart */ 503414d97d3SAlp Dener if (tao->niter >= 1 || tao->recycle) { 504c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr); 505c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr); 506c8bcdf1eSAdam Denchfield ynorm2 = ynorm*ynorm; 507c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk);CHKERRQ(ierr); 508484c7b14SAdam Denchfield if (step*dnorm < PETSC_MACHINE_EPSILON || step*dk_yk < PETSC_MACHINE_EPSILON) { 509e2570530SAlp Dener cg_restart = PETSC_TRUE; 510484c7b14SAdam Denchfield ++cg->skipped_updates; 511484c7b14SAdam Denchfield } 51250b47da0SAdam Denchfield if (cg->spaced_restart) { 51350b47da0SAdam Denchfield ierr = VecGetSize(tao->gradient, &dim);CHKERRQ(ierr); 514e2570530SAlp Dener if (tao->niter % (dim*cg->min_restart_num)) cg_restart = PETSC_TRUE; 51550b47da0SAdam Denchfield } 51650b47da0SAdam Denchfield } 51750b47da0SAdam Denchfield /* If the user wants regular restarts, do it every 6n iterations, where n=dimension */ 51850b47da0SAdam Denchfield if (cg->spaced_restart) { 51950b47da0SAdam Denchfield ierr = VecGetSize(tao->gradient, &dim);CHKERRQ(ierr); 520e2570530SAlp Dener if (0 == tao->niter % (6*dim)) cg_restart = PETSC_TRUE; 52150b47da0SAdam Denchfield } 52250b47da0SAdam Denchfield /* Compute the diagonal scaling vector if applicable */ 52350b47da0SAdam Denchfield if (cg->diag_scaling) { 52450b47da0SAdam Denchfield ierr = MatLMVMUpdate(cg->B, tao->solution, tao->gradient);CHKERRQ(ierr); 52550b47da0SAdam Denchfield } 52650b47da0SAdam Denchfield 527484c7b14SAdam Denchfield /* A note on diagonal scaling (to be added to paper): 528484c7b14SAdam Denchfield For the FR, PR, PRP, and DY methods, the diagonally scaled versions 529484c7b14SAdam Denchfield must be derived as a preconditioned CG method rather than as 530484c7b14SAdam Denchfield a Hessian initialization like in the Broyden methods. */ 53150b47da0SAdam Denchfield 532484c7b14SAdam Denchfield /* In that case, one writes the objective function as 533484c7b14SAdam Denchfield f(x) \equiv f(Ay). Gradient evaluations yield g(x_k) = A g(Ay_k) = A g(x_k). 534484c7b14SAdam Denchfield Furthermore, the direction d_k \equiv (x_k - x_{k-1})/step according to 535484c7b14SAdam Denchfield HZ (2006) becomes A^{-1} d_k, such that d_k^T g_k remains the 536484c7b14SAdam Denchfield same under preconditioning. Note that A is diagonal, such that A^T = A. */ 53750b47da0SAdam Denchfield 538484c7b14SAdam Denchfield /* This yields questions like what the dot product d_k^T y_k 539484c7b14SAdam Denchfield should look like. HZ mistakenly treats that as the same under 540484c7b14SAdam Denchfield preconditioning, but that is not necessarily true. */ 54150b47da0SAdam Denchfield 542484c7b14SAdam Denchfield /* Observe y_k \equiv g_k - g_{k-1}, and under the P.C. transformation, 543484c7b14SAdam 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}), 544484c7b14SAdam 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 545484c7b14SAdam Denchfield NOT the same if our preconditioning matrix is updated between iterations. 546484c7b14SAdam Denchfield This same issue is found when considering dot products of the form g_{k+1}^T y_k. */ 54750b47da0SAdam Denchfield 54850b47da0SAdam Denchfield /* Compute CG step direction */ 54950b47da0SAdam Denchfield if (cg_restart) { 55050b47da0SAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 551484c7b14SAdam Denchfield } else if (pcgd_fallback) { 552484c7b14SAdam Denchfield /* Just like preconditioned CG */ 553484c7b14SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 554484c7b14SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work);CHKERRQ(ierr); 55550b47da0SAdam Denchfield } else if (ynorm2 > PETSC_MACHINE_EPSILON) { 55650b47da0SAdam Denchfield switch (cg->cg_type) { 557484c7b14SAdam Denchfield case CG_PCGradientDescent: 55850b47da0SAdam Denchfield if (!cg->diag_scaling) { 559484c7b14SAdam Denchfield if (!cg->no_scaling) { 56050b47da0SAdam Denchfield cg->sts = step*step*dnorm*dnorm; 56150b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 562484c7b14SAdam Denchfield } else { 563484c7b14SAdam Denchfield tau_k = 1.0; 564484c7b14SAdam Denchfield ++cg->pure_gd_steps; 565484c7b14SAdam Denchfield } 56650b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, 0.0, tao->gradient);CHKERRQ(ierr); 56750b47da0SAdam Denchfield } else { 56850b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 56950b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work);CHKERRQ(ierr); 57050b47da0SAdam Denchfield } 57150b47da0SAdam Denchfield break; 572484c7b14SAdam Denchfield 57350b47da0SAdam Denchfield case CG_HestenesStiefel: 57450b47da0SAdam Denchfield /* Classic Hestenes-Stiefel method, modified with scalar and diagonal preconditioning. */ 57550b47da0SAdam Denchfield if (!cg->diag_scaling) { 57650b47da0SAdam Denchfield cg->sts = step*step*dnorm*dnorm; 57750b47da0SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 57850b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 57950b47da0SAdam Denchfield beta = tau_k*gkp1_yk/dk_yk; 58050b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 58150b47da0SAdam Denchfield } else { 58250b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 58350b47da0SAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr); 58450b47da0SAdam Denchfield beta = gkp1_yk/dk_yk; 58550b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 58650b47da0SAdam Denchfield } 587c8bcdf1eSAdam Denchfield break; 588484c7b14SAdam Denchfield 589c8bcdf1eSAdam Denchfield case CG_FletcherReeves: 59050b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr); 59150b47da0SAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr); 59250b47da0SAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr); 59350b47da0SAdam Denchfield ynorm2 = ynorm*ynorm; 59450b47da0SAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk);CHKERRQ(ierr); 59550b47da0SAdam Denchfield if (!cg->diag_scaling) { 59650b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, step*step*dnorm*dnorm, &tau_k, cg->alpha);CHKERRQ(ierr); 59750b47da0SAdam Denchfield beta = tau_k*gnorm2/gnorm2_old; 59850b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 59950b47da0SAdam Denchfield } else { 60050b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->g_work, &gnorm2_old);CHKERRQ(ierr); /* Before it's updated */ 60150b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 60250b47da0SAdam Denchfield ierr = VecDot(tao->gradient, cg->g_work, &tmp);CHKERRQ(ierr); 60350b47da0SAdam Denchfield beta = tmp/gnorm2_old; 60450b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 60550b47da0SAdam Denchfield } 606c8bcdf1eSAdam Denchfield break; 607484c7b14SAdam Denchfield 60850b47da0SAdam Denchfield case CG_PolakRibierePolyak: 60950b47da0SAdam Denchfield snorm = step*dnorm; 61050b47da0SAdam Denchfield if (!cg->diag_scaling) { 61150b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr); 61250b47da0SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 61350b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 61450b47da0SAdam Denchfield beta = tau_k*gkp1_yk/gnorm2_old; 61550b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 61650b47da0SAdam Denchfield } else { 61750b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->g_work, &gnorm2_old);CHKERRQ(ierr); 61850b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 61950b47da0SAdam Denchfield ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 62050b47da0SAdam Denchfield beta = gkp1_yk/gnorm2_old; 62150b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 62250b47da0SAdam Denchfield } 623c8bcdf1eSAdam Denchfield break; 624484c7b14SAdam Denchfield 625c8bcdf1eSAdam Denchfield case CG_PolakRibierePlus: 62650b47da0SAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr); 62750b47da0SAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr); 62850b47da0SAdam Denchfield ynorm2 = ynorm*ynorm; 62950b47da0SAdam Denchfield if (!cg->diag_scaling) { 63050b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr); 63150b47da0SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 63250b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 63350b47da0SAdam Denchfield beta = tau_k*gkp1_yk/gnorm2_old; 63450b47da0SAdam Denchfield beta = PetscMax(beta, 0.0); 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); /* Old gtDg */ 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 beta = PetscMax(beta, 0.0); 64250b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 64350b47da0SAdam Denchfield } 644c8bcdf1eSAdam Denchfield break; 645484c7b14SAdam Denchfield 646484c7b14SAdam Denchfield case CG_DaiYuan: 647484c7b14SAdam Denchfield /* Dai, Yu-Hong, and Yaxiang Yuan. "A nonlinear conjugate gradient method with a strong global convergence property." 648484c7b14SAdam Denchfield SIAM Journal on optimization 10, no. 1 (1999): 177-182. */ 64950b47da0SAdam Denchfield if (!cg->diag_scaling) { 65050b47da0SAdam Denchfield ierr = VecDot(tao->stepdirection, tao->gradient, &gd);CHKERRQ(ierr); 651c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 65250b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->yts, &tau_k, cg->alpha);CHKERRQ(ierr); 65350b47da0SAdam Denchfield beta = tau_k*gnorm2/(gd - gd_old); 65450b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 65550b47da0SAdam Denchfield } else { 656484c7b14SAdam Denchfield ierr = MatMult(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr); 657484c7b14SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 65850b47da0SAdam Denchfield ierr = VecDot(cg->g_work, tao->gradient, >Dg);CHKERRQ(ierr); 65950b47da0SAdam Denchfield ierr = VecDot(tao->stepdirection, cg->G_old, &gd_old);CHKERRQ(ierr); 66050b47da0SAdam Denchfield ierr = VecDot(cg->d_work, cg->g_work, &dk_yk);CHKERRQ(ierr); 66150b47da0SAdam Denchfield dk_yk = dk_yk - gd_old; 66250b47da0SAdam Denchfield beta = gtDg/dk_yk; 663c624ebd3SAlp Dener ierr = VecScale(cg->d_work, beta);CHKERRQ(ierr); 66450b47da0SAdam Denchfield ierr = VecWAXPY(tao->stepdirection, -1.0, cg->g_work, cg->d_work);CHKERRQ(ierr); 66550b47da0SAdam Denchfield } 666c8bcdf1eSAdam Denchfield break; 667484c7b14SAdam Denchfield 668c8bcdf1eSAdam Denchfield case CG_HagerZhang: 669484c7b14SAdam Denchfield /* Hager, William W., and Hongchao Zhang. "Algorithm 851: CG_DESCENT, a conjugate gradient method with guaranteed descent." 670484c7b14SAdam Denchfield ACM Transactions on Mathematical Software (TOMS) 32, no. 1 (2006): 113-137. */ 671c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 672c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 673c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 67450b47da0SAdam Denchfield snorm = dnorm*step; 67550b47da0SAdam Denchfield cg->yts = step*dk_yk; 676c8bcdf1eSAdam Denchfield if (cg->use_dynamic_restart) { 677c8bcdf1eSAdam Denchfield ierr = TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold);CHKERRQ(ierr); 678c8bcdf1eSAdam Denchfield } 67950b47da0SAdam Denchfield if (cg->dynamic_restart) { 680c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 681c8bcdf1eSAdam Denchfield } else { 682c8bcdf1eSAdam Denchfield if (!cg->diag_scaling) { 683c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 684c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 685c8bcdf1eSAdam Denchfield /* Supplying cg->alpha = -1.0 will give the CG_DESCENT 5.3 special case of tau_k = 1.0 */ 686c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 687c8bcdf1eSAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk)); 688c8bcdf1eSAdam Denchfield /* Bound beta as in CG_DESCENT 5.3, as implemented, with the third comparison from DK 2013 */ 68950b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta*tau_k*gd_old/(dnorm*dnorm)), cg->dk_eta*tau_k*gd/(dnorm*dnorm)); 690c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d */ 69150b47da0SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 692c8bcdf1eSAdam Denchfield } else { 693c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 694c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 695c8bcdf1eSAdam Denchfield cg->sts = snorm*snorm; 69650b47da0SAdam Denchfield /* Apply the diagonal scaling to all my vectors */ 69750b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 69850b47da0SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 69950b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr); 700c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 701c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr); 702c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */ 703c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 70450b47da0SAdam Denchfield ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr); 705c8bcdf1eSAdam Denchfield tau_k = -tau_k*gd/(dk_yk*dk_yk); 706c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */ 707484c7b14SAdam Denchfield beta = gkp1_yk/dk_yk + cg->hz_theta*tau_k; /* HZ; (1.15) from DK 2013 */ 708c8bcdf1eSAdam Denchfield /* From HZ2013, modified to account for diagonal scaling*/ 70950b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->d_work, &gd_old);CHKERRQ(ierr); 71050b47da0SAdam Denchfield ierr = VecDot(tao->stepdirection, cg->g_work, &gd);CHKERRQ(ierr); 71150b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta*gd_old/(dnorm*dnorm)), cg->dk_eta*gd/(dnorm*dnorm)); 712c8bcdf1eSAdam Denchfield /* Do the update */ 713484c7b14SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 71450b47da0SAdam Denchfield } 71550b47da0SAdam Denchfield } 716c8bcdf1eSAdam Denchfield break; 717484c7b14SAdam Denchfield 718c8bcdf1eSAdam Denchfield case CG_DaiKou: 719484c7b14SAdam Denchfield /* Dai, Yu-Hong, and Cai-Xia Kou. "A nonlinear conjugate gradient algorithm with an optimal property and an improved Wolfe line search." 720484c7b14SAdam Denchfield SIAM Journal on Optimization 23, no. 1 (2013): 296-320. */ 721c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 722c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 723c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 72450b47da0SAdam Denchfield snorm = step*dnorm; 72550b47da0SAdam Denchfield cg->yts = dk_yk*step; 726c8bcdf1eSAdam Denchfield if (!cg->diag_scaling) { 727c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 72850b47da0SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 729c8bcdf1eSAdam Denchfield /* Use cg->alpha = -1.0 to get tau_k = 1.0 as in CG_DESCENT 5.3 */ 730c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 73150b47da0SAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk) + gd/(dnorm*dnorm)) - step*gd/dk_yk; 73250b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta*tau_k*gd_old/(dnorm*dnorm)), cg->dk_eta*tau_k*gd/(dnorm*dnorm)); 733c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d */ 734c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, 0.0, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 735c8bcdf1eSAdam Denchfield } else { 736c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 737c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 738c8bcdf1eSAdam Denchfield cg->sts = snorm*snorm; 73950b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 74050b47da0SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 74150b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr); 742c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 743c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr); 74450b47da0SAdam Denchfield ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr); 745c8bcdf1eSAdam Denchfield tau_k = tau_k*gd/(dk_yk*dk_yk); 746c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 747c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */ 748484c7b14SAdam Denchfield beta = gkp1_yk/dk_yk - step*tmp - tau_k; 749c8bcdf1eSAdam Denchfield /* Update this for the last term in beta */ 750c8bcdf1eSAdam Denchfield ierr = VecDot(cg->y_work, tao->stepdirection, &dk_yk);CHKERRQ(ierr); 751c8bcdf1eSAdam Denchfield beta += tmp*dk_yk/(dnorm*dnorm); /* projection of y_work onto dk */ 75250b47da0SAdam Denchfield ierr = VecDot(tao->stepdirection, cg->g_work, &gd);CHKERRQ(ierr); 75350b47da0SAdam Denchfield ierr = VecDot(cg->G_old, cg->d_work, &gd_old);CHKERRQ(ierr); 75450b47da0SAdam Denchfield beta = PetscMax(PetscMax(beta, cg->hz_eta*gd_old/(dnorm*dnorm)), cg->dk_eta*gd/(dnorm*dnorm)); 755c8bcdf1eSAdam Denchfield /* Do the update */ 756484c7b14SAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 75750b47da0SAdam Denchfield } 758c8bcdf1eSAdam Denchfield break; 759484c7b14SAdam Denchfield 760c8bcdf1eSAdam Denchfield case CG_KouDai: 761110fc3b0SBarry Smith /* Kou, Cai-Xia, and Yu-Hong Dai. "A modified self-scaling memoryless Broyden-Fletcher-Goldfarb-Shanno method for unconstrained optimization." 762484c7b14SAdam Denchfield Journal of Optimization Theory and Applications 165, no. 1 (2015): 209-224. */ 763c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 764c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 765c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 76650b47da0SAdam Denchfield snorm = step*dnorm; 76750b47da0SAdam Denchfield cg->yts = dk_yk*step; 768c8bcdf1eSAdam Denchfield if (cg->use_dynamic_restart) { 769c8bcdf1eSAdam Denchfield ierr = TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold);CHKERRQ(ierr); 770c8bcdf1eSAdam Denchfield } 77150b47da0SAdam Denchfield if (cg->dynamic_restart) { 772c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 773c8bcdf1eSAdam Denchfield } else { 774c8bcdf1eSAdam Denchfield if (!cg->diag_scaling) { 775c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 776c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 777c8bcdf1eSAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk)) - step*gd/dk_yk; 778c8bcdf1eSAdam Denchfield if (beta < cg->zeta*tau_k*gd/(dnorm*dnorm)) /* 0.1 is KD's zeta parameter */ 779c8bcdf1eSAdam Denchfield { 780c8bcdf1eSAdam Denchfield beta = cg->zeta*tau_k*gd/(dnorm*dnorm); 781c8bcdf1eSAdam Denchfield gamma = 0.0; 782c8bcdf1eSAdam Denchfield } else { 783c8bcdf1eSAdam Denchfield if (gkp1_yk < 0 && cg->neg_xi) gamma = -1.0*gd/dk_yk; 784484c7b14SAdam Denchfield /* This seems to be very effective when there's no tau_k scaling. 785484c7b14SAdam Denchfield This guarantees a large descent step every iteration, going through DK 2015 Lemma 3.1's proof but allowing for negative xi */ 78650b47da0SAdam Denchfield else { 78750b47da0SAdam Denchfield gamma = cg->xi*gd/dk_yk; 78850b47da0SAdam Denchfield } 789c8bcdf1eSAdam Denchfield } 790c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d + t*tmp*yk */ 791c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, gamma*tau_k, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 792c8bcdf1eSAdam Denchfield } else { 793c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 794c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 795c8bcdf1eSAdam Denchfield cg->sts = snorm*snorm; 79650b47da0SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 79750b47da0SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 798c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 799c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1D_yk);CHKERRQ(ierr); 800c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */ 801c8bcdf1eSAdam Denchfield gamma = gd/dk_yk; 802c8bcdf1eSAdam Denchfield /* tau_k = -ytDy/(ytd)^2 * gd */ 80350b47da0SAdam Denchfield ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr); 804c8bcdf1eSAdam Denchfield tau_k = tau_k*gd/(dk_yk*dk_yk); 805c8bcdf1eSAdam Denchfield /* beta is the constant which adds the d_k contribution */ 806c8bcdf1eSAdam Denchfield beta = gkp1D_yk/dk_yk - step*gamma - tau_k; 807c8bcdf1eSAdam Denchfield /* Here is the requisite check */ 80850b47da0SAdam Denchfield ierr = VecDot(tao->stepdirection, cg->g_work, &tmp);CHKERRQ(ierr); 809c8bcdf1eSAdam Denchfield if (cg->neg_xi) { 810c8bcdf1eSAdam Denchfield /* modified KD implementation */ 811c8bcdf1eSAdam Denchfield if (gkp1D_yk/dk_yk < 0) gamma = -1.0*gd/dk_yk; 81250b47da0SAdam Denchfield else { 81350b47da0SAdam Denchfield gamma = cg->xi*gd/dk_yk; 81450b47da0SAdam Denchfield } 815c8bcdf1eSAdam Denchfield if (beta < cg->zeta*tmp/(dnorm*dnorm)) { 816c8bcdf1eSAdam Denchfield beta = cg->zeta*tmp/(dnorm*dnorm); 817c8bcdf1eSAdam Denchfield gamma = 0.0; 818c8bcdf1eSAdam Denchfield } 819c8bcdf1eSAdam Denchfield } else { /* original KD 2015 implementation */ 820c8bcdf1eSAdam Denchfield if (beta < cg->zeta*tmp/(dnorm*dnorm)) { 821c8bcdf1eSAdam Denchfield beta = cg->zeta*tmp/(dnorm*dnorm); 822c8bcdf1eSAdam Denchfield gamma = 0.0; 823c8bcdf1eSAdam Denchfield } else { 824c8bcdf1eSAdam Denchfield gamma = cg->xi*gd/dk_yk; 825c8bcdf1eSAdam Denchfield } 826c8bcdf1eSAdam Denchfield } 827c8bcdf1eSAdam Denchfield /* Do the update in two steps */ 828c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 829c8bcdf1eSAdam Denchfield ierr = VecAXPY(tao->stepdirection, gamma, cg->y_work);CHKERRQ(ierr); 83050b47da0SAdam Denchfield } 83150b47da0SAdam Denchfield } 832c8bcdf1eSAdam Denchfield break; 833484c7b14SAdam Denchfield 834484c7b14SAdam Denchfield case CG_SSML_BFGS: 835484c7b14SAdam Denchfield /* Perry, J. M. "A class of conjugate gradient algorithms with a two-step variable-metric memory." 836484c7b14SAdam Denchfield Discussion Papers 269 (1977). */ 837484c7b14SAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 838484c7b14SAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 839484c7b14SAdam Denchfield snorm = step*dnorm; 840484c7b14SAdam Denchfield cg->yts = dk_yk*step; 841484c7b14SAdam Denchfield cg->yty = ynorm2; 842484c7b14SAdam Denchfield cg->sts = snorm*snorm; 843484c7b14SAdam Denchfield if (!cg->diag_scaling) { 844484c7b14SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 845484c7b14SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 846484c7b14SAdam Denchfield tmp = gd/dk_yk; 847484c7b14SAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - cg->yty*gd/(dk_yk*dk_yk)) - step*tmp; 848484c7b14SAdam Denchfield /* d <- -t*g + beta*t*d + t*tmp*yk */ 849484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp*tau_k, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 850484c7b14SAdam Denchfield } else { 851484c7b14SAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 852484c7b14SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 853484c7b14SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 854484c7b14SAdam Denchfield /* compute scalar gamma */ 855484c7b14SAdam Denchfield ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 856484c7b14SAdam Denchfield ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr); 857484c7b14SAdam Denchfield gamma = gd/dk_yk; 858484c7b14SAdam Denchfield /* Compute scalar beta */ 859484c7b14SAdam Denchfield beta = (gkp1_yk/dk_yk - gd*tmp/(dk_yk*dk_yk)) - step*gd/dk_yk; 860484c7b14SAdam Denchfield /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 861484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr); 862484c7b14SAdam Denchfield } 863484c7b14SAdam Denchfield break; 864484c7b14SAdam Denchfield 865484c7b14SAdam Denchfield case CG_SSML_DFP: 866484c7b14SAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 867484c7b14SAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 868484c7b14SAdam Denchfield snorm = step*dnorm; 869484c7b14SAdam Denchfield cg->yts = dk_yk*step; 870484c7b14SAdam Denchfield cg->yty = ynorm2; 871484c7b14SAdam Denchfield cg->sts = snorm*snorm; 872484c7b14SAdam Denchfield if (!cg->diag_scaling) { 873484c7b14SAdam Denchfield /* Instead of a regular convex combination, we will solve a quadratic formula. */ 874484c7b14SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 875484c7b14SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 876484c7b14SAdam Denchfield tau_k = cg->dfp_scale*tau_k; 877484c7b14SAdam Denchfield tmp = tau_k*gkp1_yk/cg->yty; 878484c7b14SAdam Denchfield beta = -step*gd/dk_yk; 879484c7b14SAdam Denchfield /* d <- -t*g + beta*d + tmp*yk */ 880484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 881484c7b14SAdam Denchfield } else { 882484c7b14SAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless DFP step */ 883484c7b14SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 884484c7b14SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 885484c7b14SAdam Denchfield /* compute scalar gamma */ 886484c7b14SAdam Denchfield ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 887484c7b14SAdam Denchfield ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr); 888484c7b14SAdam Denchfield gamma = (gkp1_yk/tmp); 889484c7b14SAdam Denchfield /* Compute scalar beta */ 890484c7b14SAdam Denchfield beta = -step*gd/dk_yk; 891484c7b14SAdam Denchfield /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 892484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr); 893484c7b14SAdam Denchfield } 894484c7b14SAdam Denchfield break; 895484c7b14SAdam Denchfield 896484c7b14SAdam Denchfield case CG_SSML_BROYDEN: 897484c7b14SAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 898484c7b14SAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 899484c7b14SAdam Denchfield snorm = step*dnorm; 900484c7b14SAdam Denchfield cg->yts = step*dk_yk; 901484c7b14SAdam Denchfield cg->yty = ynorm2; 902484c7b14SAdam Denchfield cg->sts = snorm*snorm; 903484c7b14SAdam Denchfield if (!cg->diag_scaling) { 904484c7b14SAdam Denchfield /* Instead of a regular convex combination, we will solve a quadratic formula. */ 905484c7b14SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(cg->yty, step*dk_yk, snorm*snorm, &tau_bfgs, cg->bfgs_scale);CHKERRQ(ierr); 906484c7b14SAdam Denchfield ierr = TaoBNCGComputeScalarScaling(cg->yty, step*dk_yk, snorm*snorm, &tau_dfp, cg->dfp_scale);CHKERRQ(ierr); 907484c7b14SAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 908484c7b14SAdam Denchfield tau_k = cg->theta*tau_bfgs + (1.0-cg->theta)*tau_dfp; 909484c7b14SAdam Denchfield /* If bfgs_scale = 1.0, it should reproduce the bfgs tau_bfgs. If bfgs_scale = 0.0, 910484c7b14SAdam Denchfield it should reproduce the tau_dfp scaling. Same with dfp_scale. */ 911484c7b14SAdam Denchfield tmp = cg->theta*tau_bfgs*gd/dk_yk + (1-cg->theta)*tau_dfp*gkp1_yk/cg->yty; 912484c7b14SAdam Denchfield beta = cg->theta*tau_bfgs*(gkp1_yk/dk_yk - cg->yty*gd/(dk_yk*dk_yk)) - step*gd/dk_yk; 913484c7b14SAdam Denchfield /* d <- -t*g + beta*d + tmp*yk */ 914484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 915484c7b14SAdam Denchfield } else { 916484c7b14SAdam Denchfield /* We have diagonal scaling enabled */ 917484c7b14SAdam Denchfield ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 918484c7b14SAdam Denchfield ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 919484c7b14SAdam Denchfield /* compute scalar gamma */ 920484c7b14SAdam Denchfield ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 921484c7b14SAdam Denchfield ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr); 922484c7b14SAdam Denchfield gamma = cg->theta*gd/dk_yk + (1-cg->theta)*(gkp1_yk/tmp); 923484c7b14SAdam Denchfield /* Compute scalar beta */ 924484c7b14SAdam Denchfield beta = cg->theta*(gkp1_yk/dk_yk - gd*tmp/(dk_yk*dk_yk)) - step*gd/dk_yk; 925484c7b14SAdam Denchfield /* Compute stepdirection dkp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 926484c7b14SAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr); 927484c7b14SAdam Denchfield } 928484c7b14SAdam Denchfield break; 929484c7b14SAdam Denchfield 930c8bcdf1eSAdam Denchfield default: 931c8bcdf1eSAdam Denchfield break; 932484c7b14SAdam Denchfield 933c8bcdf1eSAdam Denchfield } 934c8bcdf1eSAdam Denchfield } 935c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 936c8bcdf1eSAdam Denchfield } 937c8bcdf1eSAdam Denchfield 938c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao tao, PetscReal gnorm) 939c8bcdf1eSAdam Denchfield { 940c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 941c8bcdf1eSAdam Denchfield PetscErrorCode ierr; 942c8bcdf1eSAdam Denchfield TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 9438ca2df50S PetscReal step=1.0,gnorm2,gd,dnorm=0.0; 944c8bcdf1eSAdam Denchfield PetscReal gnorm2_old,f_old,resnorm, gnorm_old; 945c624ebd3SAlp Dener PetscBool pcgd_fallback = PETSC_FALSE; 946c8bcdf1eSAdam Denchfield 947c8bcdf1eSAdam Denchfield PetscFunctionBegin; 948c8bcdf1eSAdam Denchfield /* We are now going to perform a line search along the direction. */ 949c8bcdf1eSAdam Denchfield /* Store solution and gradient info before it changes */ 950c8bcdf1eSAdam Denchfield ierr = VecCopy(tao->solution, cg->X_old);CHKERRQ(ierr); 951c8bcdf1eSAdam Denchfield ierr = VecCopy(tao->gradient, cg->G_old);CHKERRQ(ierr); 952c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->unprojected_gradient, cg->unprojected_gradient_old);CHKERRQ(ierr); 953c8bcdf1eSAdam Denchfield 954c8bcdf1eSAdam Denchfield gnorm_old = gnorm; 955c8bcdf1eSAdam Denchfield gnorm2_old = gnorm_old*gnorm_old; 956c8bcdf1eSAdam Denchfield f_old = cg->f; 957484c7b14SAdam Denchfield /* Perform bounded line search. If we are recycling a solution from a previous */ 958484c7b14SAdam Denchfield /* TaoSolve, then we want to immediately skip to calculating a new direction rather than performing a linesearch */ 959414d97d3SAlp Dener if (!(tao->recycle && 0 == tao->niter)) { 960484c7b14SAdam Denchfield /* Above logic: the below code happens every iteration, except for the first iteration of a recycled TaoSolve */ 961c8bcdf1eSAdam Denchfield ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 962c8bcdf1eSAdam Denchfield ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 963c8bcdf1eSAdam Denchfield ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 964c8bcdf1eSAdam Denchfield 965c8bcdf1eSAdam Denchfield /* Check linesearch failure */ 966c8bcdf1eSAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 967c8bcdf1eSAdam Denchfield ++cg->ls_fails; 968c624ebd3SAlp Dener if (cg->cg_type == CG_GradientDescent) { 969c8bcdf1eSAdam Denchfield /* Nothing left to do but fail out of the optimization */ 970c8bcdf1eSAdam Denchfield step = 0.0; 971c8bcdf1eSAdam Denchfield tao->reason = TAO_DIVERGED_LS_FAILURE; 972c8bcdf1eSAdam Denchfield } else { 973484c7b14SAdam Denchfield /* Restore previous point, perform preconditioned GD and regular GD steps at the last good point */ 974c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->X_old, tao->solution);CHKERRQ(ierr); 975c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->G_old, tao->gradient);CHKERRQ(ierr); 976c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->unprojected_gradient_old, cg->unprojected_gradient);CHKERRQ(ierr); 977c8bcdf1eSAdam Denchfield gnorm = gnorm_old; 978c8bcdf1eSAdam Denchfield gnorm2 = gnorm2_old; 979c8bcdf1eSAdam Denchfield cg->f = f_old; 980c8bcdf1eSAdam Denchfield 981484c7b14SAdam Denchfield /* Fall back on preconditioned CG (so long as you're not already using it) */ 982484c7b14SAdam Denchfield if (cg->cg_type != CG_PCGradientDescent && cg->diag_scaling) { 983e2570530SAlp Dener pcgd_fallback = PETSC_TRUE; 9848ca2df50S ierr = TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback);CHKERRQ(ierr); 985484c7b14SAdam Denchfield 98650b47da0SAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 987c8bcdf1eSAdam Denchfield ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 988c8bcdf1eSAdam Denchfield 989c8bcdf1eSAdam Denchfield ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 990c8bcdf1eSAdam Denchfield ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 991c8bcdf1eSAdam Denchfield ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 992c8bcdf1eSAdam Denchfield 993484c7b14SAdam Denchfield pcgd_fallback = PETSC_FALSE; 994484c7b14SAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 995484c7b14SAdam Denchfield /* Going to perform a regular gradient descent step. */ 996484c7b14SAdam Denchfield ++cg->ls_fails; 997484c7b14SAdam Denchfield step = 0.0; 998484c7b14SAdam Denchfield } 999484c7b14SAdam Denchfield } 1000484c7b14SAdam Denchfield /* Fall back on the scaled gradient step */ 1001484c7b14SAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 1002484c7b14SAdam Denchfield ++cg->ls_fails; 1003484c7b14SAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 1004484c7b14SAdam Denchfield ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1005484c7b14SAdam Denchfield ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 1006484c7b14SAdam Denchfield ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 1007484c7b14SAdam Denchfield ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 1008484c7b14SAdam Denchfield } 1009484c7b14SAdam Denchfield 1010c8bcdf1eSAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 1011c8bcdf1eSAdam Denchfield /* Nothing left to do but fail out of the optimization */ 101250b47da0SAdam Denchfield ++cg->ls_fails; 1013c8bcdf1eSAdam Denchfield step = 0.0; 1014c8bcdf1eSAdam Denchfield tao->reason = TAO_DIVERGED_LS_FAILURE; 1015484c7b14SAdam Denchfield } else { 1016484c7b14SAdam Denchfield /* One of the fallbacks worked. Set them both back equal to false. */ 1017484c7b14SAdam Denchfield pcgd_fallback = PETSC_FALSE; 1018c8bcdf1eSAdam Denchfield } 1019c8bcdf1eSAdam Denchfield } 1020c8bcdf1eSAdam Denchfield } 1021c8bcdf1eSAdam Denchfield /* Convergence test for line search failure */ 1022c8bcdf1eSAdam Denchfield if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 1023c8bcdf1eSAdam Denchfield 1024c8bcdf1eSAdam Denchfield /* Standard convergence test */ 1025c8bcdf1eSAdam Denchfield ierr = VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W);CHKERRQ(ierr); 1026c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr); 1027*3c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(resnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 1028c8bcdf1eSAdam Denchfield ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 1029c8bcdf1eSAdam Denchfield ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr); 1030c8bcdf1eSAdam Denchfield ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 1031c8bcdf1eSAdam Denchfield if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 1032484c7b14SAdam Denchfield } 1033c8bcdf1eSAdam Denchfield /* Assert we have an updated step and we need at least one more iteration. */ 1034c8bcdf1eSAdam Denchfield /* Calculate the next direction */ 1035c8bcdf1eSAdam Denchfield /* Estimate the active set at the new solution */ 1036c8bcdf1eSAdam Denchfield ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr); 1037c8bcdf1eSAdam Denchfield /* Compute the projected gradient and its norm */ 1038c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 1039c8bcdf1eSAdam Denchfield ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr); 1040c8bcdf1eSAdam Denchfield ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 1041c8bcdf1eSAdam Denchfield gnorm2 = gnorm*gnorm; 1042c8bcdf1eSAdam Denchfield 1043484c7b14SAdam Denchfield /* Calculate some quantities used in the StepDirectionUpdate. */ 1044c8bcdf1eSAdam Denchfield ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm);CHKERRQ(ierr); 1045484c7b14SAdam Denchfield /* Update the step direction. */ 10468ca2df50S ierr = TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback);CHKERRQ(ierr); 1047484c7b14SAdam Denchfield ++tao->niter; 1048c8bcdf1eSAdam Denchfield ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1049c8bcdf1eSAdam Denchfield 1050c8bcdf1eSAdam Denchfield if (cg->cg_type != CG_GradientDescent) { 1051c8bcdf1eSAdam Denchfield /* Figure out which previously active variables became inactive this iteration */ 1052c8bcdf1eSAdam Denchfield ierr = ISDestroy(&cg->new_inactives);CHKERRQ(ierr); 1053c8bcdf1eSAdam Denchfield if (cg->inactive_idx && cg->inactive_old) { 1054c8bcdf1eSAdam Denchfield ierr = ISDifference(cg->inactive_idx, cg->inactive_old, &cg->new_inactives);CHKERRQ(ierr); 1055c8bcdf1eSAdam Denchfield } 1056c8bcdf1eSAdam Denchfield /* Selectively reset the CG step those freshly inactive variables */ 1057c8bcdf1eSAdam Denchfield if (cg->new_inactives) { 1058c8bcdf1eSAdam Denchfield ierr = VecGetSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr); 1059c8bcdf1eSAdam Denchfield ierr = VecGetSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr); 1060c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->inactive_grad, cg->inactive_step);CHKERRQ(ierr); 1061c8bcdf1eSAdam Denchfield ierr = VecScale(cg->inactive_step, -1.0);CHKERRQ(ierr); 1062c8bcdf1eSAdam Denchfield ierr = VecRestoreSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr); 1063c8bcdf1eSAdam Denchfield ierr = VecRestoreSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr); 1064c8bcdf1eSAdam Denchfield } 1065c8bcdf1eSAdam Denchfield /* Verify that this is a descent direction */ 1066c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 106750b47da0SAdam Denchfield ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm);CHKERRQ(ierr); 106850b47da0SAdam Denchfield if (PetscIsInfOrNanReal(gd) || (gd/(dnorm*dnorm) <= -1e10 || gd/(dnorm*dnorm) >= -1e-10)) { 1069c8bcdf1eSAdam Denchfield /* Not a descent direction, so we reset back to projected gradient descent */ 107050b47da0SAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 1071c8bcdf1eSAdam Denchfield ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1072c8bcdf1eSAdam Denchfield ++cg->descent_error; 1073c8bcdf1eSAdam Denchfield } else { 1074c8bcdf1eSAdam Denchfield } 1075c8bcdf1eSAdam Denchfield } 1076ac9112b8SAlp Dener PetscFunctionReturn(0); 1077ac9112b8SAlp Dener } 1078484c7b14SAdam Denchfield 1079484c7b14SAdam Denchfield PetscErrorCode TaoBNCGSetH0(Tao tao, Mat H0) 1080484c7b14SAdam Denchfield { 1081484c7b14SAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 1082484c7b14SAdam Denchfield PetscErrorCode ierr; 1083484c7b14SAdam Denchfield 1084484c7b14SAdam Denchfield PetscFunctionBegin; 1085484c7b14SAdam Denchfield ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 1086484c7b14SAdam Denchfield cg->pc = H0; 1087484c7b14SAdam Denchfield PetscFunctionReturn(0); 1088484c7b14SAdam Denchfield } 1089