1fed79b8eSAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h> 2fed79b8eSAlp Dener #include <petscksp.h> 3fed79b8eSAlp Dener 4fed79b8eSAlp Dener /* 5fed79b8eSAlp Dener Implements Newton's Method with a trust region approach for solving 6fed79b8eSAlp Dener bound constrained minimization problems. 7fed79b8eSAlp Dener 8198282dbSAlp Dener ------------------------------------------------------------ 9198282dbSAlp Dener 10198282dbSAlp Dener x_0 = VecMedian(x_0) 11198282dbSAlp Dener f_0, g_0= TaoComputeObjectiveAndGradient(x_0) 12c4b75bccSAlp Dener pg_0 = project(g_0) 13198282dbSAlp Dener check convergence at pg_0 14c4b75bccSAlp Dener needH = TaoBNKInitialize(default:BNK_INIT_INTERPOLATION) 15198282dbSAlp Dener niter = 0 16c4b75bccSAlp Dener step_accepted = false 17198282dbSAlp Dener 18198282dbSAlp Dener while niter <= max_it 19198282dbSAlp Dener niter += 1 20c4b75bccSAlp Dener 21c4b75bccSAlp Dener if needH 22c4b75bccSAlp Dener If max_cg_steps > 0 23c4b75bccSAlp Dener x_k, g_k, pg_k = TaoSolve(BNCG) 24c4b75bccSAlp Dener end 25c4b75bccSAlp Dener 26198282dbSAlp Dener H_k = TaoComputeHessian(x_k) 27198282dbSAlp Dener if pc_type == BNK_PC_BFGS 28198282dbSAlp Dener add correction to BFGS approx 29198282dbSAlp Dener if scale_type == BNK_SCALE_AHESS 30198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(H_k)), 1e6) 31198282dbSAlp Dener scale BFGS with VecReciprocal(D) 32198282dbSAlp Dener end 33198282dbSAlp Dener end 34c4b75bccSAlp Dener needH = False 35198282dbSAlp Dener end 36198282dbSAlp Dener 37198282dbSAlp Dener if pc_type = BNK_PC_BFGS 38198282dbSAlp Dener B_k = BFGS 39198282dbSAlp Dener else 40198282dbSAlp Dener B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6) 41198282dbSAlp Dener B_k = VecReciprocal(B_k) 42198282dbSAlp Dener end 43198282dbSAlp Dener w = x_k - VecMedian(x_k - 0.001*B_k*g_k) 44198282dbSAlp Dener eps = min(eps, norm2(w)) 45198282dbSAlp Dener determine the active and inactive index sets such that 46198282dbSAlp Dener L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0} 47198282dbSAlp Dener U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0} 48198282dbSAlp Dener F = {i : l_i = (x_k)_i = u_i} 49198282dbSAlp Dener A = {L + U + F} 50c4b75bccSAlp Dener IA = {i : i not in A} 51198282dbSAlp Dener 52c4b75bccSAlp Dener generate the reduced system Hr_k dr_k = -gr_k for variables in IA 53198282dbSAlp Dener if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS 54198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6) 55198282dbSAlp Dener scale BFGS with VecReciprocal(D) 56198282dbSAlp Dener end 57c4b75bccSAlp Dener 58c4b75bccSAlp Dener while !stepAccepted 59198282dbSAlp Dener solve Hr_k dr_k = -gr_k 60198282dbSAlp Dener set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F 61198282dbSAlp Dener 62198282dbSAlp Dener x_{k+1} = VecMedian(x_k + d_k) 63198282dbSAlp Dener s = x_{k+1} - x_k 64198282dbSAlp Dener prered = dot(s, 0.5*gr_k - Hr_k*s) 65198282dbSAlp Dener f_{k+1} = TaoComputeObjective(x_{k+1}) 66198282dbSAlp Dener actred = f_k - f_{k+1} 67198282dbSAlp Dener 68198282dbSAlp Dener oldTrust = trust 69198282dbSAlp Dener step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION) 70198282dbSAlp Dener if step_accepted 71198282dbSAlp Dener g_{k+1} = TaoComputeGradient(x_{k+1}) 72c4b75bccSAlp Dener pg_{k+1} = project(g_{k+1}) 73198282dbSAlp Dener count the accepted Newton step 74c4b75bccSAlp Dener needH = True 75198282dbSAlp Dener else 76198282dbSAlp Dener f_{k+1} = f_k 77198282dbSAlp Dener x_{k+1} = x_k 78198282dbSAlp Dener g_{k+1} = g_k 79198282dbSAlp Dener pg_{k+1} = pg_k 80198282dbSAlp Dener if trust == oldTrust 81198282dbSAlp Dener terminate because we cannot shrink the radius any further 82198282dbSAlp Dener end 83198282dbSAlp Dener end 84198282dbSAlp Dener 85198282dbSAlp Dener check convergence at pg_{k+1} 86198282dbSAlp Dener end 87c4b75bccSAlp Dener 88c4b75bccSAlp Dener end 89fed79b8eSAlp Dener */ 90fed79b8eSAlp Dener 91e0ed867bSAlp Dener PetscErrorCode TaoSolve_BNTR(Tao tao) 92fed79b8eSAlp Dener { 93fed79b8eSAlp Dener PetscErrorCode ierr; 94fed79b8eSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 95e465cd6fSAlp Dener KSPConvergedReason ksp_reason; 96fed79b8eSAlp Dener 9789da521bSAlp Dener PetscReal oldTrust, prered, actred, steplen, resnorm; 98937a31a1SAlp Dener PetscBool cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE; 996b591159SAlp Dener PetscInt stepType, nDiff; 100fed79b8eSAlp Dener 101fed79b8eSAlp Dener PetscFunctionBegin; 10228017e9fSAlp Dener /* Initialize the preconditioner, KSP solver and trust radius/line search */ 103fed79b8eSAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 104937a31a1SAlp Dener ierr = TaoBNKInitialize(tao, bnk->init_type, &needH);CHKERRQ(ierr); 10528017e9fSAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 106fed79b8eSAlp Dener 107fed79b8eSAlp Dener /* Have not converged; continue with Newton method */ 108fed79b8eSAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 109e1e80dc8SAlp Dener /* Call general purpose update function */ 110e1e80dc8SAlp Dener if (tao->ops->update) { 111*8fcddce6SStefano Zampini ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr); 112e1e80dc8SAlp Dener } 113c4b75bccSAlp Dener ++tao->niter; 114e031d6f5SAlp Dener 11589da521bSAlp Dener if (needH && bnk->inactive_idx) { 116e031d6f5SAlp Dener /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */ 117e031d6f5SAlp Dener ierr = TaoBNKTakeCGSteps(tao, &cgTerminate);CHKERRQ(ierr); 118e031d6f5SAlp Dener if (cgTerminate) { 119e031d6f5SAlp Dener tao->reason = bnk->bncg->reason; 120e031d6f5SAlp Dener PetscFunctionReturn(0); 121fed79b8eSAlp Dener } 122937a31a1SAlp Dener /* Compute the hessian and update the BFGS preconditioner at the new iterate */ 123f7bf01afSAlp Dener ierr = (*bnk->computehessian)(tao);CHKERRQ(ierr); 124937a31a1SAlp Dener needH = PETSC_FALSE; 125937a31a1SAlp Dener } 126fed79b8eSAlp Dener 127fed79b8eSAlp Dener /* Store current solution before it changes */ 128fed79b8eSAlp Dener bnk->fold = bnk->f; 129fed79b8eSAlp Dener ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 130fed79b8eSAlp Dener ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr); 131fed79b8eSAlp Dener ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr); 132fed79b8eSAlp Dener 133937a31a1SAlp Dener /* Enter into trust region loops */ 134937a31a1SAlp Dener stepAccepted = PETSC_FALSE; 135937a31a1SAlp Dener while (!stepAccepted && tao->reason == TAO_CONTINUE_ITERATING) { 136937a31a1SAlp Dener tao->ksp_its=0; 137937a31a1SAlp Dener 138937a31a1SAlp Dener /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */ 139f7bf01afSAlp Dener ierr = (*bnk->computestep)(tao, shift, &ksp_reason, &stepType);CHKERRQ(ierr); 140937a31a1SAlp Dener 141b1c2d0e3SAlp Dener /* Temporarily accept the step and project it into the bounds */ 142fed79b8eSAlp Dener ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr); 1433b063059SAlp Dener ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr); 144b1c2d0e3SAlp Dener 145b1c2d0e3SAlp Dener /* Check if the projection changed the step direction */ 146c4b75bccSAlp Dener if (nDiff > 0) { 147c4b75bccSAlp Dener /* Projection changed the step, so we have to recompute the step and 148c4b75bccSAlp Dener the predicted reduction. Leave the trust radius unchanged. */ 149b1c2d0e3SAlp Dener ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr); 1508d5ead36SAlp Dener ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr); 1515e9b73cbSAlp Dener ierr = TaoBNKRecomputePred(tao, tao->stepdirection, &prered);CHKERRQ(ierr); 152b1c2d0e3SAlp Dener } else { 153b1c2d0e3SAlp Dener /* Step did not change, so we can just recover the pre-computed prediction */ 154b1c2d0e3SAlp Dener ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr); 155b1c2d0e3SAlp Dener } 156b1c2d0e3SAlp Dener prered = -prered; 157b1c2d0e3SAlp Dener 158b1c2d0e3SAlp Dener /* Compute the actual reduction and update the trust radius */ 159fed79b8eSAlp Dener ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr); 160b4a30f08SAlp Dener if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 161b1c2d0e3SAlp Dener actred = bnk->fold - bnk->f; 162e761ccfdSAlp Dener oldTrust = tao->trust; 16328017e9fSAlp Dener ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr); 164fed79b8eSAlp Dener 165fed79b8eSAlp Dener if (stepAccepted) { 166937a31a1SAlp Dener /* Step is good, evaluate the gradient and flip the need-Hessian switch */ 1678d5ead36SAlp Dener steplen = 1.0; 168937a31a1SAlp Dener needH = PETSC_TRUE; 169e465cd6fSAlp Dener ++bnk->newt; 170fed79b8eSAlp Dener ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 17161be54a6SAlp Dener ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr); 17261be54a6SAlp Dener ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 17361be54a6SAlp Dener ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr); 174f5766c09SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr); 175fed79b8eSAlp Dener } else { 176fed79b8eSAlp Dener /* Step is bad, revert old solution and re-solve with new radius*/ 1778d5ead36SAlp Dener steplen = 0.0; 178937a31a1SAlp Dener needH = PETSC_FALSE; 179fed79b8eSAlp Dener bnk->f = bnk->fold; 180fed79b8eSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 181fed79b8eSAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 182fed79b8eSAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 18373e4db90SAlp Dener if (oldTrust == tao->trust) { 18473e4db90SAlp Dener /* Can't change the radius anymore so just terminate */ 185fed79b8eSAlp Dener tao->reason = TAO_DIVERGED_TR_REDUCTION; 186fed79b8eSAlp Dener } 187fed79b8eSAlp Dener } 188fed79b8eSAlp Dener 189fed79b8eSAlp Dener /* Check for termination */ 1900b7db9bbSAlp Dener ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr); 1910b7db9bbSAlp Dener ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr); 192b4a30f08SAlp Dener if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 1939b6ef848SAlp Dener ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 1949b6ef848SAlp Dener ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr); 195fed79b8eSAlp Dener ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr); 196fed79b8eSAlp Dener } 197937a31a1SAlp Dener } 198fed79b8eSAlp Dener PetscFunctionReturn(0); 199fed79b8eSAlp Dener } 200fed79b8eSAlp Dener 201df278d8fSAlp Dener /*------------------------------------------------------------*/ 202df278d8fSAlp Dener 203e0ed867bSAlp Dener static PetscErrorCode TaoSetFromOptions_BNTR(PetscOptionItems *PetscOptionsObject,Tao tao) 2049b6ef848SAlp Dener { 2059b6ef848SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 2069b6ef848SAlp Dener PetscErrorCode ierr; 2079b6ef848SAlp Dener 2089b6ef848SAlp Dener PetscFunctionBegin; 209e0ed867bSAlp Dener ierr = TaoSetFromOptions_BNK(PetscOptionsObject, tao);CHKERRQ(ierr); 210e0ed867bSAlp Dener if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION; 2119b6ef848SAlp Dener if (!bnk->is_nash && !bnk->is_stcg && !bnk->is_gltr) SETERRQ(PETSC_COMM_SELF,1,"Must use a trust-region CG method for KSP (KSPNASH, KSPSTCG, KSPGLTR)"); 2129b6ef848SAlp Dener PetscFunctionReturn(0); 2139b6ef848SAlp Dener } 2149b6ef848SAlp Dener 2159b6ef848SAlp Dener /*------------------------------------------------------------*/ 2163850be85SAlp Dener /*MC 2173850be85SAlp Dener TAOBNTR - Bounded Newton Trust Region for nonlinear minimization with bound constraints. 2189b6ef848SAlp Dener 2193850be85SAlp Dener Options Database Keys: 2203850be85SAlp Dener + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop 2213850be85SAlp Dener . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation") 2223850be85SAlp Dener . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation") 2233850be85SAlp Dener - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas") 2243850be85SAlp Dener 2253850be85SAlp Dener Level: beginner 2263850be85SAlp Dener M*/ 227e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao) 228fed79b8eSAlp Dener { 229fed79b8eSAlp Dener TAO_BNK *bnk; 230fed79b8eSAlp Dener PetscErrorCode ierr; 231fed79b8eSAlp Dener 232fed79b8eSAlp Dener PetscFunctionBegin; 233fed79b8eSAlp Dener ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 234fed79b8eSAlp Dener tao->ops->solve=TaoSolve_BNTR; 235e0ed867bSAlp Dener tao->ops->setfromoptions=TaoSetFromOptions_BNTR; 236fed79b8eSAlp Dener 237fed79b8eSAlp Dener bnk = (TAO_BNK *)tao->data; 23866ed3702SAlp Dener bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */ 239fed79b8eSAlp Dener PetscFunctionReturn(0); 240fed79b8eSAlp Dener } 241