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