xref: /petsc/src/tao/bound/impls/bnk/bntl.c (revision e031d6f587cbb9f14a00ce52e5e14087387d741b)
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  initialize trust radius (default: BNK_INIT_INTERPOLATION)
11  x_0 = VecMedian(x_0)
12  f_0, g_0 = TaoComputeObjectiveAndGradient(x_0)
13  pg_0 = VecBoundGradientProjection(g_0)
14  check convergence at pg_0
15  niter = 0
16  step_accepted = true
17 
18  while niter <= max_it
19     niter += 1
20     H_k = TaoComputeHessian(x_k)
21     if pc_type == BNK_PC_BFGS
22       add correction to BFGS approx
23       if scale_type == BNK_SCALE_AHESS
24         D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
25         scale BFGS with VecReciprocal(D)
26       end
27     end
28 
29     if pc_type = BNK_PC_BFGS
30       B_k = BFGS
31     else
32       B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
33       B_k = VecReciprocal(B_k)
34     end
35     w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
36     eps = min(eps, norm2(w))
37     determine the active and inactive index sets such that
38       L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
39       U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
40       F = {i : l_i = (x_k)_i = u_i}
41       A = {L + U + F}
42       I = {i : i not in A}
43 
44     generate the reduced system Hr_k dr_k = -gr_k for variables in I
45     if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
46       D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
47       scale BFGS with VecReciprocal(D)
48     end
49     solve Hr_k dr_k = -gr_k
50     set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
51 
52     x_{k+1} = VecMedian(x_k + d_k)
53     s = x_{k+1} - x_k
54     prered = dot(s, 0.5*gr_k - Hr_k*s)
55     f_{k+1} = TaoComputeObjective(x_{k+1})
56     actred = f_k - f_{k+1}
57 
58     oldTrust = trust
59     step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION)
60     if step_accepted
61       g_{k+1} = TaoComputeGradient(x_{k+1})
62       pg_{k+1} = VecBoundGradientProjection(g_{k+1})
63       count the accepted Newton step
64     else
65       if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
66         dr_k = -BFGS*gr_k for variables in I
67         if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
68           reset the BFGS preconditioner
69           calculate scale delta and apply it to BFGS
70           dr_k = -BFGS*gr_k for variables in I
71           if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
72             dr_k = -gr_k for variables in I
73           end
74         end
75       end
76 
77       x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch()
78       if ls_failed
79         f_{k+1} = f_k
80         x_{k+1} = x_k
81         g_{k+1} = g_k
82         pg_{k+1} = pg_k
83         terminate
84       else
85         pg_{k+1} = VecBoundGradientProjection(g_{k+1})
86         trust = oldTrust
87         trust = TaoBNKUpdateTrustRadius(BNK_UPDATE_STEP)
88         count the accepted step type (Newton, BFGS, scaled grad or grad)
89       end
90     end
91 
92     check convergence at pg_{k+1}
93  end
94 */
95 
96 static PetscErrorCode TaoSolve_BNTL(Tao tao)
97 {
98   PetscErrorCode               ierr;
99   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
100   KSPConvergedReason           ksp_reason;
101   TaoLineSearchConvergedReason ls_reason;
102 
103   PetscReal                    resnorm, oldTrust, prered, actred, stepNorm, steplen;
104   PetscBool                    cgTerminate, stepAccepted = PETSC_TRUE, shift = PETSC_FALSE;
105   PetscInt                     stepType;
106 
107   PetscFunctionBegin;
108   /* Initialize the preconditioner, KSP solver and trust radius/line search */
109   tao->reason = TAO_CONTINUE_ITERATING;
110   ierr = TaoBNKInitialize(tao, bnk->init_type, &stepAccepted);CHKERRQ(ierr);
111   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
112 
113   /* Have not converged; continue with Newton method */
114   while (tao->reason == TAO_CONTINUE_ITERATING) {
115     tao->niter++;
116     tao->ksp_its=0;
117 
118     /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
119     ierr = TaoBNKTakeCGSteps(tao, &cgTerminate);CHKERRQ(ierr);
120     if (cgTerminate) {
121       tao->reason = bnk->bncg->reason;
122       PetscFunctionReturn(0);
123     }
124 
125     /* Compute the hessian, update the BFGS preconditioner and estimate the active-set at the new iterate */
126     if (stepAccepted) {
127       ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
128       ierr = TaoBNKEstimateActiveSet(tao);CHKERRQ(ierr);
129     }
130 
131     /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
132     ierr = TaoBNKComputeStep(tao, shift, &ksp_reason);CHKERRQ(ierr);
133 
134     /* Store current solution before it changes */
135     oldTrust = tao->trust;
136     bnk->fold = bnk->f;
137     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
138     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
139     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);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 = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
144 
145     /* Check if the projection changed the step direction */
146     ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
147     ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
148     ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr);
149     if (stepNorm != bnk->dnorm) {
150       /* Projection changed the step, so we have to recompute predicted reduction.
151          However, we deliberately do not change the step norm and the trust radius
152          in order for the safeguard to more closely mimic a piece-wise linesearch
153          along the bounds. */
154       ierr = TaoBNKRecomputePred(tao, tao->stepdirection, &prered);CHKERRQ(ierr);
155     } else {
156       /* Step did not change, so we can just recover the pre-computed prediction */
157       ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
158     }
159     prered = -prered;
160 
161     /* Compute the actual reduction and update the trust radius */
162     ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
163     actred = bnk->fold - bnk->f;
164     ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr);
165 
166     if (stepAccepted) {
167       /* Step is good, evaluate the gradient and the hessian */
168       steplen = 1.0;
169       ++bnk->newt;
170       ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
171       ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
172     } else {
173       /* Trust-region rejected the step. Revert the solution. */
174       bnk->f = bnk->fold;
175       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
176 
177       /* Trigger the line search */
178       ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr);
179       ierr = TaoBNKPerformLineSearch(tao, stepType, &steplen, &ls_reason);CHKERRQ(ierr);
180       if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
181         /* Line search failed, revert solution and terminate */
182         stepAccepted = PETSC_FALSE;
183         bnk->f = bnk->fold;
184         ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
185         ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
186         ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
187         tao->trust = 0.0;
188         tao->reason = TAO_DIVERGED_LS_FAILURE;
189       } else {
190         /* compute the projected gradient */
191         ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
192         ierr = VecNorm(tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
193         if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
194         /* Line search succeeded so we should update the trust radius based on the LS step length */
195         tao->trust = oldTrust;
196         ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, BNK_UPDATE_STEP, stepType, &stepAccepted);CHKERRQ(ierr);
197         /* count the accepted step type */
198         ierr = TaoBNKAddStepCounts(tao, stepType);CHKERRQ(ierr);
199       }
200     }
201 
202     /*  Check for termination */
203     ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr);
204     ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr);
205     ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
206     ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr);
207     ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr);
208   }
209   PetscFunctionReturn(0);
210 }
211 
212 /*------------------------------------------------------------*/
213 
214 PETSC_INTERN PetscErrorCode TaoSetUp_BNTL(Tao tao)
215 {
216   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
217   PetscErrorCode ierr;
218 
219   PetscFunctionBegin;
220   ierr = TaoSetUp_BNK(tao);CHKERRQ(ierr);
221   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)");
222   PetscFunctionReturn(0);
223 }
224 
225 /*------------------------------------------------------------*/
226 
227 PETSC_INTERN PetscErrorCode TaoCreate_BNTL(Tao tao)
228 {
229   TAO_BNK        *bnk;
230   PetscErrorCode ierr;
231 
232   PetscFunctionBegin;
233   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
234   tao->ops->solve=TaoSolve_BNTL;
235   tao->ops->setup=TaoSetUp_BNTL;
236 
237   bnk = (TAO_BNK *)tao->data;
238   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
239   PetscFunctionReturn(0);
240 }