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 19c4b75bccSAlp Dener if needH 20c4b75bccSAlp Dener If max_cg_steps > 0 21c4b75bccSAlp Dener x_k, g_k, pg_k = TaoSolve(BNCG) 22c4b75bccSAlp Dener end 23c4b75bccSAlp Dener 24198282dbSAlp Dener H_k = TaoComputeHessian(x_k) 25198282dbSAlp Dener if pc_type == BNK_PC_BFGS 26198282dbSAlp Dener add correction to BFGS approx 27198282dbSAlp Dener if scale_type == BNK_SCALE_AHESS 28198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(H_k)), 1e6) 29198282dbSAlp Dener scale BFGS with VecReciprocal(D) 30198282dbSAlp Dener end 31198282dbSAlp Dener end 32c4b75bccSAlp Dener needH = False 33c4b75bccSAlp Dener end 34198282dbSAlp Dener 35198282dbSAlp Dener if pc_type = BNK_PC_BFGS 36198282dbSAlp Dener B_k = BFGS 37198282dbSAlp Dener else 38198282dbSAlp Dener B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6) 39198282dbSAlp Dener B_k = VecReciprocal(B_k) 40198282dbSAlp Dener end 41198282dbSAlp Dener w = x_k - VecMedian(x_k - 0.001*B_k*g_k) 42198282dbSAlp Dener eps = min(eps, norm2(w)) 43198282dbSAlp Dener determine the active and inactive index sets such that 44198282dbSAlp Dener L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0} 45198282dbSAlp Dener U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0} 46198282dbSAlp Dener F = {i : l_i = (x_k)_i = u_i} 47198282dbSAlp Dener A = {L + U + F} 48c4b75bccSAlp Dener IA = {i : i not in A} 49198282dbSAlp Dener 50c4b75bccSAlp Dener generate the reduced system Hr_k dr_k = -gr_k for variables in IA 51198282dbSAlp Dener if p > 0 52c4b75bccSAlp Dener Hr_k += p* 53198282dbSAlp Dener end 54198282dbSAlp Dener if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS 55198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6) 56198282dbSAlp Dener scale BFGS with VecReciprocal(D) 57198282dbSAlp Dener end 58198282dbSAlp Dener solve Hr_k dr_k = -gr_k 59198282dbSAlp Dener set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F 60198282dbSAlp Dener 61198282dbSAlp Dener if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf 62198282dbSAlp Dener dr_k = -BFGS*gr_k for variables in I 63198282dbSAlp Dener if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf 64198282dbSAlp Dener reset the BFGS preconditioner 65198282dbSAlp Dener calculate scale delta and apply it to BFGS 66198282dbSAlp Dener dr_k = -BFGS*gr_k for variables in I 67198282dbSAlp Dener if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf 68198282dbSAlp Dener dr_k = -gr_k for variables in I 69198282dbSAlp Dener end 70198282dbSAlp Dener end 71198282dbSAlp Dener end 72198282dbSAlp Dener 73198282dbSAlp Dener x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch() 74198282dbSAlp Dener if ls_failed 75198282dbSAlp Dener f_{k+1} = f_k 76198282dbSAlp Dener x_{k+1} = x_k 77198282dbSAlp Dener g_{k+1} = g_k 78198282dbSAlp Dener pg_{k+1} = pg_k 79198282dbSAlp Dener terminate 80198282dbSAlp Dener else 81c4b75bccSAlp Dener pg_{k+1} = project(g_{k+1}) 82198282dbSAlp Dener count the accepted step type (Newton, BFGS, scaled grad or grad) 83198282dbSAlp Dener end 84198282dbSAlp Dener 850279bc1bSStefano Zampini niter += 1 86198282dbSAlp Dener check convergence at pg_{k+1} 87198282dbSAlp Dener end 88eb910715SAlp Dener */ 89eb910715SAlp Dener 90d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSolve_BNLS(Tao tao) 91d71ae5a4SJacob Faibussowitsch { 92eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 93e465cd6fSAlp Dener KSPConvergedReason ksp_reason; 94eb910715SAlp Dener TaoLineSearchConvergedReason ls_reason; 9589da521bSAlp Dener PetscReal steplen = 1.0, resnorm; 96937a31a1SAlp Dener PetscBool cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_TRUE; 97eb910715SAlp Dener PetscInt stepType; 98eb910715SAlp Dener 99eb910715SAlp Dener PetscFunctionBegin; 10028017e9fSAlp Dener /* Initialize the preconditioner, KSP solver and trust radius/line search */ 101eb910715SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 1029566063dSJacob Faibussowitsch PetscCall(TaoBNKInitialize(tao, bnk->init_type, &needH)); 1033ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 104eb910715SAlp Dener 105eb910715SAlp Dener /* Have not converged; continue with Newton method */ 106eb910715SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 107e1e80dc8SAlp Dener /* Call general purpose update function */ 108e1e80dc8SAlp Dener if (tao->ops->update) { 109dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, update, tao->niter, tao->user_update); 1107494f0b1SStefano Zampini PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient)); 111e1e80dc8SAlp Dener } 112eb910715SAlp Dener 11389da521bSAlp Dener if (needH && bnk->inactive_idx) { 114c0f10754SAlp Dener /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */ 1159566063dSJacob Faibussowitsch PetscCall(TaoBNKTakeCGSteps(tao, &cgTerminate)); 116c0f10754SAlp Dener if (cgTerminate) { 117c0f10754SAlp Dener tao->reason = bnk->bncg->reason; 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119c0f10754SAlp Dener } 12008752603SAlp Dener /* Compute the hessian and update the BFGS preconditioner at the new iterate */ 1219566063dSJacob Faibussowitsch PetscCall((*bnk->computehessian)(tao)); 122937a31a1SAlp Dener needH = PETSC_FALSE; 123937a31a1SAlp Dener } 124fed79b8eSAlp Dener 1258d5ead36SAlp Dener /* Use the common BNK kernel to compute the safeguarded Newton step (for inactive variables only) */ 1269566063dSJacob Faibussowitsch PetscCall((*bnk->computestep)(tao, shift, &ksp_reason, &stepType)); 1279566063dSJacob Faibussowitsch PetscCall(TaoBNKSafeguardStep(tao, ksp_reason, &stepType)); 128eb910715SAlp Dener 129080d2917SAlp Dener /* Store current solution before it changes */ 130080d2917SAlp Dener bnk->fold = bnk->f; 1319566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, bnk->Xold)); 1329566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, bnk->Gold)); 1339566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old)); 134eb910715SAlp Dener 135c14b763aSAlp Dener /* Trigger the line search */ 1369566063dSJacob Faibussowitsch PetscCall(TaoBNKPerformLineSearch(tao, &stepType, &steplen, &ls_reason)); 137eb910715SAlp Dener 138eb910715SAlp Dener if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 139eb910715SAlp Dener /* Failed to find an improving point */ 140937a31a1SAlp Dener needH = PETSC_FALSE; 141080d2917SAlp Dener bnk->f = bnk->fold; 1429566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->Xold, tao->solution)); 1439566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->Gold, tao->gradient)); 1449566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient)); 145c14b763aSAlp Dener steplen = 0.0; 146eb910715SAlp Dener tao->reason = TAO_DIVERGED_LS_FAILURE; 147e465cd6fSAlp Dener } else { 148937a31a1SAlp Dener /* new iterate so we need to recompute the Hessian */ 149937a31a1SAlp Dener needH = PETSC_TRUE; 150198282dbSAlp Dener /* compute the projected gradient */ 1519566063dSJacob Faibussowitsch PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type)); 1529566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient)); 153*976ed0a4SStefano Zampini if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0)); 1549566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm)); 1559b6ef848SAlp Dener /* update the trust radius based on the step length */ 1569566063dSJacob Faibussowitsch PetscCall(TaoBNKUpdateTrustRadius(tao, 0.0, 0.0, BNK_UPDATE_STEP, stepType, &stepAccepted)); 15762675beeSAlp Dener /* count the accepted step type */ 1589566063dSJacob Faibussowitsch PetscCall(TaoBNKAddStepCounts(tao, stepType)); 159937a31a1SAlp Dener /* active BNCG recycling for next iteration */ 1609566063dSJacob Faibussowitsch PetscCall(TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE)); 161eb910715SAlp Dener } 162eb910715SAlp Dener 163eb910715SAlp Dener /* Check for termination */ 1649566063dSJacob Faibussowitsch PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W)); 1659566063dSJacob Faibussowitsch PetscCall(VecNorm(bnk->W, NORM_2, &resnorm)); 1663c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 1670279bc1bSStefano Zampini ++tao->niter; 1689566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its)); 1699566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen)); 170dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 171eb910715SAlp Dener } 1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173eb910715SAlp Dener } 174eb910715SAlp Dener 175df278d8fSAlp Dener /*------------------------------------------------------------*/ 1763850be85SAlp Dener /*MC 1773850be85SAlp Dener TAOBNLS - Bounded Newton Line Search for nonlinear minimization with bound constraints. 178df278d8fSAlp Dener 1793850be85SAlp Dener Options Database Keys: 1803850be85SAlp Dener + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop 1813850be85SAlp Dener . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation") 1823850be85SAlp Dener . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation") 1833850be85SAlp Dener - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas") 1843850be85SAlp Dener 1853850be85SAlp Dener Level: beginner 1863850be85SAlp Dener M*/ 187d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao) 188d71ae5a4SJacob Faibussowitsch { 189fed79b8eSAlp Dener TAO_BNK *bnk; 190eb910715SAlp Dener 191eb910715SAlp Dener PetscFunctionBegin; 1929566063dSJacob Faibussowitsch PetscCall(TaoCreate_BNK(tao)); 193eb910715SAlp Dener tao->ops->solve = TaoSolve_BNLS; 194fed79b8eSAlp Dener 195fed79b8eSAlp Dener bnk = (TAO_BNK *)tao->data; 196e031d6f5SAlp Dener bnk->init_type = BNK_INIT_DIRECTION; 19766ed3702SAlp Dener bnk->update_type = BNK_UPDATE_STEP; /* trust region updates based on line search step length */ 1983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 199eb910715SAlp Dener } 200