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 x_0 = VecMedian(x_0)
9198282dbSAlp Dener f_0, g_0= TaoComputeObjectiveAndGradient(x_0)
10c4b75bccSAlp Dener pg_0 = project(g_0)
11198282dbSAlp Dener check convergence at pg_0
12c4b75bccSAlp Dener needH = TaoBNKInitialize(default:BNK_INIT_INTERPOLATION)
13198282dbSAlp Dener niter = 0
14c4b75bccSAlp Dener step_accepted = false
15198282dbSAlp Dener
16198282dbSAlp Dener while niter <= max_it
17c4b75bccSAlp Dener
18c4b75bccSAlp Dener if needH
19c4b75bccSAlp Dener If max_cg_steps > 0
20c4b75bccSAlp Dener x_k, g_k, pg_k = TaoSolve(BNCG)
21c4b75bccSAlp Dener end
22c4b75bccSAlp Dener
23198282dbSAlp Dener H_k = TaoComputeHessian(x_k)
24198282dbSAlp Dener if pc_type == BNK_PC_BFGS
25198282dbSAlp Dener add correction to BFGS approx
26198282dbSAlp Dener if scale_type == BNK_SCALE_AHESS
27198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
28198282dbSAlp Dener scale BFGS with VecReciprocal(D)
29198282dbSAlp Dener end
30198282dbSAlp Dener end
31c4b75bccSAlp Dener needH = False
32198282dbSAlp Dener end
33198282dbSAlp Dener
34198282dbSAlp Dener if pc_type = BNK_PC_BFGS
35198282dbSAlp Dener B_k = BFGS
36198282dbSAlp Dener else
37198282dbSAlp Dener B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
38198282dbSAlp Dener B_k = VecReciprocal(B_k)
39198282dbSAlp Dener end
40198282dbSAlp Dener w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
41198282dbSAlp Dener eps = min(eps, norm2(w))
42198282dbSAlp Dener determine the active and inactive index sets such that
43198282dbSAlp Dener L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
44198282dbSAlp Dener U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
45198282dbSAlp Dener F = {i : l_i = (x_k)_i = u_i}
46198282dbSAlp Dener A = {L + U + F}
47c4b75bccSAlp Dener IA = {i : i not in A}
48198282dbSAlp Dener
49c4b75bccSAlp Dener generate the reduced system Hr_k dr_k = -gr_k for variables in IA
50198282dbSAlp Dener if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
51198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
52198282dbSAlp Dener scale BFGS with VecReciprocal(D)
53198282dbSAlp Dener end
54c4b75bccSAlp Dener
55c4b75bccSAlp Dener while !stepAccepted
56198282dbSAlp Dener solve Hr_k dr_k = -gr_k
57198282dbSAlp Dener set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
58198282dbSAlp Dener
59198282dbSAlp Dener x_{k+1} = VecMedian(x_k + d_k)
60198282dbSAlp Dener s = x_{k+1} - x_k
61198282dbSAlp Dener prered = dot(s, 0.5*gr_k - Hr_k*s)
62198282dbSAlp Dener f_{k+1} = TaoComputeObjective(x_{k+1})
63198282dbSAlp Dener actred = f_k - f_{k+1}
64198282dbSAlp Dener
65198282dbSAlp Dener oldTrust = trust
66198282dbSAlp Dener step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION)
67198282dbSAlp Dener if step_accepted
68198282dbSAlp Dener g_{k+1} = TaoComputeGradient(x_{k+1})
69c4b75bccSAlp Dener pg_{k+1} = project(g_{k+1})
70198282dbSAlp Dener count the accepted Newton step
71c4b75bccSAlp Dener needH = True
72198282dbSAlp Dener else
73198282dbSAlp Dener f_{k+1} = f_k
74198282dbSAlp Dener x_{k+1} = x_k
75198282dbSAlp Dener g_{k+1} = g_k
76198282dbSAlp Dener pg_{k+1} = pg_k
77198282dbSAlp Dener if trust == oldTrust
78198282dbSAlp Dener terminate because we cannot shrink the radius any further
79198282dbSAlp Dener end
80198282dbSAlp Dener end
81198282dbSAlp Dener
82198282dbSAlp Dener end
83e84e3fd2SStefano Zampini check convergence at pg_{k+1}
840279bc1bSStefano Zampini niter += 1
85c4b75bccSAlp Dener
86c4b75bccSAlp Dener end
87fed79b8eSAlp Dener */
88fed79b8eSAlp Dener
TaoSolve_BNTR(Tao tao)89d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSolve_BNTR(Tao tao)
90d71ae5a4SJacob Faibussowitsch {
91fed79b8eSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
92e465cd6fSAlp Dener KSPConvergedReason ksp_reason;
93fed79b8eSAlp Dener
94e84e3fd2SStefano Zampini PetscReal oldTrust, prered, actred, steplen = 0.0, resnorm;
95937a31a1SAlp Dener PetscBool cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE;
966b591159SAlp Dener PetscInt stepType, nDiff;
97fed79b8eSAlp Dener
98fed79b8eSAlp Dener PetscFunctionBegin;
9928017e9fSAlp Dener /* Initialize the preconditioner, KSP solver and trust radius/line search */
100fed79b8eSAlp Dener tao->reason = TAO_CONTINUE_ITERATING;
1019566063dSJacob Faibussowitsch PetscCall(TaoBNKInitialize(tao, bnk->init_type, &needH));
1023ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
103fed79b8eSAlp Dener
104fed79b8eSAlp Dener /* Have not converged; continue with Newton method */
105fed79b8eSAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) {
106e1e80dc8SAlp Dener /* Call general purpose update function */
107e1e80dc8SAlp Dener if (tao->ops->update) {
108dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
109270bebe6SStefano Zampini PetscCall(TaoComputeObjective(tao, tao->solution, &bnk->f));
110e1e80dc8SAlp Dener }
111e031d6f5SAlp Dener
11289da521bSAlp Dener if (needH && bnk->inactive_idx) {
113e031d6f5SAlp Dener /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
1149566063dSJacob Faibussowitsch PetscCall(TaoBNKTakeCGSteps(tao, &cgTerminate));
115e031d6f5SAlp Dener if (cgTerminate) {
116e031d6f5SAlp Dener tao->reason = bnk->bncg->reason;
1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
118fed79b8eSAlp Dener }
119937a31a1SAlp Dener /* Compute the hessian and update the BFGS preconditioner at the new iterate */
1209566063dSJacob Faibussowitsch PetscCall((*bnk->computehessian)(tao));
121937a31a1SAlp Dener needH = PETSC_FALSE;
122937a31a1SAlp Dener }
123fed79b8eSAlp Dener
124fed79b8eSAlp Dener /* Store current solution before it changes */
125fed79b8eSAlp Dener bnk->fold = bnk->f;
1269566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, bnk->Xold));
1279566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, bnk->Gold));
1289566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old));
129fed79b8eSAlp Dener
130937a31a1SAlp Dener /* Enter into trust region loops */
131937a31a1SAlp Dener stepAccepted = PETSC_FALSE;
132937a31a1SAlp Dener while (!stepAccepted && tao->reason == TAO_CONTINUE_ITERATING) {
133937a31a1SAlp Dener tao->ksp_its = 0;
134937a31a1SAlp Dener
135937a31a1SAlp Dener /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
1369566063dSJacob Faibussowitsch PetscCall((*bnk->computestep)(tao, shift, &ksp_reason, &stepType));
137937a31a1SAlp Dener
138b1c2d0e3SAlp Dener /* Temporarily accept the step and project it into the bounds */
1399566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution, 1.0, tao->stepdirection));
1409566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
141b1c2d0e3SAlp Dener
142b1c2d0e3SAlp Dener /* Check if the projection changed the step direction */
143c4b75bccSAlp Dener if (nDiff > 0) {
144c4b75bccSAlp Dener /* Projection changed the step, so we have to recompute the step and
145c4b75bccSAlp Dener the predicted reduction. Leave the trust radius unchanged. */
1469566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, tao->stepdirection));
1479566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold));
1489566063dSJacob Faibussowitsch PetscCall(TaoBNKRecomputePred(tao, tao->stepdirection, &prered));
149b1c2d0e3SAlp Dener } else {
150b1c2d0e3SAlp Dener /* Step did not change, so we can just recover the pre-computed prediction */
1519566063dSJacob Faibussowitsch PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
152b1c2d0e3SAlp Dener }
153b1c2d0e3SAlp Dener prered = -prered;
154b1c2d0e3SAlp Dener
155b1c2d0e3SAlp Dener /* Compute the actual reduction and update the trust radius */
1569566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, tao->solution, &bnk->f));
157*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(bnk->f), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
158b1c2d0e3SAlp Dener actred = bnk->fold - bnk->f;
159e761ccfdSAlp Dener oldTrust = tao->trust;
1609566063dSJacob Faibussowitsch PetscCall(TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted));
161fed79b8eSAlp Dener
162fed79b8eSAlp Dener if (stepAccepted) {
163937a31a1SAlp Dener /* Step is good, evaluate the gradient and flip the need-Hessian switch */
1648d5ead36SAlp Dener steplen = 1.0;
165937a31a1SAlp Dener needH = PETSC_TRUE;
166e465cd6fSAlp Dener ++bnk->newt;
1679566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient));
1689566063dSJacob Faibussowitsch PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
1699566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
170976ed0a4SStefano Zampini if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
1719566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
172fed79b8eSAlp Dener } else {
173fed79b8eSAlp Dener /* Step is bad, revert old solution and re-solve with new radius*/
1748d5ead36SAlp Dener steplen = 0.0;
175937a31a1SAlp Dener needH = PETSC_FALSE;
176fed79b8eSAlp Dener bnk->f = bnk->fold;
1779566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->Xold, tao->solution));
1789566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->Gold, tao->gradient));
1799566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));
18073e4db90SAlp Dener if (oldTrust == tao->trust) {
18173e4db90SAlp Dener /* Can't change the radius anymore so just terminate */
182fed79b8eSAlp Dener tao->reason = TAO_DIVERGED_TR_REDUCTION;
183fed79b8eSAlp Dener }
184fed79b8eSAlp Dener }
185e84e3fd2SStefano Zampini }
186fed79b8eSAlp Dener /* Check for termination */
1879566063dSJacob Faibussowitsch PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
1889566063dSJacob Faibussowitsch PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
189*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
190e84e3fd2SStefano Zampini ++tao->niter;
1919566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
1929566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen));
193dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
194937a31a1SAlp Dener }
1953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
196fed79b8eSAlp Dener }
197fed79b8eSAlp Dener
TaoSetUp_BNTR(Tao tao)198d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_BNTR(Tao tao)
199d71ae5a4SJacob Faibussowitsch {
2002e6e4ca1SStefano Zampini KSP ksp;
2010cd8b6e2SStefano Zampini PetscBool valid;
2025eb5f4d6SAlp Dener
2035eb5f4d6SAlp Dener PetscFunctionBegin;
2049566063dSJacob Faibussowitsch PetscCall(TaoSetUp_BNK(tao));
2059566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp));
2060cd8b6e2SStefano Zampini PetscCall(PetscObjectHasFunction((PetscObject)ksp, "KSPCGSetRadius_C", &valid));
2073c859ba3SBarry 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);
2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2095eb5f4d6SAlp Dener }
2105eb5f4d6SAlp Dener
TaoSetFromOptions_BNTR(Tao tao,PetscOptionItems PetscOptionsObject)211ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_BNTR(Tao tao, PetscOptionItems PetscOptionsObject)
212d71ae5a4SJacob Faibussowitsch {
2139b6ef848SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
2149b6ef848SAlp Dener
2159b6ef848SAlp Dener PetscFunctionBegin;
216dbbe0bcdSBarry Smith PetscCall(TaoSetFromOptions_BNK(tao, PetscOptionsObject));
217e0ed867bSAlp Dener if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2199b6ef848SAlp Dener }
2209b6ef848SAlp Dener
2213850be85SAlp Dener /*MC
2223850be85SAlp Dener TAOBNTR - Bounded Newton Trust Region for nonlinear minimization with bound constraints.
2239b6ef848SAlp Dener
2243850be85SAlp Dener Options Database Keys:
2253850be85SAlp Dener + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
2263850be85SAlp Dener . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
2273850be85SAlp Dener . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
2283850be85SAlp Dener - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
2293850be85SAlp Dener
2303850be85SAlp Dener Level: beginner
2313850be85SAlp Dener M*/
TaoCreate_BNTR(Tao tao)232d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao)
233d71ae5a4SJacob Faibussowitsch {
234fed79b8eSAlp Dener TAO_BNK *bnk;
235fed79b8eSAlp Dener
236fed79b8eSAlp Dener PetscFunctionBegin;
2379566063dSJacob Faibussowitsch PetscCall(TaoCreate_BNK(tao));
238fed79b8eSAlp Dener tao->ops->solve = TaoSolve_BNTR;
2395eb5f4d6SAlp Dener tao->ops->setup = TaoSetUp_BNTR;
240e0ed867bSAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNTR;
241fed79b8eSAlp Dener
242fed79b8eSAlp Dener bnk = (TAO_BNK *)tao->data;
24366ed3702SAlp Dener bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
245fed79b8eSAlp Dener }
246