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  ------------------------------------------------------------
13 
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
21 
22  while niter <= max_it
23     niter += 1
24 
25     if needH
26       If max_cg_steps > 0
27         x_k, g_k, pg_k = TaoSolve(BNCG)
28       end
29 
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
40 
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}
55 
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
63 
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}
69 
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
88 
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
103 
104     check convergence at pg_{k+1}
105  end
106 */
107 
TaoSolve_BNTL(Tao tao)108 PetscErrorCode TaoSolve_BNTL(Tao tao)
109 {
110   PetscErrorCode               ierr;
111   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
112   KSPConvergedReason           ksp_reason;
113   TaoLineSearchConvergedReason ls_reason;
114 
115   PetscReal                    oldTrust, prered, actred, steplen, resnorm;
116   PetscBool                    cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE;
117   PetscInt                     stepType, nDiff;
118 
119   PetscFunctionBegin;
120   /* Initialize the preconditioner, KSP solver and trust radius/line search */
121   tao->reason = TAO_CONTINUE_ITERATING;
122   ierr = TaoBNKInitialize(tao, bnk->init_type, &needH);CHKERRQ(ierr);
123   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
124 
125   /* Have not converged; continue with Newton method */
126   while (tao->reason == TAO_CONTINUE_ITERATING) {
127     /* Call general purpose update function */
128     if (tao->ops->update) {
129       ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr);
130     }
131     ++tao->niter;
132 
133     if (needH && bnk->inactive_idx) {
134       /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
135       ierr = TaoBNKTakeCGSteps(tao, &cgTerminate);CHKERRQ(ierr);
136       if (cgTerminate) {
137         tao->reason = bnk->bncg->reason;
138         PetscFunctionReturn(0);
139       }
140       /* Compute the hessian and update the BFGS preconditioner at the new iterate */
141       ierr = (*bnk->computehessian)(tao);CHKERRQ(ierr);
142       needH = PETSC_FALSE;
143     }
144 
145     /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
146     ierr = (*bnk->computestep)(tao, shift, &ksp_reason, &stepType);CHKERRQ(ierr);
147 
148     /* Store current solution before it changes */
149     oldTrust = tao->trust;
150     bnk->fold = bnk->f;
151     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
152     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
153     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
154 
155     /* Temporarily accept the step and project it into the bounds */
156     ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
157     ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
158 
159     /* Check if the projection changed the step direction */
160     if (nDiff > 0) {
161       /* Projection changed the step, so we have to recompute the step and
162          the predicted reduction. Leave the trust radius unchanged. */
163       ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
164       ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
165       ierr = TaoBNKRecomputePred(tao, tao->stepdirection, &prered);CHKERRQ(ierr);
166     } else {
167       /* Step did not change, so we can just recover the pre-computed prediction */
168       ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
169     }
170     prered = -prered;
171 
172     /* Compute the actual reduction and update the trust radius */
173     ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
174     if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
175     actred = bnk->fold - bnk->f;
176     ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr);
177 
178     if (stepAccepted) {
179       /* Step is good, evaluate the gradient and the hessian */
180       steplen = 1.0;
181       needH = PETSC_TRUE;
182       ++bnk->newt;
183       ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
184       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
185       ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
186       ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
187       ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
188     } else {
189       /* Trust-region rejected the step. Revert the solution. */
190       bnk->f = bnk->fold;
191       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
192       /* Trigger the line search */
193       ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr);
194       ierr = TaoBNKPerformLineSearch(tao, &stepType, &steplen, &ls_reason);CHKERRQ(ierr);
195       if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
196         /* Line search failed, revert solution and terminate */
197         stepAccepted = PETSC_FALSE;
198         needH = PETSC_FALSE;
199         bnk->f = bnk->fold;
200         ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
201         ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
202         ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
203         tao->trust = 0.0;
204         tao->reason = TAO_DIVERGED_LS_FAILURE;
205       } else {
206         /* new iterate so we need to recompute the Hessian */
207         needH = PETSC_TRUE;
208         /* compute the projected gradient */
209         ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
210         ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
211         ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
212         ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
213         /* Line search succeeded so we should update the trust radius based on the LS step length */
214         tao->trust = oldTrust;
215         ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, BNK_UPDATE_STEP, stepType, &stepAccepted);CHKERRQ(ierr);
216         /* count the accepted step type */
217         ierr = TaoBNKAddStepCounts(tao, stepType);CHKERRQ(ierr);
218       }
219     }
220 
221     /*  Check for termination */
222     ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
223     ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
224     if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
225     ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
226     ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr);
227     ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr);
228   }
229   PetscFunctionReturn(0);
230 }
231 
232 /*------------------------------------------------------------*/
TaoSetFromOptions_BNTL(PetscOptionItems * PetscOptionsObject,Tao tao)233 static PetscErrorCode TaoSetFromOptions_BNTL(PetscOptionItems *PetscOptionsObject,Tao tao)
234 {
235   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
236   PetscErrorCode ierr;
237 
238   PetscFunctionBegin;
239   ierr = TaoSetFromOptions_BNK(PetscOptionsObject, tao);CHKERRQ(ierr);
240   if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
241   if (!bnk->is_nash && !bnk->is_stcg && !bnk->is_gltr) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Must use a trust-region CG method for KSP (KSPNASH, KSPSTCG, KSPGLTR)");
242   PetscFunctionReturn(0);
243 }
244 
245 /*------------------------------------------------------------*/
246 /*MC
247   TAOBNTL - Bounded Newton Trust Region method with line-search fall-back for nonlinear
248             minimization with bound constraints.
249 
250   Options Database Keys:
251   + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
252   . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
253   . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
254   - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
255 
256   Level: beginner
257 M*/
TaoCreate_BNTL(Tao tao)258 PETSC_EXTERN PetscErrorCode TaoCreate_BNTL(Tao tao)
259 {
260   TAO_BNK        *bnk;
261   PetscErrorCode ierr;
262 
263   PetscFunctionBegin;
264   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
265   tao->ops->solve=TaoSolve_BNTL;
266   tao->ops->setfromoptions=TaoSetFromOptions_BNTL;
267 
268   bnk = (TAO_BNK *)tao->data;
269   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
270   PetscFunctionReturn(0);
271 }
272