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 19c4b75bccSAlp Dener 20c4b75bccSAlp Dener if needH 21c4b75bccSAlp Dener If max_cg_steps > 0 22c4b75bccSAlp Dener x_k, g_k, pg_k = TaoSolve(BNCG) 23c4b75bccSAlp Dener end 24c4b75bccSAlp Dener 25198282dbSAlp Dener H_k = TaoComputeHessian(x_k) 26198282dbSAlp Dener if pc_type == BNK_PC_BFGS 27198282dbSAlp Dener add correction to BFGS approx 28198282dbSAlp Dener if scale_type == BNK_SCALE_AHESS 29198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(H_k)), 1e6) 30198282dbSAlp Dener scale BFGS with VecReciprocal(D) 31198282dbSAlp Dener end 32198282dbSAlp Dener end 33c4b75bccSAlp Dener needH = False 34198282dbSAlp Dener end 35198282dbSAlp Dener 36198282dbSAlp Dener if pc_type = BNK_PC_BFGS 37198282dbSAlp Dener B_k = BFGS 38198282dbSAlp Dener else 39198282dbSAlp Dener B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6) 40198282dbSAlp Dener B_k = VecReciprocal(B_k) 41198282dbSAlp Dener end 42198282dbSAlp Dener w = x_k - VecMedian(x_k - 0.001*B_k*g_k) 43198282dbSAlp Dener eps = min(eps, norm2(w)) 44198282dbSAlp Dener determine the active and inactive index sets such that 45198282dbSAlp Dener L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0} 46198282dbSAlp Dener U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0} 47198282dbSAlp Dener F = {i : l_i = (x_k)_i = u_i} 48198282dbSAlp Dener A = {L + U + F} 49c4b75bccSAlp Dener IA = {i : i not in A} 50198282dbSAlp Dener 51c4b75bccSAlp Dener generate the reduced system Hr_k dr_k = -gr_k for variables in IA 52198282dbSAlp Dener if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS 53198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6) 54198282dbSAlp Dener scale BFGS with VecReciprocal(D) 55198282dbSAlp Dener end 56c4b75bccSAlp Dener 57c4b75bccSAlp Dener while !stepAccepted 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 x_{k+1} = VecMedian(x_k + d_k) 62198282dbSAlp Dener s = x_{k+1} - x_k 63198282dbSAlp Dener prered = dot(s, 0.5*gr_k - Hr_k*s) 64198282dbSAlp Dener f_{k+1} = TaoComputeObjective(x_{k+1}) 65198282dbSAlp Dener actred = f_k - f_{k+1} 66198282dbSAlp Dener 67198282dbSAlp Dener oldTrust = trust 68198282dbSAlp Dener step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION) 69198282dbSAlp Dener if step_accepted 70198282dbSAlp Dener g_{k+1} = TaoComputeGradient(x_{k+1}) 71c4b75bccSAlp Dener pg_{k+1} = project(g_{k+1}) 72198282dbSAlp Dener count the accepted Newton step 73c4b75bccSAlp Dener needH = True 74198282dbSAlp Dener else 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 if trust == oldTrust 80198282dbSAlp Dener terminate because we cannot shrink the radius any further 81198282dbSAlp Dener end 82198282dbSAlp Dener end 83198282dbSAlp Dener 84198282dbSAlp Dener end 85e84e3fd2SStefano Zampini check convergence at pg_{k+1} 860279bc1bSStefano Zampini niter += 1 87c4b75bccSAlp Dener 88c4b75bccSAlp Dener end 89fed79b8eSAlp Dener */ 90fed79b8eSAlp Dener 91d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSolve_BNTR(Tao tao) 92d71ae5a4SJacob Faibussowitsch { 93fed79b8eSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 94e465cd6fSAlp Dener KSPConvergedReason ksp_reason; 95fed79b8eSAlp Dener 96e84e3fd2SStefano Zampini PetscReal oldTrust, prered, actred, steplen = 0.0, resnorm; 97937a31a1SAlp Dener PetscBool cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE; 986b591159SAlp Dener PetscInt stepType, nDiff; 99fed79b8eSAlp Dener 100fed79b8eSAlp Dener PetscFunctionBegin; 10128017e9fSAlp Dener /* Initialize the preconditioner, KSP solver and trust radius/line search */ 102fed79b8eSAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 1039566063dSJacob Faibussowitsch PetscCall(TaoBNKInitialize(tao, bnk->init_type, &needH)); 1043ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 105fed79b8eSAlp Dener 106fed79b8eSAlp Dener /* Have not converged; continue with Newton method */ 107fed79b8eSAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 108e1e80dc8SAlp Dener /* Call general purpose update function */ 109e1e80dc8SAlp Dener if (tao->ops->update) { 110dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, update, tao->niter, tao->user_update); 1117494f0b1SStefano Zampini PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient)); 112e1e80dc8SAlp Dener } 113e031d6f5SAlp Dener 11489da521bSAlp Dener if (needH && bnk->inactive_idx) { 115e031d6f5SAlp Dener /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */ 1169566063dSJacob Faibussowitsch PetscCall(TaoBNKTakeCGSteps(tao, &cgTerminate)); 117e031d6f5SAlp Dener if (cgTerminate) { 118e031d6f5SAlp Dener tao->reason = bnk->bncg->reason; 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 120fed79b8eSAlp Dener } 121937a31a1SAlp Dener /* Compute the hessian and update the BFGS preconditioner at the new iterate */ 1229566063dSJacob Faibussowitsch PetscCall((*bnk->computehessian)(tao)); 123937a31a1SAlp Dener needH = PETSC_FALSE; 124937a31a1SAlp Dener } 125fed79b8eSAlp Dener 126fed79b8eSAlp Dener /* Store current solution before it changes */ 127fed79b8eSAlp Dener bnk->fold = bnk->f; 1289566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, bnk->Xold)); 1299566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, bnk->Gold)); 1309566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old)); 131fed79b8eSAlp Dener 132937a31a1SAlp Dener /* Enter into trust region loops */ 133937a31a1SAlp Dener stepAccepted = PETSC_FALSE; 134937a31a1SAlp Dener while (!stepAccepted && tao->reason == TAO_CONTINUE_ITERATING) { 135937a31a1SAlp Dener tao->ksp_its = 0; 136937a31a1SAlp Dener 137937a31a1SAlp Dener /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */ 1389566063dSJacob Faibussowitsch PetscCall((*bnk->computestep)(tao, shift, &ksp_reason, &stepType)); 139937a31a1SAlp Dener 140b1c2d0e3SAlp Dener /* Temporarily accept the step and project it into the bounds */ 1419566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution, 1.0, tao->stepdirection)); 1429566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution)); 143b1c2d0e3SAlp Dener 144b1c2d0e3SAlp Dener /* Check if the projection changed the step direction */ 145c4b75bccSAlp Dener if (nDiff > 0) { 146c4b75bccSAlp Dener /* Projection changed the step, so we have to recompute the step and 147c4b75bccSAlp Dener the predicted reduction. Leave the trust radius unchanged. */ 1489566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, tao->stepdirection)); 1499566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold)); 1509566063dSJacob Faibussowitsch PetscCall(TaoBNKRecomputePred(tao, tao->stepdirection, &prered)); 151b1c2d0e3SAlp Dener } else { 152b1c2d0e3SAlp Dener /* Step did not change, so we can just recover the pre-computed prediction */ 1539566063dSJacob Faibussowitsch PetscCall(KSPCGGetObjFcn(tao->ksp, &prered)); 154b1c2d0e3SAlp Dener } 155b1c2d0e3SAlp Dener prered = -prered; 156b1c2d0e3SAlp Dener 157b1c2d0e3SAlp Dener /* Compute the actual reduction and update the trust radius */ 1589566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, tao->solution, &bnk->f)); 1593c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(bnk->f), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 160b1c2d0e3SAlp Dener actred = bnk->fold - bnk->f; 161e761ccfdSAlp Dener oldTrust = tao->trust; 1629566063dSJacob Faibussowitsch PetscCall(TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted)); 163fed79b8eSAlp Dener 164fed79b8eSAlp Dener if (stepAccepted) { 165937a31a1SAlp Dener /* Step is good, evaluate the gradient and flip the need-Hessian switch */ 1668d5ead36SAlp Dener steplen = 1.0; 167937a31a1SAlp Dener needH = PETSC_TRUE; 168e465cd6fSAlp Dener ++bnk->newt; 1699566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient)); 1709566063dSJacob Faibussowitsch PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type)); 1719566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient)); 172*976ed0a4SStefano Zampini if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0)); 1739566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm)); 174fed79b8eSAlp Dener } else { 175fed79b8eSAlp Dener /* Step is bad, revert old solution and re-solve with new radius*/ 1768d5ead36SAlp Dener steplen = 0.0; 177937a31a1SAlp Dener needH = PETSC_FALSE; 178fed79b8eSAlp Dener bnk->f = bnk->fold; 1799566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->Xold, tao->solution)); 1809566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->Gold, tao->gradient)); 1819566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient)); 18273e4db90SAlp Dener if (oldTrust == tao->trust) { 18373e4db90SAlp Dener /* Can't change the radius anymore so just terminate */ 184fed79b8eSAlp Dener tao->reason = TAO_DIVERGED_TR_REDUCTION; 185fed79b8eSAlp Dener } 186fed79b8eSAlp Dener } 187e84e3fd2SStefano Zampini } 188fed79b8eSAlp Dener /* Check for termination */ 1899566063dSJacob Faibussowitsch PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W)); 1909566063dSJacob Faibussowitsch PetscCall(VecNorm(bnk->W, NORM_2, &resnorm)); 1913c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 192e84e3fd2SStefano Zampini ++tao->niter; 1939566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its)); 1949566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen)); 195dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 196937a31a1SAlp Dener } 1973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 198fed79b8eSAlp Dener } 199fed79b8eSAlp Dener 200df278d8fSAlp Dener /*------------------------------------------------------------*/ 201d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_BNTR(Tao tao) 202d71ae5a4SJacob Faibussowitsch { 2032e6e4ca1SStefano Zampini KSP ksp; 2042e6e4ca1SStefano Zampini PetscVoidFunction valid; 2055eb5f4d6SAlp Dener 2065eb5f4d6SAlp Dener PetscFunctionBegin; 2079566063dSJacob Faibussowitsch PetscCall(TaoSetUp_BNK(tao)); 2089566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 2099566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)ksp, "KSPCGSetRadius_C", &valid)); 2103c859ba3SBarry Smith PetscCheck(valid, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Not for KSP type %s. Must use a trust-region CG method for KSP (e.g. KSPNASH, KSPSTCG, KSPGLTR)", ((PetscObject)ksp)->type_name); 2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2125eb5f4d6SAlp Dener } 2135eb5f4d6SAlp Dener 2145eb5f4d6SAlp Dener /*------------------------------------------------------------*/ 215df278d8fSAlp Dener 216d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetFromOptions_BNTR(Tao tao, PetscOptionItems *PetscOptionsObject) 217d71ae5a4SJacob Faibussowitsch { 2189b6ef848SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 2199b6ef848SAlp Dener 2209b6ef848SAlp Dener PetscFunctionBegin; 221dbbe0bcdSBarry Smith PetscCall(TaoSetFromOptions_BNK(tao, PetscOptionsObject)); 222e0ed867bSAlp Dener if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION; 2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2249b6ef848SAlp Dener } 2259b6ef848SAlp Dener 2269b6ef848SAlp Dener /*------------------------------------------------------------*/ 2273850be85SAlp Dener /*MC 2283850be85SAlp Dener TAOBNTR - Bounded Newton Trust Region for nonlinear minimization with bound constraints. 2299b6ef848SAlp Dener 2303850be85SAlp Dener Options Database Keys: 2313850be85SAlp Dener + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop 2323850be85SAlp Dener . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation") 2333850be85SAlp Dener . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation") 2343850be85SAlp Dener - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas") 2353850be85SAlp Dener 2363850be85SAlp Dener Level: beginner 2373850be85SAlp Dener M*/ 238d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao) 239d71ae5a4SJacob Faibussowitsch { 240fed79b8eSAlp Dener TAO_BNK *bnk; 241fed79b8eSAlp Dener 242fed79b8eSAlp Dener PetscFunctionBegin; 2439566063dSJacob Faibussowitsch PetscCall(TaoCreate_BNK(tao)); 244fed79b8eSAlp Dener tao->ops->solve = TaoSolve_BNTR; 2455eb5f4d6SAlp Dener tao->ops->setup = TaoSetUp_BNTR; 246e0ed867bSAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNTR; 247fed79b8eSAlp Dener 248fed79b8eSAlp Dener bnk = (TAO_BNK *)tao->data; 24966ed3702SAlp Dener bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */ 2503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 251fed79b8eSAlp Dener } 252