1eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h> 2eb910715SAlp Dener #include <petscksp.h> 3eb910715SAlp Dener 4eb910715SAlp Dener /* 5198282dbSAlp Dener Implements Newton's Method with a line search approach for 6198282dbSAlp Dener solving bound constrained minimization problems. 7eb910715SAlp Dener 8198282dbSAlp Dener ------------------------------------------------------------ 9eb910715SAlp 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_DIRECTION) 15198282dbSAlp Dener niter = 0 16c4b75bccSAlp Dener step_accepted = true 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 35c4b75bccSAlp 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 p > 0 54c4b75bccSAlp Dener Hr_k += p* 55198282dbSAlp Dener end 56198282dbSAlp Dener if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS 57198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6) 58198282dbSAlp Dener scale BFGS with VecReciprocal(D) 59198282dbSAlp Dener end 60198282dbSAlp Dener solve Hr_k dr_k = -gr_k 61198282dbSAlp Dener set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F 62198282dbSAlp Dener 63198282dbSAlp Dener if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf 64198282dbSAlp Dener dr_k = -BFGS*gr_k for variables in I 65198282dbSAlp Dener if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf 66198282dbSAlp Dener reset the BFGS preconditioner 67198282dbSAlp Dener calculate scale delta and apply it to BFGS 68198282dbSAlp Dener dr_k = -BFGS*gr_k for variables in I 69198282dbSAlp Dener if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf 70198282dbSAlp Dener dr_k = -gr_k for variables in I 71198282dbSAlp Dener end 72198282dbSAlp Dener end 73198282dbSAlp Dener end 74198282dbSAlp Dener 75198282dbSAlp Dener x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch() 76198282dbSAlp Dener if ls_failed 77198282dbSAlp Dener f_{k+1} = f_k 78198282dbSAlp Dener x_{k+1} = x_k 79198282dbSAlp Dener g_{k+1} = g_k 80198282dbSAlp Dener pg_{k+1} = pg_k 81198282dbSAlp Dener terminate 82198282dbSAlp Dener else 83c4b75bccSAlp Dener pg_{k+1} = project(g_{k+1}) 84198282dbSAlp Dener count the accepted step type (Newton, BFGS, scaled grad or grad) 85198282dbSAlp Dener end 86198282dbSAlp Dener 87198282dbSAlp Dener check convergence at pg_{k+1} 88198282dbSAlp Dener end 89eb910715SAlp Dener */ 90eb910715SAlp Dener 91e0ed867bSAlp Dener PetscErrorCode TaoSolve_BNLS(Tao tao) 92eb910715SAlp Dener { 93eb910715SAlp Dener PetscErrorCode ierr; 94eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 95e465cd6fSAlp Dener KSPConvergedReason ksp_reason; 96eb910715SAlp Dener TaoLineSearchConvergedReason ls_reason; 97eb910715SAlp Dener 9889da521bSAlp Dener PetscReal steplen = 1.0, resnorm; 99937a31a1SAlp Dener PetscBool cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_TRUE; 100eb910715SAlp Dener PetscInt stepType; 101eb910715SAlp Dener 102eb910715SAlp Dener PetscFunctionBegin; 10328017e9fSAlp Dener /* Initialize the preconditioner, KSP solver and trust radius/line search */ 104eb910715SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 105937a31a1SAlp Dener ierr = TaoBNKInitialize(tao, bnk->init_type, &needH);CHKERRQ(ierr); 10628017e9fSAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 107eb910715SAlp Dener 108eb910715SAlp Dener /* Have not converged; continue with Newton method */ 109eb910715SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 110e1e80dc8SAlp Dener /* Call general purpose update function */ 111e1e80dc8SAlp Dener if (tao->ops->update) { 112*8fcddce6SStefano Zampini ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr); 113e1e80dc8SAlp Dener } 114eb910715SAlp Dener ++tao->niter; 115eb910715SAlp Dener 11689da521bSAlp Dener if (needH && bnk->inactive_idx) { 117c0f10754SAlp Dener /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */ 118c0f10754SAlp Dener ierr = TaoBNKTakeCGSteps(tao, &cgTerminate);CHKERRQ(ierr); 119c0f10754SAlp Dener if (cgTerminate) { 120c0f10754SAlp Dener tao->reason = bnk->bncg->reason; 121c0f10754SAlp Dener PetscFunctionReturn(0); 122c0f10754SAlp Dener } 12308752603SAlp Dener /* Compute the hessian and update the BFGS preconditioner at the new iterate */ 124f7bf01afSAlp Dener ierr = (*bnk->computehessian)(tao);CHKERRQ(ierr); 125937a31a1SAlp Dener needH = PETSC_FALSE; 126937a31a1SAlp Dener } 127fed79b8eSAlp Dener 1288d5ead36SAlp Dener /* Use the common BNK kernel to compute the safeguarded Newton step (for inactive variables only) */ 129f7bf01afSAlp Dener ierr = (*bnk->computestep)(tao, shift, &ksp_reason, &stepType);CHKERRQ(ierr); 130e465cd6fSAlp Dener ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr); 131eb910715SAlp Dener 132080d2917SAlp Dener /* Store current solution before it changes */ 133080d2917SAlp Dener bnk->fold = bnk->f; 134eb910715SAlp Dener ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 135eb910715SAlp Dener ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr); 13609164190SAlp Dener ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr); 137eb910715SAlp Dener 138c14b763aSAlp Dener /* Trigger the line search */ 139937a31a1SAlp Dener ierr = TaoBNKPerformLineSearch(tao, &stepType, &steplen, &ls_reason);CHKERRQ(ierr); 140eb910715SAlp Dener 141eb910715SAlp Dener if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 142eb910715SAlp Dener /* Failed to find an improving point */ 143937a31a1SAlp Dener needH = PETSC_FALSE; 144080d2917SAlp Dener bnk->f = bnk->fold; 145eb910715SAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 146eb910715SAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 14709164190SAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 148c14b763aSAlp Dener steplen = 0.0; 149eb910715SAlp Dener tao->reason = TAO_DIVERGED_LS_FAILURE; 150e465cd6fSAlp Dener } else { 151937a31a1SAlp Dener /* new iterate so we need to recompute the Hessian */ 152937a31a1SAlp Dener needH = PETSC_TRUE; 153198282dbSAlp Dener /* compute the projected gradient */ 15433c78596SAlp Dener ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr); 15561be54a6SAlp Dener ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 15661be54a6SAlp Dener ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr); 157f5766c09SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr); 1589b6ef848SAlp Dener /* update the trust radius based on the step length */ 1599b6ef848SAlp Dener ierr = TaoBNKUpdateTrustRadius(tao, 0.0, 0.0, BNK_UPDATE_STEP, stepType, &stepAccepted);CHKERRQ(ierr); 16062675beeSAlp Dener /* count the accepted step type */ 16162675beeSAlp Dener ierr = TaoBNKAddStepCounts(tao, stepType);CHKERRQ(ierr); 162937a31a1SAlp Dener /* active BNCG recycling for next iteration */ 163937a31a1SAlp Dener ierr = TaoBNCGSetRecycleFlag(bnk->bncg, PETSC_TRUE);CHKERRQ(ierr); 164eb910715SAlp Dener } 165eb910715SAlp Dener 166eb910715SAlp Dener /* Check for termination */ 1670b7db9bbSAlp Dener ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr); 1680b7db9bbSAlp Dener ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr); 169b4a30f08SAlp Dener if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 1709b6ef848SAlp Dener ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 1719b6ef848SAlp Dener ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr); 172eb910715SAlp Dener ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr); 173eb910715SAlp Dener } 174eb910715SAlp Dener PetscFunctionReturn(0); 175eb910715SAlp Dener } 176eb910715SAlp Dener 177df278d8fSAlp Dener /*------------------------------------------------------------*/ 1783850be85SAlp Dener /*MC 1793850be85SAlp Dener TAOBNLS - Bounded Newton Line Search for nonlinear minimization with bound constraints. 180df278d8fSAlp Dener 1813850be85SAlp Dener Options Database Keys: 1823850be85SAlp Dener + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop 1833850be85SAlp Dener . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation") 1843850be85SAlp Dener . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation") 1853850be85SAlp Dener - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas") 1863850be85SAlp Dener 1873850be85SAlp Dener Level: beginner 1883850be85SAlp Dener M*/ 189e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao) 190eb910715SAlp Dener { 191fed79b8eSAlp Dener TAO_BNK *bnk; 192eb910715SAlp Dener PetscErrorCode ierr; 193eb910715SAlp Dener 194eb910715SAlp Dener PetscFunctionBegin; 195eb910715SAlp Dener ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 196eb910715SAlp Dener tao->ops->solve = TaoSolve_BNLS; 197fed79b8eSAlp Dener 198fed79b8eSAlp Dener bnk = (TAO_BNK *)tao->data; 199e031d6f5SAlp Dener bnk->init_type = BNK_INIT_DIRECTION; 20066ed3702SAlp Dener bnk->update_type = BNK_UPDATE_STEP; /* trust region updates based on line search step length */ 201eb910715SAlp Dener PetscFunctionReturn(0); 202eb910715SAlp Dener } 203