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