1ac9112b8SAlp Dener #include <petsctaolinesearch.h> 2ac9112b8SAlp Dener #include <../src/tao/bound/impls/bncg/bncg.h> 3ac9112b8SAlp Dener 4*c8bcdf1eSAdam Denchfield #define CG_GradientDescent 0 5*c8bcdf1eSAdam Denchfield #define CG_HestenesStiefel 1 6*c8bcdf1eSAdam Denchfield #define CG_FletcherReeves 2 7*c8bcdf1eSAdam Denchfield #define CG_PolakRibiere 3 8*c8bcdf1eSAdam Denchfield #define CG_PolakRibierePlus 4 9*c8bcdf1eSAdam Denchfield #define CG_DaiYuan 5 10*c8bcdf1eSAdam Denchfield #define CG_HagerZhang 6 11*c8bcdf1eSAdam Denchfield #define CG_DaiKou 7 12*c8bcdf1eSAdam Denchfield #define CG_KouDai 8 13*c8bcdf1eSAdam Denchfield #define CG_SSML_BFGS 9 14*c8bcdf1eSAdam Denchfield #define CG_SSML_DFP 10 15*c8bcdf1eSAdam Denchfield #define CG_SSML_BROYDEN 11 16*c8bcdf1eSAdam Denchfield #define CG_BFGS_Mod 12 17*c8bcdf1eSAdam Denchfield #define CG_Types 13 18ac9112b8SAlp Dener 19*c8bcdf1eSAdam Denchfield static const char *CG_Table[64] = {"gd", "hs", "fr", "pr", "prp", "dy", "hz", "dk", "kd", "ssml_bfgs", "ssml_dfp", "ssml_brdn", "ssml_bfgs_mod"}; 20ac9112b8SAlp Dener 2161be54a6SAlp Dener #define CG_AS_NONE 0 2261be54a6SAlp Dener #define CG_AS_BERTSEKAS 1 2361be54a6SAlp Dener #define CG_AS_SIZE 2 24ac9112b8SAlp Dener 2561be54a6SAlp Dener static const char *CG_AS_TYPE[64] = {"none", "bertsekas"}; 26ac9112b8SAlp Dener 27c0f10754SAlp Dener PetscErrorCode TaoBNCGSetRecycleFlag(Tao tao, PetscBool recycle) 28c0f10754SAlp Dener { 29c0f10754SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 30c0f10754SAlp Dener 31c0f10754SAlp Dener PetscFunctionBegin; 32c0f10754SAlp Dener cg->recycle = recycle; 33c0f10754SAlp Dener PetscFunctionReturn(0); 34c0f10754SAlp Dener } 35c0f10754SAlp Dener 3661be54a6SAlp Dener PetscErrorCode TaoBNCGEstimateActiveSet(Tao tao, PetscInt asType) 3761be54a6SAlp Dener { 3861be54a6SAlp Dener PetscErrorCode ierr; 3961be54a6SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 4061be54a6SAlp Dener 4161be54a6SAlp Dener PetscFunctionBegin; 4261be54a6SAlp Dener ierr = ISDestroy(&cg->inactive_old);CHKERRQ(ierr); 4361be54a6SAlp Dener if (cg->inactive_idx) { 4461be54a6SAlp Dener ierr = ISDuplicate(cg->inactive_idx, &cg->inactive_old);CHKERRQ(ierr); 4561be54a6SAlp Dener ierr = ISCopy(cg->inactive_idx, cg->inactive_old);CHKERRQ(ierr); 4661be54a6SAlp Dener } 4761be54a6SAlp Dener switch (asType) { 4861be54a6SAlp Dener case CG_AS_NONE: 4961be54a6SAlp Dener ierr = ISDestroy(&cg->inactive_idx);CHKERRQ(ierr); 5061be54a6SAlp Dener ierr = VecWhichInactive(tao->XL, tao->solution, cg->unprojected_gradient, tao->XU, PETSC_TRUE, &cg->inactive_idx);CHKERRQ(ierr); 5161be54a6SAlp Dener ierr = ISDestroy(&cg->active_idx);CHKERRQ(ierr); 5261be54a6SAlp Dener ierr = ISComplementVec(cg->inactive_idx, tao->solution, &cg->active_idx);CHKERRQ(ierr); 5361be54a6SAlp Dener break; 5461be54a6SAlp Dener 5561be54a6SAlp Dener case CG_AS_BERTSEKAS: 5661be54a6SAlp Dener /* Use gradient descent to estimate the active set */ 5761be54a6SAlp Dener ierr = VecCopy(cg->unprojected_gradient, cg->W);CHKERRQ(ierr); 5861be54a6SAlp Dener ierr = VecScale(cg->W, -1.0);CHKERRQ(ierr); 5989da521bSAlp Dener ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, cg->unprojected_gradient, cg->W, cg->work, cg->as_step, &cg->as_tol, 6089da521bSAlp Dener &cg->active_lower, &cg->active_upper, &cg->active_fixed, &cg->active_idx, &cg->inactive_idx);CHKERRQ(ierr); 61c4b75bccSAlp Dener break; 6261be54a6SAlp Dener 6361be54a6SAlp Dener default: 6461be54a6SAlp Dener break; 6561be54a6SAlp Dener } 6661be54a6SAlp Dener PetscFunctionReturn(0); 6761be54a6SAlp Dener } 6861be54a6SAlp Dener 69a1318120SAlp Dener PetscErrorCode TaoBNCGBoundStep(Tao tao, PetscInt asType, Vec step) 7061be54a6SAlp Dener { 7161be54a6SAlp Dener PetscErrorCode ierr; 7261be54a6SAlp Dener TAO_BNCG *cg = (TAO_BNCG *)tao->data; 7361be54a6SAlp Dener 7461be54a6SAlp Dener PetscFunctionBegin; 75a1318120SAlp Dener switch (asType) { 7661be54a6SAlp Dener case CG_AS_NONE: 77c4b75bccSAlp Dener ierr = VecISSet(step, cg->active_idx, 0.0);CHKERRQ(ierr); 7861be54a6SAlp Dener break; 7961be54a6SAlp Dener 8061be54a6SAlp Dener case CG_AS_BERTSEKAS: 81c4b75bccSAlp Dener ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, cg->active_lower, cg->active_upper, cg->active_fixed, 1.0, step);CHKERRQ(ierr); 8261be54a6SAlp Dener break; 8361be54a6SAlp Dener 8461be54a6SAlp Dener default: 8561be54a6SAlp Dener break; 8661be54a6SAlp Dener } 8761be54a6SAlp Dener PetscFunctionReturn(0); 8861be54a6SAlp Dener } 8961be54a6SAlp Dener 90ac9112b8SAlp Dener static PetscErrorCode TaoSolve_BNCG(Tao tao) 91ac9112b8SAlp Dener { 92ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 93ac9112b8SAlp Dener PetscErrorCode ierr; 94*c8bcdf1eSAdam Denchfield PetscReal step=1.0,gnorm,gnorm2; 95c4b75bccSAlp Dener PetscInt nDiff; 96ac9112b8SAlp Dener 97ac9112b8SAlp Dener PetscFunctionBegin; 98ac9112b8SAlp Dener /* Project the current point onto the feasible set */ 99ac9112b8SAlp Dener ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 100cd929ea3SAlp Dener ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 101ac9112b8SAlp Dener 102*c8bcdf1eSAdam Denchfield /* Project the initial point onto the feasible region */ 103*c8bcdf1eSAdam Denchfield ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr); 104c0f10754SAlp Dener ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient);CHKERRQ(ierr); 105ac9112b8SAlp Dener ierr = VecNorm(cg->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr); 106c0f10754SAlp Dener if (PetscIsInfOrNanReal(cg->f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 107ac9112b8SAlp Dener 10861be54a6SAlp Dener /* Estimate the active set and compute the projected gradient */ 10961be54a6SAlp Dener ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr); 11061be54a6SAlp Dener 111ac9112b8SAlp Dener /* Project the gradient and calculate the norm */ 11261be54a6SAlp Dener ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 11361be54a6SAlp Dener ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr); 114ac9112b8SAlp Dener ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 115ac9112b8SAlp Dener gnorm2 = gnorm*gnorm; 116ac9112b8SAlp Dener 117*c8bcdf1eSAdam Denchfield /* Initialize counters */ 118e031d6f5SAlp Dener tao->niter = 0; 119*c8bcdf1eSAdam Denchfield cg->ls_fails = cg->broken_ortho = cg->descent_error = 0; 120*c8bcdf1eSAdam Denchfield cg->resets = -1; 121*c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 122*c8bcdf1eSAdam Denchfield 123*c8bcdf1eSAdam Denchfield /* Convergence test at the starting point. */ 124ac9112b8SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 125*c8bcdf1eSAdam Denchfield ierr = TaoLogConvergenceHistory(tao, cg->f, gnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 126*c8bcdf1eSAdam Denchfield ierr = TaoMonitor(tao, tao->niter, cg->f, gnorm, 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); 129ac9112b8SAlp Dener 130*c8bcdf1eSAdam Denchfield /* Assert that we have not converged. Calculate initial direction. */ 131*c8bcdf1eSAdam Denchfield /* This is where recycling goes. The outcome of this code must be */ 132*c8bcdf1eSAdam Denchfield /* the direction that we will use. */ 133*c8bcdf1eSAdam Denchfield if (cg->recycle) { 134ac9112b8SAlp Dener } 135*c8bcdf1eSAdam Denchfield /* Initial gradient descent step. Scaling by 1.0 also does a decent job for some problems. */ 136*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); 137ac9112b8SAlp Dener 138*c8bcdf1eSAdam Denchfield while(1) { 139*c8bcdf1eSAdam Denchfield ierr = TaoBNCGConductIteration(tao, gnorm); CHKERRQ(ierr); 140*c8bcdf1eSAdam Denchfield if (cg->use_steffenson) ierr = TaoBNCGSteffensonAcceleration(tao); CHKERRQ(ierr); 141*c8bcdf1eSAdam Denchfield if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 142ac9112b8SAlp Dener } 143ac9112b8SAlp Dener PetscFunctionReturn(0); 144ac9112b8SAlp Dener } 145ac9112b8SAlp Dener 146ac9112b8SAlp Dener static PetscErrorCode TaoSetUp_BNCG(Tao tao) 147ac9112b8SAlp Dener { 148ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 149ac9112b8SAlp Dener PetscErrorCode ierr; 150ac9112b8SAlp Dener 151ac9112b8SAlp Dener PetscFunctionBegin; 152c4b75bccSAlp Dener if (!tao->gradient) { 153c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 154c4b75bccSAlp Dener } 155c4b75bccSAlp Dener if (!tao->stepdirection) { 156c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); 157c4b75bccSAlp Dener } 158c4b75bccSAlp Dener if (!cg->W) { 159c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&cg->W);CHKERRQ(ierr); 160c4b75bccSAlp Dener } 161c4b75bccSAlp Dener if (!cg->work) { 162c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&cg->work);CHKERRQ(ierr); 163c4b75bccSAlp Dener } 164*c8bcdf1eSAdam Denchfield if (!cg->sk) { 165*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->sk);CHKERRQ(ierr); 166*c8bcdf1eSAdam Denchfield } 167*c8bcdf1eSAdam Denchfield if (!cg->yk) { 168*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->gradient,&cg->yk);CHKERRQ(ierr); 169*c8bcdf1eSAdam Denchfield } 170c4b75bccSAlp Dener if (!cg->X_old) { 171c4b75bccSAlp Dener ierr = VecDuplicate(tao->solution,&cg->X_old);CHKERRQ(ierr); 172c4b75bccSAlp Dener } 173c4b75bccSAlp Dener if (!cg->G_old) { 174c4b75bccSAlp Dener ierr = VecDuplicate(tao->gradient,&cg->G_old);CHKERRQ(ierr); 175*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->g_work);CHKERRQ(ierr); 176*c8bcdf1eSAdam Denchfield } 177*c8bcdf1eSAdam Denchfield if (cg->use_steffenson){ 178*c8bcdf1eSAdam Denchfield if (!cg->steffnm1) ierr = VecDuplicate(tao->solution, &cg->steffnm1); CHKERRQ(ierr); 179*c8bcdf1eSAdam Denchfield if (!cg->steffn) ierr = VecDuplicate(tao->solution, &cg->steffn); CHKERRQ(ierr); 180*c8bcdf1eSAdam Denchfield if (!cg->steffnp1) ierr = VecDuplicate(tao->solution, &cg->steffnp1); CHKERRQ(ierr); 181*c8bcdf1eSAdam Denchfield if (!cg->steffva) ierr = VecDuplicate(tao->solution, &cg->steffva); CHKERRQ(ierr); 182*c8bcdf1eSAdam Denchfield if (!cg->steffvatmp) ierr = VecDuplicate(tao->solution, &cg->steffvatmp); CHKERRQ(ierr); 183*c8bcdf1eSAdam Denchfield } 184*c8bcdf1eSAdam Denchfield if (cg->diag_scaling){ 185*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->gradient,&cg->invD);CHKERRQ(ierr); 186*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->gradient,&cg->invDnew);CHKERRQ(ierr); 187*c8bcdf1eSAdam Denchfield ierr = VecSet(cg->invDnew, 1.0);CHKERRQ(ierr); 188*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->gradient,&cg->bfgs_work);CHKERRQ(ierr); 189*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->gradient,&cg->dfp_work);CHKERRQ(ierr); 190*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->gradient,&cg->U);CHKERRQ(ierr); 191*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->gradient,&cg->V);CHKERRQ(ierr); 192*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->gradient,&cg->W);CHKERRQ(ierr); 193*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->d_work);CHKERRQ(ierr); 194*c8bcdf1eSAdam Denchfield ierr = VecDuplicate(tao->solution,&cg->y_work);CHKERRQ(ierr); 195*c8bcdf1eSAdam Denchfield 196c4b75bccSAlp Dener } 197c4b75bccSAlp Dener if (!cg->unprojected_gradient) { 198c4b75bccSAlp Dener ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient);CHKERRQ(ierr); 199c4b75bccSAlp Dener } 200c4b75bccSAlp Dener if (!cg->unprojected_gradient_old) { 201c4b75bccSAlp Dener ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient_old);CHKERRQ(ierr); 202c4b75bccSAlp Dener } 203ac9112b8SAlp Dener PetscFunctionReturn(0); 204ac9112b8SAlp Dener } 205ac9112b8SAlp Dener 206ac9112b8SAlp Dener static PetscErrorCode TaoDestroy_BNCG(Tao tao) 207ac9112b8SAlp Dener { 208ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*) tao->data; 209ac9112b8SAlp Dener PetscErrorCode ierr; 210ac9112b8SAlp Dener 211ac9112b8SAlp Dener PetscFunctionBegin; 212ac9112b8SAlp Dener if (tao->setupcalled) { 21361be54a6SAlp Dener ierr = VecDestroy(&cg->W);CHKERRQ(ierr); 214c4b75bccSAlp Dener ierr = VecDestroy(&cg->work);CHKERRQ(ierr); 215ac9112b8SAlp Dener ierr = VecDestroy(&cg->X_old);CHKERRQ(ierr); 216ac9112b8SAlp Dener ierr = VecDestroy(&cg->G_old);CHKERRQ(ierr); 217ac9112b8SAlp Dener ierr = VecDestroy(&cg->unprojected_gradient);CHKERRQ(ierr); 218ac9112b8SAlp Dener ierr = VecDestroy(&cg->unprojected_gradient_old);CHKERRQ(ierr); 219ac9112b8SAlp Dener } 220ca964c49SAlp Dener ierr = ISDestroy(&cg->active_lower);CHKERRQ(ierr); 221ca964c49SAlp Dener ierr = ISDestroy(&cg->active_upper);CHKERRQ(ierr); 222ca964c49SAlp Dener ierr = ISDestroy(&cg->active_fixed);CHKERRQ(ierr); 223ca964c49SAlp Dener ierr = ISDestroy(&cg->active_idx);CHKERRQ(ierr); 224ca964c49SAlp Dener ierr = ISDestroy(&cg->inactive_idx);CHKERRQ(ierr); 225ca964c49SAlp Dener ierr = ISDestroy(&cg->inactive_old);CHKERRQ(ierr); 226ca964c49SAlp Dener ierr = ISDestroy(&cg->new_inactives);CHKERRQ(ierr); 227ac9112b8SAlp Dener ierr = PetscFree(tao->data);CHKERRQ(ierr); 228ac9112b8SAlp Dener PetscFunctionReturn(0); 229ac9112b8SAlp Dener } 230ac9112b8SAlp Dener 231ac9112b8SAlp Dener static PetscErrorCode TaoSetFromOptions_BNCG(PetscOptionItems *PetscOptionsObject,Tao tao) 232ac9112b8SAlp Dener { 233ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 234ac9112b8SAlp Dener PetscErrorCode ierr; 235ac9112b8SAlp Dener 236ac9112b8SAlp Dener PetscFunctionBegin; 237ac9112b8SAlp Dener ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 238ac9112b8SAlp Dener ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Conjugate Gradient method for unconstrained optimization");CHKERRQ(ierr); 239*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_eta","cutoff tolerance for HZ", "", cg->eta,&cg->eta,NULL);CHKERRQ(ierr); 240*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_xi","Parameter in KD, HZ, and DK methods", "", cg->xi,&cg->xi,NULL);CHKERRQ(ierr); 241*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_gamma", "stabilization term for multiple CG methods", "", cg->gamma, &cg->gamma, NULL); CHKERRQ(ierr); 242*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_theta", "update parameter for some CG methods", "", cg->theta, &cg->theta, NULL); CHKERRQ(ierr); 243*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_hz_theta", "parameter for the HZ (2006) method", "", cg->hz_theta, &cg->hz_theta, NULL); CHKERRQ(ierr); 244*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_rho","(developer) update limiter in the J0 scaling","",cg->rho,&cg->rho,NULL);CHKERRQ(ierr); 245*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_alpha","(developer) convex ratio in the J0 scaling","",cg->alpha,&cg->alpha,NULL);CHKERRQ(ierr); 246*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_beta","(developer) exponential factor in the diagonal J0 scaling","",cg->beta,&cg->beta,NULL);CHKERRQ(ierr); 247*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_bfgs_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->bfgs_scale, &cg->bfgs_scale, NULL); CHKERRQ(ierr); 248*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_dfp_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->dfp_scale, &cg->dfp_scale, NULL); CHKERRQ(ierr); 249*c8bcdf1eSAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_diag_scaling","Enable diagonal Broyden-like preconditioning","",cg->diag_scaling,&cg->diag_scaling,NULL);CHKERRQ(ierr); 250*c8bcdf1eSAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_inv_sig","(developer) test parameter to invert the sigma scaling","",cg->inv_sig,&cg->inv_sig,NULL);CHKERRQ(ierr); 251*c8bcdf1eSAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_dynamic_restart","use dynamic restarts as in HZ, DK, KD","",cg->use_dynamic_restart,&cg->use_dynamic_restart,NULL);CHKERRQ(ierr); 252*c8bcdf1eSAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_use_steffenson","(incomplete) use vector-Steffenson acceleration","",cg->use_steffenson,&cg->use_steffenson,NULL);CHKERRQ(ierr); 253*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_zeta", "Free parameter for the Kou-Dai method", "", cg->zeta, &cg->zeta, NULL); CHKERRQ(ierr); 254*c8bcdf1eSAdam 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); 255*c8bcdf1eSAdam Denchfield ierr = PetscOptionsInt("-tao_bncg_min_restart_num", "Number of iterations between restarts (times dimension)", "", cg->min_restart_num, &cg->min_restart_num, NULL); CHKERRQ(ierr); 256*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_rho","(developer) descent direction tolerance", "", cg->rho,&cg->rho,NULL);CHKERRQ(ierr); 257*c8bcdf1eSAdam Denchfield ierr = PetscOptionsReal("-tao_bncg_pow","(developer) descent direction exponent", "", cg->pow,&cg->pow,NULL);CHKERRQ(ierr); 25861be54a6SAlp Dener ierr = PetscOptionsEList("-tao_bncg_type","cg formula", "", CG_Table, CG_Types, CG_Table[cg->cg_type], &cg->cg_type,NULL);CHKERRQ(ierr); 259*c8bcdf1eSAdam 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); 26061be54a6SAlp Dener ierr = PetscOptionsBool("-tao_bncg_recycle","enable recycling the existing solution and gradient at the start of a new solve","",cg->recycle,&cg->recycle,NULL);CHKERRQ(ierr); 261*c8bcdf1eSAdam Denchfield ierr = PetscOptionsBool("-tao_bncg_spaced_restart","Enable regular steepest descent restarting every ixed number of iterations","",cg->spaced_restart,&cg->spaced_restart,NULL);CHKERRQ(ierr); 262*c8bcdf1eSAdam 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); 26361be54a6SAlp Dener ierr = PetscOptionsReal("-tao_bncg_as_tol", "initial tolerance used when estimating actively bounded variables","",cg->as_tol,&cg->as_tol,NULL);CHKERRQ(ierr); 26461be54a6SAlp Dener ierr = PetscOptionsReal("-tao_bncg_as_step", "step length used when estimating actively bounded variables","",cg->as_step,&cg->as_step,NULL);CHKERRQ(ierr); 265ac9112b8SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 266ac9112b8SAlp Dener PetscFunctionReturn(0); 267ac9112b8SAlp Dener } 268ac9112b8SAlp Dener 269ac9112b8SAlp Dener static PetscErrorCode TaoView_BNCG(Tao tao, PetscViewer viewer) 270ac9112b8SAlp Dener { 271ac9112b8SAlp Dener PetscBool isascii; 272ac9112b8SAlp Dener TAO_BNCG *cg = (TAO_BNCG*)tao->data; 273ac9112b8SAlp Dener PetscErrorCode ierr; 274ac9112b8SAlp Dener 275ac9112b8SAlp Dener PetscFunctionBegin; 276ac9112b8SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 277ac9112b8SAlp Dener if (isascii) { 278ac9112b8SAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 279ac9112b8SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cg->cg_type]);CHKERRQ(ierr); 280ac9112b8SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Resets: %i\n", cg->resets);CHKERRQ(ierr); 281ac9112b8SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " Broken ortho: %i\n", cg->broken_ortho);CHKERRQ(ierr); 282ac9112b8SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " Not a descent dir.: %i\n", cg->descent_error);CHKERRQ(ierr); 283ac9112b8SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "Line search fails: %i\n", cg->ls_fails);CHKERRQ(ierr); 284ac9112b8SAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 285ac9112b8SAlp Dener } 286ac9112b8SAlp Dener PetscFunctionReturn(0); 287ac9112b8SAlp Dener } 288ac9112b8SAlp Dener 289*c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal yty, PetscReal yts, PetscReal sts, PetscReal *scale, PetscReal alpha) 290*c8bcdf1eSAdam Denchfield { 291*c8bcdf1eSAdam Denchfield PetscReal a, b, c, sig1, sig2; 292*c8bcdf1eSAdam Denchfield 293*c8bcdf1eSAdam Denchfield PetscFunctionBegin; 294*c8bcdf1eSAdam Denchfield *scale = 0.0; 295*c8bcdf1eSAdam Denchfield 296*c8bcdf1eSAdam Denchfield if (alpha == 1.0){ 297*c8bcdf1eSAdam Denchfield *scale = yts/yty; 298*c8bcdf1eSAdam Denchfield } else if (alpha == 0.5) { 299*c8bcdf1eSAdam Denchfield *scale = sts/yty; 300*c8bcdf1eSAdam Denchfield } 301*c8bcdf1eSAdam Denchfield else if (alpha == 0.0){ 302*c8bcdf1eSAdam Denchfield *scale = sts/yts; 303*c8bcdf1eSAdam Denchfield } 304*c8bcdf1eSAdam Denchfield else if (alpha == -1.0) *scale = 1.0; 305*c8bcdf1eSAdam Denchfield else { 306*c8bcdf1eSAdam Denchfield a = yty; 307*c8bcdf1eSAdam Denchfield b = yts; 308*c8bcdf1eSAdam Denchfield c = sts; 309*c8bcdf1eSAdam Denchfield a *= alpha; 310*c8bcdf1eSAdam Denchfield b *= -(2.0*alpha - 1.0); 311*c8bcdf1eSAdam Denchfield c *= alpha - 1.0; 312*c8bcdf1eSAdam Denchfield sig1 = (-b + PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a); 313*c8bcdf1eSAdam Denchfield sig2 = (-b - PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a); 314*c8bcdf1eSAdam Denchfield /* accept the positive root as the scalar */ 315*c8bcdf1eSAdam Denchfield if (sig1 > 0.0) { 316*c8bcdf1eSAdam Denchfield *scale = sig1; 317*c8bcdf1eSAdam Denchfield } else if (sig2 > 0.0) { 318*c8bcdf1eSAdam Denchfield *scale = sig2; 319*c8bcdf1eSAdam Denchfield } else { 320*c8bcdf1eSAdam Denchfield SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Cannot find positive scalar"); 321*c8bcdf1eSAdam Denchfield } 322*c8bcdf1eSAdam Denchfield } 323*c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 324*c8bcdf1eSAdam Denchfield } 325*c8bcdf1eSAdam Denchfield 326ac9112b8SAlp Dener /*MC 327ac9112b8SAlp Dener TAOBNCG - Bound-constrained Nonlinear Conjugate Gradient method. 328ac9112b8SAlp Dener 329ac9112b8SAlp Dener Options Database Keys: 330c4b75bccSAlp Dener + -tao_bncg_recycle - enable recycling the latest calculated gradient vector in subsequent TaoSolve() calls 331c4b75bccSAlp Dener . -tao_bncg_eta <r> - restart tolerance 33261be54a6SAlp Dener . -tao_bncg_type <taocg_type> - cg formula 333c4b75bccSAlp Dener . -tao_bncg_as_type <none,bertsekas> - active set estimation method 334c4b75bccSAlp Dener . -tao_bncg_as_tol <r> - tolerance used in Bertsekas active-set estimation 335c4b75bccSAlp Dener . -tao_bncg_as_step <r> - trial step length used in Bertsekas active-set estimation 336ac9112b8SAlp Dener 337ac9112b8SAlp Dener Notes: 338ac9112b8SAlp Dener CG formulas are: 339ac9112b8SAlp Dener "fr" - Fletcher-Reeves 340ac9112b8SAlp Dener "pr" - Polak-Ribiere 341ac9112b8SAlp Dener "prp" - Polak-Ribiere-Plus 342ac9112b8SAlp Dener "hs" - Hestenes-Steifel 343ac9112b8SAlp Dener "dy" - Dai-Yuan 344560169d0SAlp Dener "gd" - Gradient Descent 345ac9112b8SAlp Dener Level: beginner 346ac9112b8SAlp Dener M*/ 347ac9112b8SAlp Dener 348ac9112b8SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNCG(Tao tao) 349ac9112b8SAlp Dener { 350ac9112b8SAlp Dener TAO_BNCG *cg; 351ac9112b8SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 352ac9112b8SAlp Dener PetscErrorCode ierr; 353ac9112b8SAlp Dener 354ac9112b8SAlp Dener PetscFunctionBegin; 355ac9112b8SAlp Dener tao->ops->setup = TaoSetUp_BNCG; 356ac9112b8SAlp Dener tao->ops->solve = TaoSolve_BNCG; 357ac9112b8SAlp Dener tao->ops->view = TaoView_BNCG; 358ac9112b8SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNCG; 359ac9112b8SAlp Dener tao->ops->destroy = TaoDestroy_BNCG; 360ac9112b8SAlp Dener 361ac9112b8SAlp Dener /* Override default settings (unless already changed) */ 362ac9112b8SAlp Dener if (!tao->max_it_changed) tao->max_it = 2000; 363ac9112b8SAlp Dener if (!tao->max_funcs_changed) tao->max_funcs = 4000; 364ac9112b8SAlp Dener 365ac9112b8SAlp Dener /* Note: nondefault values should be used for nonlinear conjugate gradient */ 366ac9112b8SAlp Dener /* method. In particular, gtol should be less that 0.5; the value used in */ 367ac9112b8SAlp Dener /* Nocedal and Wright is 0.10. We use the default values for the */ 368ac9112b8SAlp Dener /* linesearch because it seems to work better. */ 369ac9112b8SAlp Dener ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 370ac9112b8SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 371ac9112b8SAlp Dener ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 372ac9112b8SAlp Dener ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr); 373ac9112b8SAlp Dener ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 374ac9112b8SAlp Dener 375ac9112b8SAlp Dener ierr = PetscNewLog(tao,&cg);CHKERRQ(ierr); 376ac9112b8SAlp Dener tao->data = (void*)cg; 377*c8bcdf1eSAdam Denchfield 378ac9112b8SAlp Dener cg->pow = 2.1; 379ac9112b8SAlp Dener cg->eta = 0.5; 380*c8bcdf1eSAdam Denchfield cg->dynamic_restart = PETSC_FALSE; 381*c8bcdf1eSAdam Denchfield cg->unscaled_restart = PETSC_FALSE; 382*c8bcdf1eSAdam Denchfield cg->theta = 1.0; 383*c8bcdf1eSAdam Denchfield cg->hz_theta = 1.0; 384*c8bcdf1eSAdam Denchfield cg->dfp_scale = 1.0; 385*c8bcdf1eSAdam Denchfield cg->bfgs_scale = 1.0; 386*c8bcdf1eSAdam Denchfield cg->gamma = 0.4; 387*c8bcdf1eSAdam Denchfield cg->zeta = 0.5; 388*c8bcdf1eSAdam Denchfield cg->min_quad = 3; 389*c8bcdf1eSAdam Denchfield cg->min_restart_num = 6; /* As in CG_DESCENT and KD2015*/ 390*c8bcdf1eSAdam Denchfield cg->xi = 1.0; 391*c8bcdf1eSAdam Denchfield cg->neg_xi = PETSC_FALSE; 392*c8bcdf1eSAdam Denchfield cg->spaced_restart = PETSC_FALSE; 393*c8bcdf1eSAdam Denchfield cg->tol_quad = 1e-8; 39461be54a6SAlp Dener cg->as_step = 0.001; 39561be54a6SAlp Dener cg->as_tol = 0.001; 396*c8bcdf1eSAdam Denchfield cg->epsilon = PETSC_MACHINE_EPSILON; 397*c8bcdf1eSAdam Denchfield cg->eps_23 = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); 39861be54a6SAlp Dener cg->as_type = CG_AS_BERTSEKAS; 399*c8bcdf1eSAdam Denchfield cg->cg_type = CG_SSML_BFGS; 400c0f10754SAlp Dener cg->recycle = PETSC_FALSE; 401*c8bcdf1eSAdam Denchfield cg->alpha = 1.0; 402*c8bcdf1eSAdam Denchfield cg->rho = 1.0; 403*c8bcdf1eSAdam Denchfield cg->beta = 0.5; /* Default a la Alp */ 404*c8bcdf1eSAdam Denchfield cg->diag_scaling = PETSC_TRUE; 405*c8bcdf1eSAdam Denchfield cg->inv_sig = PETSC_FALSE; 406*c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 407*c8bcdf1eSAdam Denchfield } 408*c8bcdf1eSAdam Denchfield 409*c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGComputeDiagScaling(Tao tao, PetscReal yts, PetscReal yty){ 410*c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 411*c8bcdf1eSAdam Denchfield PetscErrorCode ierr; 412*c8bcdf1eSAdam Denchfield PetscScalar ytDy, ytDs, stDs; 413*c8bcdf1eSAdam Denchfield PetscReal sigma; 414*c8bcdf1eSAdam Denchfield 415*c8bcdf1eSAdam Denchfield PetscFunctionBegin; 416*c8bcdf1eSAdam Denchfield /* W = Hy */ 417*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->W, cg->invD, cg->yk);CHKERRQ(ierr); 418*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->W, cg->yk, &ytDy);CHKERRQ(ierr); 419*c8bcdf1eSAdam Denchfield /* Compute Hadamard product V = s*s^T */ 420*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->V, cg->sk, cg->sk);CHKERRQ(ierr); 421*c8bcdf1eSAdam Denchfield /* Compute the BFGS contribution bfgs_work */ 422*c8bcdf1eSAdam Denchfield /* first construct sDy - Hadamard product */ 423*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->bfgs_work, cg->sk, cg->W); CHKERRQ(ierr); 424*c8bcdf1eSAdam Denchfield /* Now assemble the BFGS component in there - denom of yts added later */ 425*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(cg->bfgs_work, ytDy/(yts), -2.0, cg->V); CHKERRQ(ierr); 426*c8bcdf1eSAdam Denchfield /* Start assembling the new inverse diagonal, the pure DFP component - denom of ytDy added later */ 427*c8bcdf1eSAdam Denchfield /* Compute Hadamard product U = (Hy)*(Hy)^T */ 428*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->U, cg->W, cg->W); CHKERRQ(ierr); 429*c8bcdf1eSAdam Denchfield /* D_{k+1} = D_k + V/yts + (1-theta)*BFGS + theta*DFP */ 430*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->invD, cg->invDnew); CHKERRQ(ierr); 431*c8bcdf1eSAdam Denchfield /* The factors I left out in BFGS and DFP */ 432*c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(cg->invDnew, (1-cg->theta)/yts, -cg->theta/(ytDy), 1.0, cg->bfgs_work, cg->U); CHKERRQ(ierr); 433*c8bcdf1eSAdam Denchfield ierr = VecAXPY(cg->invDnew, 1/cg->yts, cg->V); 434*c8bcdf1eSAdam Denchfield /* Ensure positive definite */ 435*c8bcdf1eSAdam Denchfield ierr = VecAbs(cg->invDnew); CHKERRQ(ierr); 436*c8bcdf1eSAdam Denchfield /* Start with re-scaling on the newly computed diagonal */ 437*c8bcdf1eSAdam Denchfield /* Compute y^T H^{2*beta} y */ 438*c8bcdf1eSAdam Denchfield /* Note that VecPow has special cases tabulated for +-1.0, +-0.5, 0.0, and 2.0 */ 439*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->invDnew, cg->work);CHKERRQ(ierr); 440*c8bcdf1eSAdam Denchfield ierr = VecPow(cg->work, 2.0*cg->beta);CHKERRQ(ierr); 441*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->work, cg->work, cg->yk);CHKERRQ(ierr); 442*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->work, &ytDy);CHKERRQ(ierr); 443*c8bcdf1eSAdam Denchfield /* Compute y^T H^{2*beta - 1} s */ 444*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->invDnew, cg->work);CHKERRQ(ierr); 445*c8bcdf1eSAdam Denchfield ierr = VecPow(cg->work, 2.0*cg->beta - 1.0);CHKERRQ(ierr); 446*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->work, cg->work, cg->sk);CHKERRQ(ierr); 447*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->work, &ytDs);CHKERRQ(ierr); 448*c8bcdf1eSAdam Denchfield /* Compute s^T H^{2*beta - 2} s */ 449*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->invDnew, cg->work);CHKERRQ(ierr); 450*c8bcdf1eSAdam Denchfield ierr = VecPow(cg->work, 2.0*cg->beta - 2.0);CHKERRQ(ierr); 451*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->work, cg->work, cg->sk);CHKERRQ(ierr); 452*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->sk, cg->work, &stDs);CHKERRQ(ierr); 453*c8bcdf1eSAdam Denchfield /* Compute the diagonal scaling */ 454*c8bcdf1eSAdam Denchfield sigma = 0.0; 455*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ytDy, ytDs, stDs, &sigma, cg->alpha); 456*c8bcdf1eSAdam Denchfield 457*c8bcdf1eSAdam Denchfield /* If Q has small values, then Q^(r_beta - 1) */ 458*c8bcdf1eSAdam Denchfield /* can have very large values. Hence, ys_sum */ 459*c8bcdf1eSAdam Denchfield /* and ss_sum can be infinity. In this case, */ 460*c8bcdf1eSAdam Denchfield /* sigma can either be not-a-number or infinity. */ 461*c8bcdf1eSAdam Denchfield 462*c8bcdf1eSAdam Denchfield if (PetscIsInfOrNanReal(sigma)) { 463*c8bcdf1eSAdam Denchfield /* sigma is not-a-number; skip rescaling */ 464*c8bcdf1eSAdam Denchfield } else if (!sigma) { 465*c8bcdf1eSAdam Denchfield /* sigma is zero; this is a bad case; skip rescaling */ 466*c8bcdf1eSAdam Denchfield } else { 467*c8bcdf1eSAdam Denchfield /* sigma is positive */ 468*c8bcdf1eSAdam Denchfield ierr = VecScale(cg->invDnew, sigma);CHKERRQ(ierr); 469*c8bcdf1eSAdam Denchfield } 470*c8bcdf1eSAdam Denchfield 471*c8bcdf1eSAdam Denchfield /* Combine the old diagonal and the new diagonal using a convex limiter */ 472*c8bcdf1eSAdam Denchfield if (cg->rho == 1.0) { 473*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->invDnew, cg->invD);CHKERRQ(ierr); 474*c8bcdf1eSAdam Denchfield } else if (cg->rho) { 475*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(cg->invD, 1.0-cg->rho, cg->rho, cg->invDnew);CHKERRQ(ierr); 476*c8bcdf1eSAdam Denchfield } 477*c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 478*c8bcdf1eSAdam Denchfield } 479*c8bcdf1eSAdam Denchfield 480*c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGResetUpdate(Tao tao, PetscReal gnormsq) 481*c8bcdf1eSAdam Denchfield { 482*c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 483*c8bcdf1eSAdam Denchfield PetscErrorCode ierr; 484*c8bcdf1eSAdam Denchfield PetscReal scaling; 485*c8bcdf1eSAdam Denchfield 486*c8bcdf1eSAdam Denchfield PetscFunctionBegin; 487*c8bcdf1eSAdam Denchfield ++cg->resets; 488*c8bcdf1eSAdam Denchfield scaling = 2.0 * PetscMax(1.0, PetscAbsScalar(cg->f)) / PetscMax(gnormsq, cg->eps_23); 489*c8bcdf1eSAdam Denchfield scaling = PetscMin(100.0, PetscMax(1e-7, scaling)); 490*c8bcdf1eSAdam Denchfield if (cg->unscaled_restart) scaling = 1.0; 491*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -scaling, 0.0, tao->gradient);CHKERRQ(ierr); 492*c8bcdf1eSAdam Denchfield /* Also want to reset our diagonal scaling with each restart */ 493*c8bcdf1eSAdam Denchfield if (cg->diag_scaling) { 494*c8bcdf1eSAdam Denchfield ierr = VecSet(cg->invD, 1.0);CHKERRQ(ierr); 495*c8bcdf1eSAdam Denchfield } 496*c8bcdf1eSAdam Denchfield 497*c8bcdf1eSAdam Denchfield 498*c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 499*c8bcdf1eSAdam Denchfield } 500*c8bcdf1eSAdam Denchfield 501*c8bcdf1eSAdam Denchfield PetscErrorCode TaoBNCGCheckDynamicRestart(Tao tao, PetscReal stepsize, PetscReal gd, PetscReal gd_old, PetscBool *dynrestart, PetscReal fold) 502*c8bcdf1eSAdam Denchfield { 503*c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 504*c8bcdf1eSAdam Denchfield PetscReal quadinterp; 505*c8bcdf1eSAdam Denchfield 506*c8bcdf1eSAdam Denchfield PetscFunctionBegin; 507*c8bcdf1eSAdam Denchfield if (cg->f < cg->min_quad/10) {*dynrestart = PETSC_FALSE; PetscFunctionReturn(0);} /* just skip this since this strategy doesn't work well for functions near zero */ 508*c8bcdf1eSAdam Denchfield quadinterp = 2*(cg->f - fold)/(stepsize*(gd + gd_old)); 509*c8bcdf1eSAdam Denchfield if (PetscAbs(quadinterp - 1.0) < cg->tol_quad) cg->iter_quad++; 510*c8bcdf1eSAdam Denchfield else { 511*c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 512*c8bcdf1eSAdam Denchfield *dynrestart = PETSC_FALSE; 513*c8bcdf1eSAdam Denchfield } 514*c8bcdf1eSAdam Denchfield 515*c8bcdf1eSAdam Denchfield if (cg->iter_quad >= cg->min_quad) { 516*c8bcdf1eSAdam Denchfield cg->iter_quad = 0; 517*c8bcdf1eSAdam Denchfield *dynrestart = PETSC_TRUE; 518*c8bcdf1eSAdam Denchfield } 519*c8bcdf1eSAdam Denchfield 520*c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 521*c8bcdf1eSAdam Denchfield } 522*c8bcdf1eSAdam Denchfield 523*c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao tao, PetscReal gnorm2, PetscReal step, PetscReal fold, PetscReal gnorm2_old, PetscBool cg_restart, PetscReal dnorm, PetscReal ginner){ 524*c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 525*c8bcdf1eSAdam Denchfield PetscErrorCode ierr; 526*c8bcdf1eSAdam Denchfield PetscReal gamma, tau_k, beta; 527*c8bcdf1eSAdam Denchfield PetscReal tmp, ynorm, ynorm2, snorm, dk_yk, gd; 528*c8bcdf1eSAdam Denchfield PetscReal gkp1_yk, gd_old, tau_bfgs, tau_dfp, gkp1D_yk; 529*c8bcdf1eSAdam Denchfield PetscInt dim; 530*c8bcdf1eSAdam Denchfield 531*c8bcdf1eSAdam Denchfield PetscFunctionBegin; 532*c8bcdf1eSAdam Denchfield 533*c8bcdf1eSAdam Denchfield /* Want to implement P.C. versions eventually */ 534*c8bcdf1eSAdam Denchfield /* Compute CG step */ 535*c8bcdf1eSAdam Denchfield if (cg_restart) { 536*c8bcdf1eSAdam Denchfield beta = 0.0; 537*c8bcdf1eSAdam Denchfield ++cg->resets; 538*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient);CHKERRQ(ierr); 539*c8bcdf1eSAdam Denchfield } else { 540*c8bcdf1eSAdam Denchfield switch (cg->cg_type) { 541*c8bcdf1eSAdam Denchfield case CG_GradientDescent: 542*c8bcdf1eSAdam Denchfield beta = 0.0; 543*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, 0.0, tao->gradient);CHKERRQ(ierr); 544*c8bcdf1eSAdam Denchfield break; 545*c8bcdf1eSAdam Denchfield 546*c8bcdf1eSAdam Denchfield case CG_HestenesStiefel: 547*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 548*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 549*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient); CHKERRQ(ierr); 550*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm); CHKERRQ(ierr); 551*c8bcdf1eSAdam Denchfield ynorm2 = ynorm*ynorm; 552*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk); CHKERRQ(ierr); 553*c8bcdf1eSAdam Denchfield //tau_k = step*dk_yk/ynorm2; 554*c8bcdf1eSAdam Denchfield //beta = tau_k*(gnorm2 - ginner) / (gd - gd_old); 555*c8bcdf1eSAdam Denchfield beta = (gnorm2 - ginner) / (gd - gd_old); 556*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient);CHKERRQ(ierr); 557*c8bcdf1eSAdam Denchfield break; 558*c8bcdf1eSAdam Denchfield case CG_FletcherReeves: 559*c8bcdf1eSAdam Denchfield beta = gnorm2 / gnorm2_old; 560*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient);CHKERRQ(ierr); 561*c8bcdf1eSAdam Denchfield break; 562*c8bcdf1eSAdam Denchfield case CG_PolakRibiere: 563*c8bcdf1eSAdam Denchfield beta = (gnorm2 - ginner) / gnorm2_old; 564*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient);CHKERRQ(ierr); 565*c8bcdf1eSAdam Denchfield break; 566*c8bcdf1eSAdam Denchfield case CG_PolakRibierePlus: 567*c8bcdf1eSAdam Denchfield beta = PetscMax((gnorm2-ginner)/gnorm2_old, 0.0); 568*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient);CHKERRQ(ierr); 569*c8bcdf1eSAdam Denchfield break; 570*c8bcdf1eSAdam Denchfield case CG_DaiYuan: 571*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 572*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 573*c8bcdf1eSAdam Denchfield beta = gnorm2 / (gd - gd_old); 574*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient);CHKERRQ(ierr); 575*c8bcdf1eSAdam Denchfield break; 576*c8bcdf1eSAdam Denchfield case CG_SSML_BFGS: 577*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 578*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient); CHKERRQ(ierr); 579*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution); CHKERRQ(ierr); 580*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm); CHKERRQ(ierr); 581*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->sk, NORM_2, &snorm); CHKERRQ(ierr); 582*c8bcdf1eSAdam Denchfield ynorm2 = ynorm*ynorm; 583*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk); CHKERRQ(ierr); 584*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->sk, &cg->yts); CHKERRQ(ierr); 585*c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 586*c8bcdf1eSAdam Denchfield cg->sts = snorm*snorm; 587*c8bcdf1eSAdam Denchfield /* TODO: Should only need one of dk_yk and yts */ 588*c8bcdf1eSAdam Denchfield if (ynorm2 < cg->epsilon){ 589*c8bcdf1eSAdam Denchfield /* The gradient hasn't changed, so we should stay in the same direction as before. Don't update it or anything.*/ 590*c8bcdf1eSAdam Denchfield if (snorm < cg->eps_23){ 591*c8bcdf1eSAdam Denchfield /* We're making no progress. Scaled gradient descent step.*/ 592*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 593*c8bcdf1eSAdam Denchfield } 594*c8bcdf1eSAdam Denchfield } else { 595*c8bcdf1eSAdam Denchfield if (!cg->diag_scaling){ 596*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk); CHKERRQ(ierr); 597*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha); CHKERRQ(ierr); 598*c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 599*c8bcdf1eSAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk) - (step/tau_k)*gd/dk_yk); 600*c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d + t*tmp*yk */ 601*c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp*tau_k, beta, tao->gradient, cg->yk); CHKERRQ(ierr); 602*c8bcdf1eSAdam Denchfield } 603*c8bcdf1eSAdam Denchfield else if (snorm < cg->epsilon || cg->yts < cg->epsilon) { 604*c8bcdf1eSAdam Denchfield /* We're making no progress. Scaled gradient descent step and reset diagonal scaling */ 605*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 606*c8bcdf1eSAdam Denchfield } else { 607*c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 608*c8bcdf1eSAdam Denchfield /* Compute the invD vector */ 609*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeDiagScaling(tao, cg->yts, ynorm2); CHKERRQ(ierr); 610*c8bcdf1eSAdam Denchfield /* Apply the invD scaling to all my vectors */ 611*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->g_work, cg->invD, tao->gradient); CHKERRQ(ierr); 612*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->d_work, cg->invD, tao->stepdirection); CHKERRQ(ierr); 613*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->y_work, cg->invD, cg->yk); CHKERRQ(ierr); 614*c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 615*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk); CHKERRQ(ierr); 616*c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */ 617*c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 618*c8bcdf1eSAdam Denchfield /* tau_bfgs can be a parameter we tune via command line in future versions. For now, just setting to one. May instead put this inside the ComputeDiagonalScaling function... */ 619*c8bcdf1eSAdam Denchfield cg->tau_bfgs = 1.0; 620*c8bcdf1eSAdam Denchfield /* tau_k = -ytDy/(ytd)^2 * gd */ 621*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->y_work, &tau_k); 622*c8bcdf1eSAdam Denchfield tau_k = -tau_k*gd/(dk_yk*dk_yk); 623*c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */ 624*c8bcdf1eSAdam Denchfield beta = cg->tau_bfgs*gkp1_yk/dk_yk - step*tmp + tau_k; 625*c8bcdf1eSAdam Denchfield /* Do the update in two steps */ 626*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -cg->tau_bfgs, beta, cg->g_work); CHKERRQ(ierr); 627*c8bcdf1eSAdam Denchfield ierr = VecAXPY(tao->stepdirection, tmp, cg->y_work); CHKERRQ(ierr); 628*c8bcdf1eSAdam Denchfield }} 629*c8bcdf1eSAdam Denchfield break; 630*c8bcdf1eSAdam Denchfield case CG_SSML_DFP: 631*c8bcdf1eSAdam Denchfield 632*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 633*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient); CHKERRQ(ierr); 634*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution); CHKERRQ(ierr); 635*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm); CHKERRQ(ierr); 636*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->sk, NORM_2, &snorm); CHKERRQ(ierr); 637*c8bcdf1eSAdam Denchfield ynorm2 = ynorm*ynorm; 638*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk); CHKERRQ(ierr); 639*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk); CHKERRQ(ierr); 640*c8bcdf1eSAdam Denchfield if (ynorm < cg->epsilon){ 641*c8bcdf1eSAdam Denchfield /* The gradient hasn't changed, so we should stay in the same direction as before. Don't update it or anything.*/ 642*c8bcdf1eSAdam Denchfield if (snorm < cg->epsilon){ 643*c8bcdf1eSAdam Denchfield /* We're making no progress. Gradient descent step. Not scaled by tau_k since it's' is 0 or NaN in this case.*/ 644*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); 645*c8bcdf1eSAdam Denchfield } 646*c8bcdf1eSAdam Denchfield } else { 647*c8bcdf1eSAdam Denchfield 648*c8bcdf1eSAdam Denchfield tau_k = cg->dfp_scale*snorm*snorm/(step*dk_yk); 649*c8bcdf1eSAdam Denchfield tmp = tau_k*gkp1_yk/ynorm2; 650*c8bcdf1eSAdam Denchfield beta = -step*gd/dk_yk; 651*c8bcdf1eSAdam Denchfield 652*c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*d + tmp*yk */ 653*c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk); CHKERRQ(ierr); 654*c8bcdf1eSAdam Denchfield } 655*c8bcdf1eSAdam Denchfield break; 656*c8bcdf1eSAdam Denchfield case CG_SSML_BROYDEN: 657*c8bcdf1eSAdam Denchfield 658*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 659*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient); CHKERRQ(ierr); 660*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution); CHKERRQ(ierr); 661*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm); CHKERRQ(ierr); 662*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->sk, NORM_2, &snorm); CHKERRQ(ierr); 663*c8bcdf1eSAdam Denchfield ynorm2 = ynorm*ynorm; 664*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk); CHKERRQ(ierr); 665*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk); CHKERRQ(ierr); 666*c8bcdf1eSAdam Denchfield 667*c8bcdf1eSAdam Denchfield if (ynorm < cg->epsilon){ 668*c8bcdf1eSAdam Denchfield /* The gradient hasn't changed, so we should stay in the same direction as before. Don't update it or anything.*/ 669*c8bcdf1eSAdam Denchfield if (snorm < cg->epsilon){ 670*c8bcdf1eSAdam Denchfield /* We're making no progress. Gradient descent step.*/ 671*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); 672*c8bcdf1eSAdam Denchfield } 673*c8bcdf1eSAdam Denchfield } else { 674*c8bcdf1eSAdam Denchfield /* Instead of a regular convex combination, we will solve a quadratic formula. */ 675*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_bfgs, cg->bfgs_scale); 676*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_dfp, cg->dfp_scale); 677*c8bcdf1eSAdam Denchfield tau_k = cg->theta*tau_bfgs + (1.0-cg->theta)*tau_dfp; 678*c8bcdf1eSAdam Denchfield 679*c8bcdf1eSAdam Denchfield /* Used for the gradient */ 680*c8bcdf1eSAdam Denchfield /* If bfgs_scale = 1.0, it should reproduce the bfgs tau_bfgs. If bfgs_scale = 0.0, it should reproduce the tau_dfp scaling. Same with dfp_scale. */ 681*c8bcdf1eSAdam Denchfield tmp = cg->theta*tau_bfgs*gd/dk_yk + (1-cg->theta)*tau_dfp*gkp1_yk/ynorm2; 682*c8bcdf1eSAdam Denchfield beta = cg->theta*tau_bfgs*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk)) - step*gd/dk_yk; 683*c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*d + tmp*yk */ 684*c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk); CHKERRQ(ierr); 685*c8bcdf1eSAdam Denchfield } 686*c8bcdf1eSAdam Denchfield break; 687*c8bcdf1eSAdam Denchfield 688*c8bcdf1eSAdam Denchfield case CG_HagerZhang: 689*c8bcdf1eSAdam Denchfield /* Their 2006 paper. Comes from deleting the y_k term from SSML_BFGS, introducing a theta parameter, and using a cutoff for beta. See DK 2013 pg. 315 for a review of CG_DESCENT 5.3 */ 690*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 691*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 692*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient); CHKERRQ(ierr); 693*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution); CHKERRQ(ierr); 694*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm); CHKERRQ(ierr); 695*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->sk, NORM_2, &snorm); CHKERRQ(ierr); 696*c8bcdf1eSAdam Denchfield ynorm2 = ynorm*ynorm; 697*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk); CHKERRQ(ierr); 698*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->sk, &cg->yts); CHKERRQ(ierr); 699*c8bcdf1eSAdam Denchfield if (cg->use_dynamic_restart){ 700*c8bcdf1eSAdam Denchfield ierr = TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold); CHKERRQ(ierr); 701*c8bcdf1eSAdam Denchfield } 702*c8bcdf1eSAdam Denchfield if (ynorm2 < cg->epsilon){ 703*c8bcdf1eSAdam Denchfield /* The gradient hasn't changed, so we should stay in the same direction as before. Don't update it or anything.*/ 704*c8bcdf1eSAdam Denchfield if (snorm < cg->eps_23){ 705*c8bcdf1eSAdam Denchfield /* We're making no progress. Scaled gradient descent step.*/ 706*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 707*c8bcdf1eSAdam Denchfield } 708*c8bcdf1eSAdam Denchfield } else if (cg->dynamic_restart){ 709*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 710*c8bcdf1eSAdam Denchfield } else if (tao->niter % (cg->min_restart_num*dim) == 0 && cg->spaced_restart){ 711*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 712*c8bcdf1eSAdam Denchfield } else { 713*c8bcdf1eSAdam Denchfield if (!cg->diag_scaling){ 714*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk); CHKERRQ(ierr); 715*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha); CHKERRQ(ierr); 716*c8bcdf1eSAdam Denchfield /* Supplying cg->alpha = -1.0 will give the CG_DESCENT 5.3 special case of tau_k = 1.0 */ 717*c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 718*c8bcdf1eSAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk)); 719*c8bcdf1eSAdam Denchfield /* Bound beta as in CG_DESCENT 5.3, as implemented, with the third comparison from DK 2013 */ 720*c8bcdf1eSAdam Denchfield beta = PetscMax(PetscMax(beta, 0.4*tau_k*gd_old/(dnorm*dnorm)), 0.5*tau_k*gd/(dnorm*dnorm)); 721*c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d */ 722*c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, 0.0, beta, tao->gradient, cg->yk); CHKERRQ(ierr); 723*c8bcdf1eSAdam Denchfield } 724*c8bcdf1eSAdam Denchfield else if (snorm < cg->epsilon || cg->yts < cg->epsilon) { 725*c8bcdf1eSAdam Denchfield /* We're making no progress. Scaled gradient descent step and reset diagonal scaling */ 726*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 727*c8bcdf1eSAdam Denchfield } else { 728*c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 729*c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 730*c8bcdf1eSAdam Denchfield cg->sts = snorm*snorm; 731*c8bcdf1eSAdam Denchfield /* Compute the invD vector */ 732*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeDiagScaling(tao, cg->yts, ynorm2); CHKERRQ(ierr); 733*c8bcdf1eSAdam Denchfield /* Apply the invD scaling to all my vectors */ 734*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->g_work, cg->invD, tao->gradient); CHKERRQ(ierr); 735*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->d_work, cg->invD, tao->stepdirection); CHKERRQ(ierr); 736*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->y_work, cg->invD, cg->yk); CHKERRQ(ierr); 737*c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 738*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk); CHKERRQ(ierr); 739*c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */ 740*c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 741*c8bcdf1eSAdam Denchfield /* tau_bfgs can be a parameter we tune via command line in future versions. For now, just setting to one. May instead put this inside the ComputeDiagonalScaling function... */ 742*c8bcdf1eSAdam Denchfield cg->tau_bfgs = 1.0; 743*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->y_work, &tau_k); 744*c8bcdf1eSAdam Denchfield tau_k = -tau_k*gd/(dk_yk*dk_yk); 745*c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */ 746*c8bcdf1eSAdam Denchfield beta = cg->tau_bfgs*gkp1_yk/dk_yk + cg->hz_theta*tau_k; /* HZ; (1.15) from DK 2013 */ 747*c8bcdf1eSAdam Denchfield /* From HZ2013, modified to account for diagonal scaling*/ 748*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, cg->d_work, &gd_old); 749*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->stepdirection, cg->g_work, &gd); 750*c8bcdf1eSAdam Denchfield beta = PetscMax(PetscMax(beta, 0.4*gd_old/(dnorm*dnorm)), 0.5*gd/(dnorm*dnorm)); 751*c8bcdf1eSAdam Denchfield /* Do the update */ 752*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -cg->tau_bfgs, beta, cg->g_work); CHKERRQ(ierr); 753*c8bcdf1eSAdam Denchfield }} 754*c8bcdf1eSAdam Denchfield break; 755*c8bcdf1eSAdam Denchfield case CG_DaiKou: 756*c8bcdf1eSAdam Denchfield /* 2013 paper. */ 757*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 758*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 759*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient); CHKERRQ(ierr); 760*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution); CHKERRQ(ierr); 761*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm); CHKERRQ(ierr); 762*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->sk, NORM_2, &snorm); CHKERRQ(ierr); 763*c8bcdf1eSAdam Denchfield ynorm2 = ynorm*ynorm; 764*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk); CHKERRQ(ierr); 765*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->sk, &cg->yts); CHKERRQ(ierr); 766*c8bcdf1eSAdam Denchfield /* TODO: Should only need one of dk_yk and yts */ 767*c8bcdf1eSAdam Denchfield if (ynorm2 < cg->epsilon){ 768*c8bcdf1eSAdam Denchfield /* The gradient hasn't changed, so we should stay in the same direction as before. Don't update it or anything.*/ 769*c8bcdf1eSAdam Denchfield if (snorm < cg->eps_23){ 770*c8bcdf1eSAdam Denchfield /* We're making no progress. Scaled gradient descent step.*/ 771*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 772*c8bcdf1eSAdam Denchfield } 773*c8bcdf1eSAdam Denchfield } else { 774*c8bcdf1eSAdam Denchfield if (!cg->diag_scaling){ 775*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk); CHKERRQ(ierr); 776*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha); 777*c8bcdf1eSAdam Denchfield /* Use cg->alpha = -1.0 to get tau_k = 1.0 as in CG_DESCENT 5.3 */ 778*c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 779*c8bcdf1eSAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk) - (step/tau_k)*gd/dk_yk + gd/(dnorm*dnorm)); 780*c8bcdf1eSAdam Denchfield beta = PetscMax(PetscMax(beta, 0.4*tau_k*gd_old/(dnorm*dnorm)), 0.5*tau_k*gd/(dnorm*dnorm)); 781*c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d */ 782*c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, 0.0, beta, tao->gradient, cg->yk); CHKERRQ(ierr); 783*c8bcdf1eSAdam Denchfield } 784*c8bcdf1eSAdam Denchfield else if (snorm < cg->epsilon || cg->yts < cg->epsilon) { 785*c8bcdf1eSAdam Denchfield /* We're making no progress. Scaled gradient descent step and reset diagonal scaling */ 786*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 787*c8bcdf1eSAdam Denchfield } else { 788*c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 789*c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 790*c8bcdf1eSAdam Denchfield cg->sts = snorm*snorm; 791*c8bcdf1eSAdam Denchfield /* Compute the invD vector */ 792*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeDiagScaling(tao, cg->yts, ynorm2); CHKERRQ(ierr); 793*c8bcdf1eSAdam Denchfield /* Apply the invD scaling to all my vectors */ 794*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->g_work, cg->invD, tao->gradient); CHKERRQ(ierr); 795*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->d_work, cg->invD, tao->stepdirection); CHKERRQ(ierr); 796*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->y_work, cg->invD, cg->yk); CHKERRQ(ierr); 797*c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 798*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk); CHKERRQ(ierr); 799*c8bcdf1eSAdam Denchfield /* tau_bfgs can be a parameter we tune via command line in future versions. For now, just setting to one. May instead put this inside the ComputeDiagonalScaling function... */ 800*c8bcdf1eSAdam Denchfield cg->tau_bfgs = 1.0; 801*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->y_work, &tau_k); 802*c8bcdf1eSAdam Denchfield tau_k = tau_k*gd/(dk_yk*dk_yk); 803*c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 804*c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */ 805*c8bcdf1eSAdam Denchfield beta = cg->tau_bfgs*gkp1_yk/dk_yk - step*tmp - tau_k; 806*c8bcdf1eSAdam Denchfield 807*c8bcdf1eSAdam Denchfield /* Update this for the last term in beta */ 808*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->y_work, tao->stepdirection, &dk_yk); CHKERRQ(ierr); 809*c8bcdf1eSAdam Denchfield beta += tmp*dk_yk/(dnorm*dnorm); /* projection of y_work onto dk */ 810*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->stepdirection, cg->g_work, &gd); 811*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, cg->d_work, &gd_old); 812*c8bcdf1eSAdam Denchfield beta = PetscMax(PetscMax(beta, 0.4*gd_old/(dnorm*dnorm)), 0.5*gd/(dnorm*dnorm)); 813*c8bcdf1eSAdam Denchfield /* TODO: need to change these hardcoded constants into user parameters */ 814*c8bcdf1eSAdam Denchfield /* Do the update */ 815*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -cg->tau_bfgs, beta, cg->g_work); CHKERRQ(ierr); 816*c8bcdf1eSAdam Denchfield }} 817*c8bcdf1eSAdam Denchfield break; 818*c8bcdf1eSAdam Denchfield 819*c8bcdf1eSAdam Denchfield case CG_KouDai: 820*c8bcdf1eSAdam Denchfield 821*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 822*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 823*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient); CHKERRQ(ierr); 824*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution); CHKERRQ(ierr); 825*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm); CHKERRQ(ierr); 826*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->sk, NORM_2, &snorm); CHKERRQ(ierr); 827*c8bcdf1eSAdam Denchfield ynorm2 = ynorm*ynorm; 828*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk); CHKERRQ(ierr); 829*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->sk, &cg->yts); CHKERRQ(ierr); 830*c8bcdf1eSAdam Denchfield ierr = VecGetSize(tao->gradient, &dim); CHKERRQ(ierr); 831*c8bcdf1eSAdam Denchfield /* TODO: Should only need one of dk_yk and yts */ 832*c8bcdf1eSAdam Denchfield if (cg->use_dynamic_restart){ 833*c8bcdf1eSAdam Denchfield ierr = TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold); CHKERRQ(ierr); 834*c8bcdf1eSAdam Denchfield } 835*c8bcdf1eSAdam Denchfield if (ynorm2 < cg->epsilon){ 836*c8bcdf1eSAdam Denchfield /* The gradient hasn't changed, so we should stay in the same direction as before. Don't update it or anything.*/ 837*c8bcdf1eSAdam Denchfield if (snorm < cg->eps_23){ 838*c8bcdf1eSAdam Denchfield /* We're making no progress. Scaled gradient descent step.*/ 839*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 840*c8bcdf1eSAdam Denchfield } 841*c8bcdf1eSAdam Denchfield } else if (cg->dynamic_restart){ 842*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 843*c8bcdf1eSAdam Denchfield } else if (tao->niter % (cg->min_restart_num*dim) == 0 && cg->spaced_restart){ 844*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 845*c8bcdf1eSAdam Denchfield } else { 846*c8bcdf1eSAdam Denchfield if (!cg->diag_scaling){ 847*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk); CHKERRQ(ierr); 848*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha); CHKERRQ(ierr); 849*c8bcdf1eSAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk)) - step*gd/dk_yk; 850*c8bcdf1eSAdam Denchfield if (beta < cg->zeta*tau_k*gd/(dnorm*dnorm)) /* 0.1 is KD's zeta parameter */ 851*c8bcdf1eSAdam Denchfield { 852*c8bcdf1eSAdam Denchfield beta = cg->zeta*tau_k*gd/(dnorm*dnorm); 853*c8bcdf1eSAdam Denchfield gamma = 0.0; 854*c8bcdf1eSAdam Denchfield } else { 855*c8bcdf1eSAdam Denchfield if (gkp1_yk < 0 && cg->neg_xi) gamma = -1.0*gd/dk_yk; 856*c8bcdf1eSAdam Denchfield /* This seems to be very effective when there's no tau_k scaling... this guarantees a large descent step every iteration, going through DK 2015 Lemma 3.1's proof but allowing for negative xi */ 857*c8bcdf1eSAdam Denchfield else gamma = cg->xi*gd/dk_yk; 858*c8bcdf1eSAdam Denchfield } 859*c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d + t*tmp*yk */ 860*c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, gamma*tau_k, beta, tao->gradient, cg->yk); CHKERRQ(ierr); 861*c8bcdf1eSAdam Denchfield } 862*c8bcdf1eSAdam Denchfield else if (snorm < cg->epsilon || cg->yts < cg->epsilon) { 863*c8bcdf1eSAdam Denchfield /* We're making no progress. Scaled gradient descent step and reset diagonal scaling */ 864*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 865*c8bcdf1eSAdam Denchfield } else { 866*c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 867*c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 868*c8bcdf1eSAdam Denchfield cg->sts = snorm*snorm; 869*c8bcdf1eSAdam Denchfield /* Compute the invD vector */ 870*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeDiagScaling(tao, cg->yts, ynorm2); CHKERRQ(ierr); 871*c8bcdf1eSAdam Denchfield /* Apply the invD scaling to all my vectors */ 872*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->g_work, cg->invD, tao->gradient); CHKERRQ(ierr); 873*c8bcdf1eSAdam Denchfield /* ierr = VecPointwiseMult(cg->d_work, cg->invD, tao->stepdirection); CHKERRQ(ierr); */ 874*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->y_work, cg->invD, cg->yk); CHKERRQ(ierr); 875*c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 876*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1D_yk); CHKERRQ(ierr); 877*c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */ 878*c8bcdf1eSAdam Denchfield gamma = gd/dk_yk; 879*c8bcdf1eSAdam Denchfield /* tau_k = -ytDy/(ytd)^2 * gd */ 880*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->y_work, &tau_k); 881*c8bcdf1eSAdam Denchfield tau_k = tau_k*gd/(dk_yk*dk_yk); 882*c8bcdf1eSAdam Denchfield /* beta is the constant which adds the d_k contribution */ 883*c8bcdf1eSAdam Denchfield beta = gkp1D_yk/dk_yk - step*gamma - tau_k; 884*c8bcdf1eSAdam Denchfield /* Here is the requisite check */ 885*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->stepdirection, cg->g_work, &tmp); 886*c8bcdf1eSAdam Denchfield if (cg->neg_xi){ 887*c8bcdf1eSAdam Denchfield /* modified KD implementation */ 888*c8bcdf1eSAdam Denchfield if (gkp1D_yk/dk_yk < 0) gamma = -1.0*gd/dk_yk; 889*c8bcdf1eSAdam Denchfield else gamma = cg->xi*gd/dk_yk; 890*c8bcdf1eSAdam Denchfield if (beta < cg->zeta*tmp/(dnorm*dnorm)){ 891*c8bcdf1eSAdam Denchfield beta = cg->zeta*tmp/(dnorm*dnorm); 892*c8bcdf1eSAdam Denchfield gamma = 0.0; 893*c8bcdf1eSAdam Denchfield } 894*c8bcdf1eSAdam Denchfield } else { /* original KD 2015 implementation */ 895*c8bcdf1eSAdam Denchfield if (beta < cg->zeta*tmp/(dnorm*dnorm)) { 896*c8bcdf1eSAdam Denchfield beta = cg->zeta*tmp/(dnorm*dnorm); 897*c8bcdf1eSAdam Denchfield gamma = 0.0; 898*c8bcdf1eSAdam Denchfield } else { 899*c8bcdf1eSAdam Denchfield gamma = cg->xi*gd/dk_yk; 900*c8bcdf1eSAdam Denchfield } 901*c8bcdf1eSAdam Denchfield } 902*c8bcdf1eSAdam Denchfield /* Do the update in two steps */ 903*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work); CHKERRQ(ierr); 904*c8bcdf1eSAdam Denchfield ierr = VecAXPY(tao->stepdirection, gamma, cg->y_work); CHKERRQ(ierr); 905*c8bcdf1eSAdam Denchfield }} 906*c8bcdf1eSAdam Denchfield break; 907*c8bcdf1eSAdam Denchfield case CG_BFGS_Mod: 908*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 909*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient); CHKERRQ(ierr); 910*c8bcdf1eSAdam Denchfield ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution); CHKERRQ(ierr); 911*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->yk, NORM_2, &ynorm); CHKERRQ(ierr); 912*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->sk, NORM_2, &snorm); CHKERRQ(ierr); 913*c8bcdf1eSAdam Denchfield ynorm2 = ynorm*ynorm; 914*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk); CHKERRQ(ierr); 915*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->sk, &cg->yts); CHKERRQ(ierr); 916*c8bcdf1eSAdam Denchfield /* TODO: Should only need one of dk_yk and yts */ 917*c8bcdf1eSAdam Denchfield if (ynorm2 < cg->epsilon){ 918*c8bcdf1eSAdam Denchfield /* The gradient hasn't changed, so we should stay in the same direction as before. Don't update it or anything.*/ 919*c8bcdf1eSAdam Denchfield if (snorm < cg->eps_23){ 920*c8bcdf1eSAdam Denchfield /* We're making no progress. Scaled gradient descent step.*/ 921*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 922*c8bcdf1eSAdam Denchfield } 923*c8bcdf1eSAdam Denchfield } else { 924*c8bcdf1eSAdam Denchfield if (!cg->diag_scaling){ 925*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk); CHKERRQ(ierr); 926*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha); CHKERRQ(ierr); 927*c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 928*c8bcdf1eSAdam Denchfield beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk) - (step/tau_k)*gd/dk_yk); 929*c8bcdf1eSAdam Denchfield /* d <- -t*g + beta*t*d + t*tmp*yk */ 930*c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp*tau_k, beta, tao->gradient, cg->yk); CHKERRQ(ierr); 931*c8bcdf1eSAdam Denchfield } 932*c8bcdf1eSAdam Denchfield else if (snorm < cg->epsilon || cg->yts < cg->epsilon) { 933*c8bcdf1eSAdam Denchfield /* We're making no progress. Scaled gradient descent step and reset diagonal scaling */ 934*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); CHKERRQ(ierr); 935*c8bcdf1eSAdam Denchfield } else { 936*c8bcdf1eSAdam Denchfield /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 937*c8bcdf1eSAdam Denchfield cg->yty = ynorm2; 938*c8bcdf1eSAdam Denchfield cg->sts = snorm*snorm; 939*c8bcdf1eSAdam Denchfield /* Compute the invD vector */ 940*c8bcdf1eSAdam Denchfield ierr = TaoBNCGComputeDiagScaling(tao, cg->yts, ynorm2); CHKERRQ(ierr); 941*c8bcdf1eSAdam Denchfield /* Apply the invD scaling to all my vectors */ 942*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->g_work, cg->invD, tao->gradient); CHKERRQ(ierr); 943*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->d_work, cg->invD, tao->stepdirection); CHKERRQ(ierr); 944*c8bcdf1eSAdam Denchfield ierr = VecPointwiseMult(cg->y_work, cg->invD, cg->yk); CHKERRQ(ierr); 945*c8bcdf1eSAdam Denchfield /* Construct the constant ytDgkp1 */ 946*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk); CHKERRQ(ierr); 947*c8bcdf1eSAdam Denchfield /* Construct the constant for scaling Dkyk in the update */ 948*c8bcdf1eSAdam Denchfield tmp = gd/dk_yk; 949*c8bcdf1eSAdam Denchfield /* tau_bfgs can be a parameter we tune via command line in future versions. For now, just setting to one. May instead put this inside the ComputeDiagonalScaling function... */ 950*c8bcdf1eSAdam Denchfield cg->tau_bfgs = 1.0; 951*c8bcdf1eSAdam Denchfield /* tau_k = -ytDy/(ytd)^2 * gd */ 952*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->yk, cg->y_work, &tau_k); 953*c8bcdf1eSAdam Denchfield tau_k = -tau_k*gd/(dk_yk*dk_yk); 954*c8bcdf1eSAdam Denchfield /* beta is the constant which adds the dk contribution */ 955*c8bcdf1eSAdam Denchfield beta = cg->tau_bfgs*gkp1_yk/dk_yk - step*tmp + tau_k; 956*c8bcdf1eSAdam Denchfield /* Do the update in two steps */ 957*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(tao->stepdirection, -cg->tau_bfgs, beta, cg->g_work); CHKERRQ(ierr); 958*c8bcdf1eSAdam Denchfield ierr = VecAXPY(tao->stepdirection, tmp, cg->y_work); CHKERRQ(ierr); 959*c8bcdf1eSAdam Denchfield }} 960*c8bcdf1eSAdam Denchfield default: 961*c8bcdf1eSAdam Denchfield beta = 0.0; 962*c8bcdf1eSAdam Denchfield break; 963*c8bcdf1eSAdam Denchfield } 964*c8bcdf1eSAdam Denchfield } 965*c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 966*c8bcdf1eSAdam Denchfield } 967*c8bcdf1eSAdam Denchfield 968*c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGSteffensonAcceleration(Tao tao){ 969*c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 970*c8bcdf1eSAdam Denchfield PetscErrorCode ierr; 971*c8bcdf1eSAdam Denchfield PetscReal mag1, mag2; 972*c8bcdf1eSAdam Denchfield PetscReal resnorm; 973*c8bcdf1eSAdam Denchfield PetscReal steff_f; 974*c8bcdf1eSAdam Denchfield PetscFunctionBegin; 975*c8bcdf1eSAdam Denchfield if (tao->niter > 2 && tao->niter % 2 == 0){ 976*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->steffnp1, cg->steffnm1); CHKERRQ(ierr); /* X_np1 to X_nm1 since it's been two iterations*/ 977*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->X_old, cg->steffn); CHKERRQ(ierr); /* Get X_n */ 978*c8bcdf1eSAdam Denchfield ierr = VecCopy(tao->solution, cg->steffnp1); CHKERRQ(ierr); 979*c8bcdf1eSAdam Denchfield 980*c8bcdf1eSAdam Denchfield /* Begin step 4 */ 981*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->steffnm1, cg->W); CHKERRQ(ierr); 982*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(cg->W, 1.0, -1.0, cg->steffn); CHKERRQ(ierr); 983*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->W, cg->W, &mag1); 984*c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(cg->W, -1.0, 1.0, -1.0, cg->steffn, cg->steffnp1); CHKERRQ(ierr); 985*c8bcdf1eSAdam Denchfield ierr = VecDot(cg->W, cg->W, &mag2); 986*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->steffnm1, cg->steffva); CHKERRQ(ierr); 987*c8bcdf1eSAdam Denchfield ierr = VecAXPY(cg->steffva, -mag1/mag2, cg->W); CHKERRQ(ierr); 988*c8bcdf1eSAdam Denchfield 989*c8bcdf1eSAdam Denchfield ierr = TaoComputeObjectiveAndGradient(tao, cg->steffva, &steff_f, cg->g_work); CHKERRQ(ierr); 990*c8bcdf1eSAdam Denchfield 991*c8bcdf1eSAdam Denchfield /* Check if the accelerated point has converged*/ 992*c8bcdf1eSAdam Denchfield ierr = VecFischer(cg->steffva, cg->g_work, tao->XL, tao->XU, cg->W);CHKERRQ(ierr); 993*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr); 994*c8bcdf1eSAdam Denchfield if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 995*c8bcdf1eSAdam Denchfield //ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 996*c8bcdf1eSAdam Denchfield //ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr); 997*c8bcdf1eSAdam Denchfield ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 998*c8bcdf1eSAdam Denchfield ierr = PetscPrintf(PETSC_COMM_SELF, "va f: %20.19e\t va gnorm: %20.19e\n", steff_f, resnorm); 999*c8bcdf1eSAdam Denchfield 1000*c8bcdf1eSAdam Denchfield } else if (tao->niter == 2){ 1001*c8bcdf1eSAdam Denchfield ierr = VecCopy(tao->solution, cg->steffnp1); CHKERRQ(ierr); 1002*c8bcdf1eSAdam Denchfield mag1 = cg->sts; /* = |x1 - x0|^2 */ 1003*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->steffnm1, cg->steffvatmp); CHKERRQ(ierr); 1004*c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(cg->steffva, 1.0, -1.0, 0.0, cg->steffnp1, cg->steffn); 1005*c8bcdf1eSAdam Denchfield ierr = VecAXPBYPCZ(cg->steffva, 1.0, -1.0, 1.0, cg->steffnm1, cg->steffn); 1006*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->steffva, NORM_2, &mag2); CHKERRQ(ierr); 1007*c8bcdf1eSAdam Denchfield mag2 = mag2*mag2; 1008*c8bcdf1eSAdam Denchfield ierr = VecAXPBY(cg->steffva, 1.0, -mag1/mag2, cg->steffvatmp); CHKERRQ(ierr); 1009*c8bcdf1eSAdam Denchfield // finished step 2 1010*c8bcdf1eSAdam Denchfield } else if (tao->niter == 1){ 1011*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->X_old, cg->steffnm1); CHKERRQ(ierr); 1012*c8bcdf1eSAdam Denchfield ierr = VecCopy(tao->solution, cg->steffn); CHKERRQ(ierr); 1013*c8bcdf1eSAdam Denchfield } 1014*c8bcdf1eSAdam Denchfield /* Now have step 2 done of method */ 1015*c8bcdf1eSAdam Denchfield 1016*c8bcdf1eSAdam Denchfield PetscFunctionReturn(0); 1017*c8bcdf1eSAdam Denchfield } 1018*c8bcdf1eSAdam Denchfield 1019*c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao tao, PetscReal gnorm) 1020*c8bcdf1eSAdam Denchfield { 1021*c8bcdf1eSAdam Denchfield TAO_BNCG *cg = (TAO_BNCG*)tao->data; 1022*c8bcdf1eSAdam Denchfield PetscErrorCode ierr; 1023*c8bcdf1eSAdam Denchfield TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 1024*c8bcdf1eSAdam Denchfield PetscReal step=1.0,gnorm2,gd,ginner,dnorm; 1025*c8bcdf1eSAdam Denchfield PetscReal gnorm2_old,f_old,resnorm, gnorm_old; 1026*c8bcdf1eSAdam Denchfield PetscBool cg_restart, gd_fallback = PETSC_FALSE; 1027*c8bcdf1eSAdam Denchfield 1028*c8bcdf1eSAdam Denchfield PetscFunctionBegin; 1029*c8bcdf1eSAdam Denchfield /* We are now going to perform a line search along the direction. */ 1030*c8bcdf1eSAdam Denchfield ++tao->niter; 1031*c8bcdf1eSAdam Denchfield 1032*c8bcdf1eSAdam Denchfield /* Store solution and gradient info before it changes */ 1033*c8bcdf1eSAdam Denchfield ierr = VecCopy(tao->solution, cg->X_old);CHKERRQ(ierr); 1034*c8bcdf1eSAdam Denchfield ierr = VecCopy(tao->gradient, cg->G_old);CHKERRQ(ierr); 1035*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->unprojected_gradient, cg->unprojected_gradient_old);CHKERRQ(ierr); 1036*c8bcdf1eSAdam Denchfield 1037*c8bcdf1eSAdam Denchfield gnorm_old = gnorm; 1038*c8bcdf1eSAdam Denchfield gnorm2_old = gnorm_old*gnorm_old; 1039*c8bcdf1eSAdam Denchfield f_old = cg->f; 1040*c8bcdf1eSAdam Denchfield /* Perform bounded line search */ 1041*c8bcdf1eSAdam Denchfield ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 1042*c8bcdf1eSAdam Denchfield ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 1043*c8bcdf1eSAdam Denchfield ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 1044*c8bcdf1eSAdam Denchfield 1045*c8bcdf1eSAdam Denchfield /* Check linesearch failure */ 1046*c8bcdf1eSAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 1047*c8bcdf1eSAdam Denchfield ++cg->ls_fails; 1048*c8bcdf1eSAdam Denchfield 1049*c8bcdf1eSAdam Denchfield if (cg->cg_type == CG_GradientDescent || gd_fallback){ 1050*c8bcdf1eSAdam Denchfield /* Nothing left to do but fail out of the optimization */ 1051*c8bcdf1eSAdam Denchfield step = 0.0; 1052*c8bcdf1eSAdam Denchfield tao->reason = TAO_DIVERGED_LS_FAILURE; 1053*c8bcdf1eSAdam Denchfield } else { 1054*c8bcdf1eSAdam Denchfield /* Restore previous point */ 1055*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->X_old, tao->solution);CHKERRQ(ierr); 1056*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->G_old, tao->gradient);CHKERRQ(ierr); 1057*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->unprojected_gradient_old, cg->unprojected_gradient);CHKERRQ(ierr); 1058*c8bcdf1eSAdam Denchfield 1059*c8bcdf1eSAdam Denchfield gnorm = gnorm_old; 1060*c8bcdf1eSAdam Denchfield gnorm2 = gnorm2_old; 1061*c8bcdf1eSAdam Denchfield cg->f = f_old; 1062*c8bcdf1eSAdam Denchfield 1063*c8bcdf1eSAdam Denchfield /* Fall back on the scaled gradient step */ 1064*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); 1065*c8bcdf1eSAdam Denchfield ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1066*c8bcdf1eSAdam Denchfield 1067*c8bcdf1eSAdam Denchfield ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 1068*c8bcdf1eSAdam Denchfield ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 1069*c8bcdf1eSAdam Denchfield ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 1070*c8bcdf1eSAdam Denchfield 1071*c8bcdf1eSAdam Denchfield if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER){ 1072*c8bcdf1eSAdam Denchfield /* Nothing left to do but fail out of the optimization */ 1073*c8bcdf1eSAdam Denchfield cg->ls_fails++; 1074*c8bcdf1eSAdam Denchfield step = 0.0; 1075*c8bcdf1eSAdam Denchfield tao->reason = TAO_DIVERGED_LS_FAILURE; 1076*c8bcdf1eSAdam Denchfield } 1077*c8bcdf1eSAdam Denchfield } 1078*c8bcdf1eSAdam Denchfield } 1079*c8bcdf1eSAdam Denchfield /* Convergence test for line search failure */ 1080*c8bcdf1eSAdam Denchfield if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 1081*c8bcdf1eSAdam Denchfield 1082*c8bcdf1eSAdam Denchfield /* Standard convergence test */ 1083*c8bcdf1eSAdam Denchfield /* Make sure convergence test is the same. */ 1084*c8bcdf1eSAdam Denchfield ierr = VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W);CHKERRQ(ierr); 1085*c8bcdf1eSAdam Denchfield ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr); 1086*c8bcdf1eSAdam Denchfield if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 1087*c8bcdf1eSAdam Denchfield ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 1088*c8bcdf1eSAdam Denchfield ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr); 1089*c8bcdf1eSAdam Denchfield ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 1090*c8bcdf1eSAdam Denchfield if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 1091*c8bcdf1eSAdam Denchfield 1092*c8bcdf1eSAdam Denchfield /* Assert we have an updated step and we need at least one more iteration. */ 1093*c8bcdf1eSAdam Denchfield /* Calculate the next direction */ 1094*c8bcdf1eSAdam Denchfield /* Estimate the active set at the new solution */ 1095*c8bcdf1eSAdam Denchfield ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr); 1096*c8bcdf1eSAdam Denchfield /* Compute the projected gradient and its norm */ 1097*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 1098*c8bcdf1eSAdam Denchfield ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr); 1099*c8bcdf1eSAdam Denchfield ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 1100*c8bcdf1eSAdam Denchfield gnorm2 = gnorm*gnorm; 1101*c8bcdf1eSAdam Denchfield 1102*c8bcdf1eSAdam Denchfield /* Check restart conditions for using steepest descent */ 1103*c8bcdf1eSAdam Denchfield cg_restart = PETSC_FALSE; 1104*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, cg->G_old, &ginner);CHKERRQ(ierr); 1105*c8bcdf1eSAdam Denchfield ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm);CHKERRQ(ierr); 1106*c8bcdf1eSAdam Denchfield 1107*c8bcdf1eSAdam Denchfield ierr = TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, cg_restart, dnorm, ginner); CHKERRQ(ierr); 1108*c8bcdf1eSAdam Denchfield 1109*c8bcdf1eSAdam Denchfield ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1110*c8bcdf1eSAdam Denchfield 1111*c8bcdf1eSAdam Denchfield if (cg->cg_type != CG_GradientDescent) { 1112*c8bcdf1eSAdam Denchfield /* Figure out which previously active variables became inactive this iteration */ 1113*c8bcdf1eSAdam Denchfield ierr = ISDestroy(&cg->new_inactives);CHKERRQ(ierr); 1114*c8bcdf1eSAdam Denchfield if (cg->inactive_idx && cg->inactive_old) { 1115*c8bcdf1eSAdam Denchfield ierr = ISDifference(cg->inactive_idx, cg->inactive_old, &cg->new_inactives);CHKERRQ(ierr); 1116*c8bcdf1eSAdam Denchfield } 1117*c8bcdf1eSAdam Denchfield /* Selectively reset the CG step those freshly inactive variables */ 1118*c8bcdf1eSAdam Denchfield if (cg->new_inactives) { 1119*c8bcdf1eSAdam Denchfield ierr = VecGetSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr); 1120*c8bcdf1eSAdam Denchfield ierr = VecGetSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr); 1121*c8bcdf1eSAdam Denchfield ierr = VecCopy(cg->inactive_grad, cg->inactive_step);CHKERRQ(ierr); 1122*c8bcdf1eSAdam Denchfield ierr = VecScale(cg->inactive_step, -1.0);CHKERRQ(ierr); 1123*c8bcdf1eSAdam Denchfield ierr = VecRestoreSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr); 1124*c8bcdf1eSAdam Denchfield ierr = VecRestoreSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr); 1125*c8bcdf1eSAdam Denchfield } 1126*c8bcdf1eSAdam Denchfield 1127*c8bcdf1eSAdam Denchfield /* Verify that this is a descent direction */ 1128*c8bcdf1eSAdam Denchfield ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 1129*c8bcdf1eSAdam Denchfield ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm); 1130*c8bcdf1eSAdam Denchfield /* TODO: potentially remove the third check */ 1131*c8bcdf1eSAdam Denchfield if (gd >= 0 || PetscIsInfOrNanReal(gd) || gd >= -cg->epsilon) { 1132*c8bcdf1eSAdam Denchfield /* Not a descent direction, so we reset back to projected gradient descent */ 1133*c8bcdf1eSAdam Denchfield PetscPrintf(PETSC_COMM_SELF, "gd is small or positive: %20.19e\n", gd); 1134*c8bcdf1eSAdam Denchfield ierr = TaoBNCGResetUpdate(tao, gnorm2); 1135*c8bcdf1eSAdam Denchfield ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1136*c8bcdf1eSAdam Denchfield ++cg->descent_error; 1137*c8bcdf1eSAdam Denchfield gd_fallback = PETSC_TRUE; 1138*c8bcdf1eSAdam Denchfield } else { 1139*c8bcdf1eSAdam Denchfield gd_fallback = PETSC_FALSE; 1140*c8bcdf1eSAdam Denchfield } 1141*c8bcdf1eSAdam Denchfield } 1142*c8bcdf1eSAdam Denchfield 1143ac9112b8SAlp Dener PetscFunctionReturn(0); 1144ac9112b8SAlp Dener } 1145