xref: /petsc/src/tao/bound/impls/bnk/bntr.c (revision c3b5f7ba6bc5ce25a01a67bb37ba5d34b02bbbd7)
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  ------------------------------------------------------------
9 
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
17 
18  while niter <= max_it
19 
20     if needH
21       If max_cg_steps > 0
22         x_k, g_k, pg_k = TaoSolve(BNCG)
23       end
24 
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
35 
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}
50 
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
56 
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
60 
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}
66 
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
83 
84       check convergence at pg_{k+1}
85     end
86     niter += 1
87 
88  end
89 */
90 
91 PetscErrorCode TaoSolve_BNTR(Tao tao)
92 {
93   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
94   KSPConvergedReason           ksp_reason;
95 
96   PetscReal                    oldTrust, prered, actred, steplen, resnorm;
97   PetscBool                    cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE;
98   PetscInt                     stepType, nDiff;
99 
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(0);
105 
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       PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update));
111     }
112 
113     if (needH && bnk->inactive_idx) {
114       /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
115       PetscCall(TaoBNKTakeCGSteps(tao, &cgTerminate));
116       if (cgTerminate) {
117         tao->reason = bnk->bncg->reason;
118         PetscFunctionReturn(0);
119       }
120       /* Compute the hessian and update the BFGS preconditioner at the new iterate */
121       PetscCall((*bnk->computehessian)(tao));
122       needH = PETSC_FALSE;
123     }
124 
125     /* Store current solution before it changes */
126     bnk->fold = bnk->f;
127     PetscCall(VecCopy(tao->solution, bnk->Xold));
128     PetscCall(VecCopy(tao->gradient, bnk->Gold));
129     PetscCall(VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old));
130 
131     /* Enter into trust region loops */
132     stepAccepted = PETSC_FALSE;
133     while (!stepAccepted && tao->reason == TAO_CONTINUE_ITERATING) {
134       tao->ksp_its=0;
135 
136       /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
137       PetscCall((*bnk->computestep)(tao, shift, &ksp_reason, &stepType));
138 
139       /* Temporarily accept the step and project it into the bounds */
140       PetscCall(VecAXPY(tao->solution, 1.0, tao->stepdirection));
141       PetscCall(TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution));
142 
143       /* Check if the projection changed the step direction */
144       if (nDiff > 0) {
145         /* Projection changed the step, so we have to recompute the step and
146            the predicted reduction. Leave the trust radius unchanged. */
147         PetscCall(VecCopy(tao->solution, tao->stepdirection));
148         PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold));
149         PetscCall(TaoBNKRecomputePred(tao, tao->stepdirection, &prered));
150       } else {
151         /* Step did not change, so we can just recover the pre-computed prediction */
152         PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
153       }
154       prered = -prered;
155 
156       /* Compute the actual reduction and update the trust radius */
157       PetscCall(TaoComputeObjective(tao, tao->solution, &bnk->f));
158       PetscCheck(!PetscIsInfOrNanReal(bnk->f),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
159       actred = bnk->fold - bnk->f;
160       oldTrust = tao->trust;
161       PetscCall(TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted));
162 
163       if (stepAccepted) {
164         /* Step is good, evaluate the gradient and flip the need-Hessian switch */
165         steplen = 1.0;
166         needH = PETSC_TRUE;
167         ++bnk->newt;
168         PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient));
169         PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
170         PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
171         PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
172         PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
173       } else {
174         /* Step is bad, revert old solution and re-solve with new radius*/
175         steplen = 0.0;
176         needH = PETSC_FALSE;
177         bnk->f = bnk->fold;
178         PetscCall(VecCopy(bnk->Xold, tao->solution));
179         PetscCall(VecCopy(bnk->Gold, tao->gradient));
180         PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));
181         if (oldTrust == tao->trust) {
182           /* Can't change the radius anymore so just terminate */
183           tao->reason = TAO_DIVERGED_TR_REDUCTION;
184         }
185       }
186 
187       /*  Check for termination */
188       PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
189       PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
190       PetscCheck(!PetscIsInfOrNanReal(resnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
191       ++tao->niter; /* temporarily increment the iteration number */
192       PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
193       PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen));
194       PetscCall((*tao->ops->convergencetest)(tao, tao->cnvP));
195       --tao->niter; /* decrement */
196     }
197     ++tao->niter; /* finished iterating */
198   }
199   PetscFunctionReturn(0);
200 }
201 
202 /*------------------------------------------------------------*/
203 static PetscErrorCode TaoSetUp_BNTR(Tao tao)
204 {
205   KSP               ksp;
206   PetscVoidFunction valid;
207 
208   PetscFunctionBegin;
209   PetscCall(TaoSetUp_BNK(tao));
210   PetscCall(TaoGetKSP(tao,&ksp));
211   PetscCall(PetscObjectQueryFunction((PetscObject)ksp,"KSPCGSetRadius_C",&valid));
212   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);
213   PetscFunctionReturn(0);
214 }
215 
216 /*------------------------------------------------------------*/
217 
218 static PetscErrorCode TaoSetFromOptions_BNTR(PetscOptionItems *PetscOptionsObject,Tao tao)
219 {
220   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
221 
222   PetscFunctionBegin;
223   PetscCall(TaoSetFromOptions_BNK(PetscOptionsObject, tao));
224   if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
225   PetscFunctionReturn(0);
226 }
227 
228 /*------------------------------------------------------------*/
229 /*MC
230   TAOBNTR - Bounded Newton Trust Region for nonlinear minimization with bound constraints.
231 
232   Options Database Keys:
233 + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
234 . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
235 . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
236 - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
237 
238   Level: beginner
239 M*/
240 PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao)
241 {
242   TAO_BNK        *bnk;
243 
244   PetscFunctionBegin;
245   PetscCall(TaoCreate_BNK(tao));
246   tao->ops->solve=TaoSolve_BNTR;
247   tao->ops->setup=TaoSetUp_BNTR;
248   tao->ops->setfromoptions=TaoSetFromOptions_BNTR;
249 
250   bnk = (TAO_BNK *)tao->data;
251   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
252   PetscFunctionReturn(0);
253 }
254