Actual source code: bnls.c
1: #include <../src/tao/bound/impls/bnk/bnk.h>
2: #include <petscksp.h>
4: /*
5: Implements Newton's Method with a line search approach for
6: solving 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_DIRECTION)
15: niter = 0
16: step_accepted = true
18: while niter < max_it
19: niter += 1
21: if needH
22: If max_cg_steps > 0
23: x_k, g_k, pg_k = TaoSolve(BNCG)
24: end
26: H_k = TaoComputeHessian(x_k)
27: if pc_type == BNK_PC_BFGS
28: add correction to BFGS approx
29: if scale_type == BNK_SCALE_AHESS
30: D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
31: scale BFGS with VecReciprocal(D)
32: end
33: end
34: needH = False
35: end
37: if pc_type = BNK_PC_BFGS
38: B_k = BFGS
39: else
40: B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
41: B_k = VecReciprocal(B_k)
42: end
43: w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
44: eps = min(eps, norm2(w))
45: determine the active and inactive index sets such that
46: L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
47: U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
48: F = {i : l_i = (x_k)_i = u_i}
49: A = {L + U + F}
50: IA = {i : i not in A}
52: generate the reduced system Hr_k dr_k = -gr_k for variables in IA
53: if p > 0
54: Hr_k += p*
55: end
56: if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
57: D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
58: scale BFGS with VecReciprocal(D)
59: end
60: solve Hr_k dr_k = -gr_k
61: set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
63: if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
64: dr_k = -BFGS*gr_k for variables in I
65: if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
66: reset the BFGS preconditioner
67: calculate scale delta and apply it to BFGS
68: dr_k = -BFGS*gr_k for variables in I
69: if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
70: dr_k = -gr_k for variables in I
71: end
72: end
73: end
75: x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch()
76: if ls_failed
77: f_{k+1} = f_k
78: x_{k+1} = x_k
79: g_{k+1} = g_k
80: pg_{k+1} = pg_k
81: terminate
82: else
83: pg_{k+1} = project(g_{k+1})
84: count the accepted step type (Newton, BFGS, scaled grad or grad)
85: end
87: check convergence at pg_{k+1}
88: end
89: */
91: PetscErrorCode TaoSolve_BNLS(Tao tao)
92: {
93: PetscErrorCode ierr;
94: TAO_BNK *bnk = (TAO_BNK *)tao->data;
95: KSPConvergedReason ksp_reason;
96: TaoLineSearchConvergedReason ls_reason;
98: PetscReal steplen = 1.0, resnorm;
99: PetscBool cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_TRUE;
100: PetscInt stepType;
103: /* Initialize the preconditioner, KSP solver and trust radius/line search */
104: tao->reason = TAO_CONTINUE_ITERATING;
105: TaoBNKInitialize(tao, bnk->init_type, &needH);
106: if (tao->reason != TAO_CONTINUE_ITERATING) return(0);
108: /* Have not converged; continue with Newton method */
109: while (tao->reason == TAO_CONTINUE_ITERATING) {
110: /* Call general purpose update function */
111: if (tao->ops->update) {
112: (*tao->ops->update)(tao, tao->niter, tao->user_update);
113: }
114: ++tao->niter;
116: if (needH && bnk->inactive_idx) {
117: /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
118: TaoBNKTakeCGSteps(tao, &cgTerminate);
119: if (cgTerminate) {
120: tao->reason = bnk->bncg->reason;
121: return(0);
122: }
123: /* Compute the hessian and update the BFGS preconditioner at the new iterate */
124: (*bnk->computehessian)(tao);
125: needH = PETSC_FALSE;
126: }
128: /* Use the common BNK kernel to compute the safeguarded Newton step (for inactive variables only) */
129: (*bnk->computestep)(tao, shift, &ksp_reason, &stepType);
130: TaoBNKSafeguardStep(tao, ksp_reason, &stepType);
132: /* Store current solution before it changes */
133: bnk->fold = bnk->f;
134: VecCopy(tao->solution, bnk->Xold);
135: VecCopy(tao->gradient, bnk->Gold);
136: VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);
138: /* Trigger the line search */
139: TaoBNKPerformLineSearch(tao, &stepType, &steplen, &ls_reason);
141: if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
142: /* Failed to find an improving point */
143: needH = PETSC_FALSE;
144: bnk->f = bnk->fold;
145: VecCopy(bnk->Xold, tao->solution);
146: VecCopy(bnk->Gold, tao->gradient);
147: VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);
148: steplen = 0.0;
149: tao->reason = TAO_DIVERGED_LS_FAILURE;
150: } else {
151: /* new iterate so we need to recompute the Hessian */
152: needH = PETSC_TRUE;
153: /* compute the projected gradient */
154: TaoBNKEstimateActiveSet(tao, bnk->as_type);
155: VecCopy(bnk->unprojected_gradient, tao->gradient);
156: VecISSet(tao->gradient, bnk->active_idx, 0.0);
157: TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);
158: /* update the trust radius based on the step length */
159: TaoBNKUpdateTrustRadius(tao, 0.0, 0.0, BNK_UPDATE_STEP, stepType, &stepAccepted);
160: /* count the accepted step type */
161: TaoBNKAddStepCounts(tao, stepType);
162: /* active BNCG recycling for next iteration */
163: TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE);
164: }
166: /* Check for termination */
167: VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);
168: VecNorm(bnk->W, NORM_2, &resnorm);
169: if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
170: TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);
171: TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);
172: (*tao->ops->convergencetest)(tao, tao->cnvP);
173: }
174: return(0);
175: }
177: /*------------------------------------------------------------*/
178: /*MC
179: TAOBNLS - Bounded Newton Line Search for nonlinear minimization with bound constraints.
181: Options Database Keys:
182: + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
183: . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
184: . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
185: - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
187: Level: beginner
188: M*/
189: PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao)
190: {
191: TAO_BNK *bnk;
195: TaoCreate_BNK(tao);
196: tao->ops->solve = TaoSolve_BNLS;
198: bnk = (TAO_BNK *)tao->data;
199: bnk->init_type = BNK_INIT_DIRECTION;
200: bnk->update_type = BNK_UPDATE_STEP; /* trust region updates based on line search step length */
201: return(0);
202: }