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