xref: /petsc/src/tao/bound/impls/bnk/bntl.c (revision 0875260307e97403c513ee2ad2ca062f01ff5a8a)
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 and update the BFGS preconditioner at the new iterate */
126     if (stepAccepted) {ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);}
127 
128     /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
129     ierr = TaoBNKComputeStep(tao, shift, &ksp_reason);CHKERRQ(ierr);
130 
131     /* Store current solution before it changes */
132     oldTrust = tao->trust;
133     bnk->fold = bnk->f;
134     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
135     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
136     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
137 
138     /* Temporarily accept the step and project it into the bounds */
139     ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
140     ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
141 
142     /* Check if the projection changed the step direction */
143     ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
144     ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
145     ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr);
146     if (stepNorm != bnk->dnorm) {
147       /* Projection changed the step, so we have to recompute predicted reduction.
148          However, we deliberately do not change the step norm and the trust radius
149          in order for the safeguard to more closely mimic a piece-wise linesearch
150          along the bounds. */
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     actred = bnk->fold - bnk->f;
161     ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr);
162 
163     if (stepAccepted) {
164       /* Step is good, evaluate the gradient and the hessian */
165       steplen = 1.0;
166       ++bnk->newt;
167       ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
168       ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
169     } else {
170       /* Trust-region rejected the step. Revert the solution. */
171       bnk->f = bnk->fold;
172       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
173 
174       /* Trigger the line search */
175       ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr);
176       ierr = TaoBNKPerformLineSearch(tao, stepType, &steplen, &ls_reason);CHKERRQ(ierr);
177       if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
178         /* Line search failed, revert solution and terminate */
179         stepAccepted = PETSC_FALSE;
180         bnk->f = bnk->fold;
181         ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
182         ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
183         ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
184         tao->trust = 0.0;
185         tao->reason = TAO_DIVERGED_LS_FAILURE;
186       } else {
187         /* compute the projected gradient */
188         ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
189         ierr = VecNorm(tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
190         if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
191         /* Line search succeeded so we should update the trust radius based on the LS step length */
192         tao->trust = oldTrust;
193         ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, BNK_UPDATE_STEP, stepType, &stepAccepted);CHKERRQ(ierr);
194         /* count the accepted step type */
195         ierr = TaoBNKAddStepCounts(tao, stepType);CHKERRQ(ierr);
196       }
197     }
198 
199     /*  Check for termination */
200     ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr);
201     ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr);
202     ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
203     ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr);
204     ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr);
205   }
206   PetscFunctionReturn(0);
207 }
208 
209 /*------------------------------------------------------------*/
210 
211 PETSC_INTERN PetscErrorCode TaoSetUp_BNTL(Tao tao)
212 {
213   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
214   PetscErrorCode ierr;
215 
216   PetscFunctionBegin;
217   ierr = TaoSetUp_BNK(tao);CHKERRQ(ierr);
218   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)");
219   PetscFunctionReturn(0);
220 }
221 
222 /*------------------------------------------------------------*/
223 
224 PETSC_INTERN PetscErrorCode TaoCreate_BNTL(Tao tao)
225 {
226   TAO_BNK        *bnk;
227   PetscErrorCode ierr;
228 
229   PetscFunctionBegin;
230   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
231   tao->ops->solve=TaoSolve_BNTL;
232   tao->ops->setup=TaoSetUp_BNTL;
233 
234   bnk = (TAO_BNK *)tao->data;
235   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
236   PetscFunctionReturn(0);
237 }