Actual source code: bntr.c
1: #include <../src/tao/bound/impls/bnk/bnk.h>
2: #include <petscksp.h>
4: /*
5: Implements Newton's Method with a trust region approach for solving
6: bound constrained minimization problems.
8: ------------------------------------------------------------
10: x_0 = VecMedian(x_0)
11: f_0, g_0= TaoComputeObjectiveAndGradient(x_0)
12: pg_0 = project(g_0)
13: check convergence at pg_0
14: needH = TaoBNKInitialize(default:BNK_INIT_INTERPOLATION)
15: niter = 0
16: step_accepted = false
18: while niter <= max_it
20: if needH
21: If max_cg_steps > 0
22: x_k, g_k, pg_k = TaoSolve(BNCG)
23: end
25: H_k = TaoComputeHessian(x_k)
26: if pc_type == BNK_PC_BFGS
27: add correction to BFGS approx
28: if scale_type == BNK_SCALE_AHESS
29: D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
30: scale BFGS with VecReciprocal(D)
31: end
32: end
33: needH = False
34: end
36: if pc_type = BNK_PC_BFGS
37: B_k = BFGS
38: else
39: B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
40: B_k = VecReciprocal(B_k)
41: end
42: w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
43: eps = min(eps, norm2(w))
44: determine the active and inactive index sets such that
45: L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
46: U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
47: F = {i : l_i = (x_k)_i = u_i}
48: A = {L + U + F}
49: IA = {i : i not in A}
51: generate the reduced system Hr_k dr_k = -gr_k for variables in IA
52: if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
53: D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
54: scale BFGS with VecReciprocal(D)
55: end
57: while !stepAccepted
58: solve Hr_k dr_k = -gr_k
59: set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
61: x_{k+1} = VecMedian(x_k + d_k)
62: s = x_{k+1} - x_k
63: prered = dot(s, 0.5*gr_k - Hr_k*s)
64: f_{k+1} = TaoComputeObjective(x_{k+1})
65: actred = f_k - f_{k+1}
67: oldTrust = trust
68: step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION)
69: if step_accepted
70: g_{k+1} = TaoComputeGradient(x_{k+1})
71: pg_{k+1} = project(g_{k+1})
72: count the accepted Newton step
73: needH = True
74: else
75: f_{k+1} = f_k
76: x_{k+1} = x_k
77: g_{k+1} = g_k
78: pg_{k+1} = pg_k
79: if trust == oldTrust
80: terminate because we cannot shrink the radius any further
81: end
82: end
84: end
85: check convergence at pg_{k+1}
86: niter += 1
88: end
89: */
91: PetscErrorCode TaoSolve_BNTR(Tao tao)
92: {
93: TAO_BNK *bnk = (TAO_BNK *)tao->data;
94: KSPConvergedReason ksp_reason;
96: PetscReal oldTrust, prered, actred, steplen = 0.0, resnorm;
97: PetscBool cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE;
98: PetscInt stepType, nDiff;
100: PetscFunctionBegin;
101: /* Initialize the preconditioner, KSP solver and trust radius/line search */
102: tao->reason = TAO_CONTINUE_ITERATING;
103: PetscCall(TaoBNKInitialize(tao, bnk->init_type, &needH));
104: if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
106: /* Have not converged; continue with Newton method */
107: while (tao->reason == TAO_CONTINUE_ITERATING) {
108: /* Call general purpose update function */
109: if (tao->ops->update) {
110: PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
111: PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient));
112: }
114: if (needH && bnk->inactive_idx) {
115: /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
116: PetscCall(TaoBNKTakeCGSteps(tao, &cgTerminate));
117: if (cgTerminate) {
118: tao->reason = bnk->bncg->reason;
119: PetscFunctionReturn(PETSC_SUCCESS);
120: }
121: /* Compute the hessian and update the BFGS preconditioner at the new iterate */
122: PetscCall((*bnk->computehessian)(tao));
123: needH = PETSC_FALSE;
124: }
126: /* Store current solution before it changes */
127: bnk->fold = bnk->f;
128: PetscCall(VecCopy(tao->solution, bnk->Xold));
129: PetscCall(VecCopy(tao->gradient, bnk->Gold));
130: PetscCall(VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old));
132: /* Enter into trust region loops */
133: stepAccepted = PETSC_FALSE;
134: while (!stepAccepted && tao->reason == TAO_CONTINUE_ITERATING) {
135: tao->ksp_its = 0;
137: /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
138: PetscCall((*bnk->computestep)(tao, shift, &ksp_reason, &stepType));
140: /* Temporarily accept the step and project it into the bounds */
141: PetscCall(VecAXPY(tao->solution, 1.0, tao->stepdirection));
142: PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
144: /* Check if the projection changed the step direction */
145: if (nDiff > 0) {
146: /* Projection changed the step, so we have to recompute the step and
147: the predicted reduction. Leave the trust radius unchanged. */
148: PetscCall(VecCopy(tao->solution, tao->stepdirection));
149: PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold));
150: PetscCall(TaoBNKRecomputePred(tao, tao->stepdirection, &prered));
151: } else {
152: /* Step did not change, so we can just recover the pre-computed prediction */
153: PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
154: }
155: prered = -prered;
157: /* Compute the actual reduction and update the trust radius */
158: PetscCall(TaoComputeObjective(tao, tao->solution, &bnk->f));
159: PetscCheck(!PetscIsInfOrNanReal(bnk->f), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
160: actred = bnk->fold - bnk->f;
161: oldTrust = tao->trust;
162: PetscCall(TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted));
164: if (stepAccepted) {
165: /* Step is good, evaluate the gradient and flip the need-Hessian switch */
166: steplen = 1.0;
167: needH = PETSC_TRUE;
168: ++bnk->newt;
169: PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient));
170: PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
171: PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
172: PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
173: PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
174: } else {
175: /* Step is bad, revert old solution and re-solve with new radius*/
176: steplen = 0.0;
177: needH = PETSC_FALSE;
178: bnk->f = bnk->fold;
179: PetscCall(VecCopy(bnk->Xold, tao->solution));
180: PetscCall(VecCopy(bnk->Gold, tao->gradient));
181: PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));
182: if (oldTrust == tao->trust) {
183: /* Can't change the radius anymore so just terminate */
184: tao->reason = TAO_DIVERGED_TR_REDUCTION;
185: }
186: }
187: }
188: /* Check for termination */
189: PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
190: PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
191: PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
192: ++tao->niter;
193: PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
194: PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen));
195: PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
196: }
197: PetscFunctionReturn(PETSC_SUCCESS);
198: }
200: /*------------------------------------------------------------*/
201: static PetscErrorCode TaoSetUp_BNTR(Tao tao)
202: {
203: KSP ksp;
204: PetscVoidFunction valid;
206: PetscFunctionBegin;
207: PetscCall(TaoSetUp_BNK(tao));
208: PetscCall(TaoGetKSP(tao, &ksp));
209: PetscCall(PetscObjectQueryFunction((PetscObject)ksp, "KSPCGSetRadius_C", &valid));
210: 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);
211: PetscFunctionReturn(PETSC_SUCCESS);
212: }
214: /*------------------------------------------------------------*/
216: static PetscErrorCode TaoSetFromOptions_BNTR(Tao tao, PetscOptionItems *PetscOptionsObject)
217: {
218: TAO_BNK *bnk = (TAO_BNK *)tao->data;
220: PetscFunctionBegin;
221: PetscCall(TaoSetFromOptions_BNK(tao, PetscOptionsObject));
222: if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
223: PetscFunctionReturn(PETSC_SUCCESS);
224: }
226: /*------------------------------------------------------------*/
227: /*MC
228: TAOBNTR - Bounded Newton Trust Region for nonlinear minimization with bound constraints.
230: Options Database Keys:
231: + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
232: . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
233: . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
234: - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
236: Level: beginner
237: M*/
238: PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao)
239: {
240: TAO_BNK *bnk;
242: PetscFunctionBegin;
243: PetscCall(TaoCreate_BNK(tao));
244: tao->ops->solve = TaoSolve_BNTR;
245: tao->ops->setup = TaoSetUp_BNTR;
246: tao->ops->setfromoptions = TaoSetFromOptions_BNTR;
248: bnk = (TAO_BNK *)tao->data;
249: bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
250: PetscFunctionReturn(PETSC_SUCCESS);
251: }