Actual source code: bntl.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: In this variant, the trust region failures trigger a line search with
9: the existing Newton step instead of re-solving the step with a
10: different radius.
12: ------------------------------------------------------------
14: x_0 = VecMedian(x_0)
15: f_0, g_0 = TaoComputeObjectiveAndGradient(x_0)
16: pg_0 = project(g_0)
17: check convergence at pg_0
18: needH = TaoBNKInitialize(default:BNK_INIT_INTERPOLATION)
19: niter = 0
20: step_accepted = true
22: while niter <= max_it
23: niter += 1
25: if needH
26: If max_cg_steps > 0
27: x_k, g_k, pg_k = TaoSolve(BNCG)
28: end
30: H_k = TaoComputeHessian(x_k)
31: if pc_type == BNK_PC_BFGS
32: add correction to BFGS approx
33: if scale_type == BNK_SCALE_AHESS
34: D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
35: scale BFGS with VecReciprocal(D)
36: end
37: end
38: needH = False
39: end
41: if pc_type = BNK_PC_BFGS
42: B_k = BFGS
43: else
44: B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
45: B_k = VecReciprocal(B_k)
46: end
47: w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
48: eps = min(eps, norm2(w))
49: determine the active and inactive index sets such that
50: L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
51: U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
52: F = {i : l_i = (x_k)_i = u_i}
53: A = {L + U + F}
54: IA = {i : i not in A}
56: generate the reduced system Hr_k dr_k = -gr_k for variables in IA
57: if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
58: D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
59: scale BFGS with VecReciprocal(D)
60: end
61: solve Hr_k dr_k = -gr_k
62: set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
64: x_{k+1} = VecMedian(x_k + d_k)
65: s = x_{k+1} - x_k
66: prered = dot(s, 0.5*gr_k - Hr_k*s)
67: f_{k+1} = TaoComputeObjective(x_{k+1})
68: actred = f_k - f_{k+1}
70: oldTrust = trust
71: step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION)
72: if step_accepted
73: g_{k+1} = TaoComputeGradient(x_{k+1})
74: pg_{k+1} = project(g_{k+1})
75: count the accepted Newton step
76: else
77: if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
78: dr_k = -BFGS*gr_k for variables in I
79: if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
80: reset the BFGS preconditioner
81: calculate scale delta and apply it to BFGS
82: dr_k = -BFGS*gr_k for variables in I
83: if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
84: dr_k = -gr_k for variables in I
85: end
86: end
87: end
89: x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch()
90: if ls_failed
91: f_{k+1} = f_k
92: x_{k+1} = x_k
93: g_{k+1} = g_k
94: pg_{k+1} = pg_k
95: terminate
96: else
97: pg_{k+1} = project(g_{k+1})
98: trust = oldTrust
99: trust = TaoBNKUpdateTrustRadius(BNK_UPDATE_STEP)
100: count the accepted step type (Newton, BFGS, scaled grad or grad)
101: end
102: end
104: check convergence at pg_{k+1}
105: end
106: */
108: PetscErrorCode TaoSolve_BNTL(Tao tao)
109: {
110: TAO_BNK *bnk = (TAO_BNK *)tao->data;
111: KSPConvergedReason ksp_reason;
112: TaoLineSearchConvergedReason ls_reason;
114: PetscReal oldTrust, prered, actred, steplen, resnorm;
115: PetscBool cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE;
116: PetscInt stepType, nDiff;
118: /* Initialize the preconditioner, KSP solver and trust radius/line search */
119: tao->reason = TAO_CONTINUE_ITERATING;
120: TaoBNKInitialize(tao, bnk->init_type, &needH);
121: if (tao->reason != TAO_CONTINUE_ITERATING) return 0;
123: /* Have not converged; continue with Newton method */
124: while (tao->reason == TAO_CONTINUE_ITERATING) {
125: /* Call general purpose update function */
126: if (tao->ops->update) {
127: PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
128: TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);
129: }
131: if (needH && bnk->inactive_idx) {
132: /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
133: TaoBNKTakeCGSteps(tao, &cgTerminate);
134: if (cgTerminate) {
135: tao->reason = bnk->bncg->reason;
136: return 0;
137: }
138: /* Compute the hessian and update the BFGS preconditioner at the new iterate */
139: (*bnk->computehessian)(tao);
140: needH = PETSC_FALSE;
141: }
143: /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
144: (*bnk->computestep)(tao, shift, &ksp_reason, &stepType);
146: /* Store current solution before it changes */
147: oldTrust = tao->trust;
148: bnk->fold = bnk->f;
149: VecCopy(tao->solution, bnk->Xold);
150: VecCopy(tao->gradient, bnk->Gold);
151: VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);
153: /* Temporarily accept the step and project it into the bounds */
154: VecAXPY(tao->solution, 1.0, tao->stepdirection);
155: TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution);
157: /* Check if the projection changed the step direction */
158: if (nDiff > 0) {
159: /* Projection changed the step, so we have to recompute the step and
160: the predicted reduction. Leave the trust radius unchanged. */
161: VecCopy(tao->solution, tao->stepdirection);
162: VecAXPY(tao->stepdirection, -1.0, bnk->Xold);
163: TaoBNKRecomputePred(tao, tao->stepdirection, &prered);
164: } else {
165: /* Step did not change, so we can just recover the pre-computed prediction */
166: KSPCGGetObjFcn(tao->ksp, &prered);
167: }
168: prered = -prered;
170: /* Compute the actual reduction and update the trust radius */
171: TaoComputeObjective(tao, tao->solution, &bnk->f);
173: actred = bnk->fold - bnk->f;
174: TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);
176: if (stepAccepted) {
177: /* Step is good, evaluate the gradient and the hessian */
178: steplen = 1.0;
179: needH = PETSC_TRUE;
180: ++bnk->newt;
181: TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);
182: TaoBNKEstimateActiveSet(tao, bnk->as_type);
183: VecCopy(bnk->unprojected_gradient, tao->gradient);
184: VecISSet(tao->gradient, bnk->active_idx, 0.0);
185: TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);
186: } else {
187: /* Trust-region rejected the step. Revert the solution. */
188: bnk->f = bnk->fold;
189: VecCopy(bnk->Xold, tao->solution);
190: /* Trigger the line search */
191: TaoBNKSafeguardStep(tao, ksp_reason, &stepType);
192: TaoBNKPerformLineSearch(tao, &stepType, &steplen, &ls_reason);
193: if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
194: /* Line search failed, revert solution and terminate */
195: stepAccepted = PETSC_FALSE;
196: needH = PETSC_FALSE;
197: bnk->f = bnk->fold;
198: VecCopy(bnk->Xold, tao->solution);
199: VecCopy(bnk->Gold, tao->gradient);
200: VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);
201: tao->trust = 0.0;
202: tao->reason = TAO_DIVERGED_LS_FAILURE;
203: } else {
204: /* new iterate so we need to recompute the Hessian */
205: needH = PETSC_TRUE;
206: /* compute the projected gradient */
207: TaoBNKEstimateActiveSet(tao, bnk->as_type);
208: VecCopy(bnk->unprojected_gradient, tao->gradient);
209: VecISSet(tao->gradient, bnk->active_idx, 0.0);
210: TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);
211: /* Line search succeeded so we should update the trust radius based on the LS step length */
212: tao->trust = oldTrust;
213: TaoBNKUpdateTrustRadius(tao, prered, actred, BNK_UPDATE_STEP, stepType, &stepAccepted);
214: /* count the accepted step type */
215: TaoBNKAddStepCounts(tao, stepType);
216: }
217: }
219: /* Check for termination */
220: VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);
221: VecNorm(bnk->W, NORM_2, &resnorm);
223: ++tao->niter;
224: TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);
225: TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);
226: PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
227: }
228: return 0;
229: }
231: /*------------------------------------------------------------*/
232: static PetscErrorCode TaoSetUp_BNTL(Tao tao)
233: {
234: KSP ksp;
235: PetscVoidFunction valid;
237: TaoSetUp_BNK(tao);
238: TaoGetKSP(tao, &ksp);
239: PetscObjectQueryFunction((PetscObject)ksp, "KSPCGSetRadius_C", &valid);
241: return 0;
242: }
244: /*------------------------------------------------------------*/
245: static PetscErrorCode TaoSetFromOptions_BNTL(Tao tao, PetscOptionItems *PetscOptionsObject)
246: {
247: TAO_BNK *bnk = (TAO_BNK *)tao->data;
249: TaoSetFromOptions_BNK(tao, PetscOptionsObject);
250: if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
251: return 0;
252: }
254: /*------------------------------------------------------------*/
255: /*MC
256: TAOBNTL - Bounded Newton Trust Region method with line-search fall-back for nonlinear
257: minimization with bound constraints.
259: Options Database Keys:
260: + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
261: . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
262: . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
263: - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
265: Level: beginner
266: M*/
267: PETSC_EXTERN PetscErrorCode TaoCreate_BNTL(Tao tao)
268: {
269: TAO_BNK *bnk;
271: TaoCreate_BNK(tao);
272: tao->ops->solve = TaoSolve_BNTL;
273: tao->ops->setup = TaoSetUp_BNTL;
274: tao->ops->setfromoptions = TaoSetFromOptions_BNTL;
276: bnk = (TAO_BNK *)tao->data;
277: bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
278: return 0;
279: }