xref: /petsc/src/tao/bound/impls/bnk/bntr.c (revision 65f8aed5f7eaa1e2ef2ddeffe666264e0669c876)
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     ++tao->niter;
110 
111     if (needH && bnk->inactive_idx) {
112       /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
113       ierr = TaoBNKTakeCGSteps(tao, &cgTerminate);CHKERRQ(ierr);
114       if (cgTerminate) {
115         tao->reason = bnk->bncg->reason;
116         PetscFunctionReturn(0);
117       }
118       /* Compute the hessian and update the BFGS preconditioner at the new iterate */
119       ierr = (*bnk->computehessian)(tao);CHKERRQ(ierr);
120       needH = PETSC_FALSE;
121     }
122 
123     /* Store current solution before it changes */
124     bnk->fold = bnk->f;
125     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
126     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
127     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
128 
129     /* Enter into trust region loops */
130     stepAccepted = PETSC_FALSE;
131     while (!stepAccepted && tao->reason == TAO_CONTINUE_ITERATING) {
132       tao->ksp_its=0;
133 
134       /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
135       ierr = (*bnk->computestep)(tao, shift, &ksp_reason, &stepType);CHKERRQ(ierr);
136 
137       /* Temporarily accept the step and project it into the bounds */
138       ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
139       ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
140 
141       /* Check if the projection changed the step direction */
142       if (nDiff > 0) {
143         /* Projection changed the step, so we have to recompute the step and
144            the predicted reduction. Leave the trust radius unchanged. */
145         ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
146         ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
147         ierr = TaoBNKRecomputePred(tao, tao->stepdirection, &prered);CHKERRQ(ierr);
148       } else {
149         /* Step did not change, so we can just recover the pre-computed prediction */
150         ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
151       }
152       prered = -prered;
153 
154       /* Compute the actual reduction and update the trust radius */
155       ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
156       if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
157       actred = bnk->fold - bnk->f;
158       oldTrust = tao->trust;
159       ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr);
160 
161       if (stepAccepted) {
162         /* Step is good, evaluate the gradient and flip the need-Hessian switch */
163         steplen = 1.0;
164         needH = PETSC_TRUE;
165         ++bnk->newt;
166         ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
167         ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
168         ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
169         ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
170         ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
171       } else {
172         /* Step is bad, revert old solution and re-solve with new radius*/
173         steplen = 0.0;
174         needH = PETSC_FALSE;
175         bnk->f = bnk->fold;
176         ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
177         ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
178         ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
179         if (oldTrust == tao->trust) {
180           /* Can't change the radius anymore so just terminate */
181           tao->reason = TAO_DIVERGED_TR_REDUCTION;
182         }
183       }
184 
185       /*  Check for termination */
186       ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
187       ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
188       if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
189       ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
190       ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr);
191       ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr);
192     }
193   }
194   PetscFunctionReturn(0);
195 }
196 
197 /*------------------------------------------------------------*/
198 
199 static PetscErrorCode TaoSetFromOptions_BNTR(PetscOptionItems *PetscOptionsObject,Tao tao)
200 {
201   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
202   PetscErrorCode ierr;
203 
204   PetscFunctionBegin;
205   ierr = TaoSetFromOptions_BNK(PetscOptionsObject, tao);CHKERRQ(ierr);
206   if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
207   if (!bnk->is_nash && !bnk->is_stcg && !bnk->is_gltr) SETERRQ(PETSC_COMM_SELF,1,"Must use a trust-region CG method for KSP (KSPNASH, KSPSTCG, KSPGLTR)");
208   PetscFunctionReturn(0);
209 }
210 
211 /*------------------------------------------------------------*/
212 
213 PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao)
214 {
215   TAO_BNK        *bnk;
216   PetscErrorCode ierr;
217 
218   PetscFunctionBegin;
219   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
220   tao->ops->solve=TaoSolve_BNTR;
221   tao->ops->setfromoptions=TaoSetFromOptions_BNTR;
222 
223   bnk = (TAO_BNK *)tao->data;
224   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
225   PetscFunctionReturn(0);
226 }
227