xref: /petsc/src/tao/bound/impls/bnk/bntl.c (revision 9b6ef8482b0295dfac0fe7cdd25b1a30d1be2d60)
1c14b763aSAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
2c14b763aSAlp Dener #include <petscksp.h>
3c14b763aSAlp Dener 
4c14b763aSAlp Dener /*
5c14b763aSAlp Dener  Implements Newton's Method with a trust region approach for solving
6198282dbSAlp Dener  bound constrained minimization problems.
7c14b763aSAlp Dener 
8198282dbSAlp Dener  ------------------------------------------------------------
9198282dbSAlp Dener 
10198282dbSAlp Dener  initialize trust radius (default: BNK_INIT_INTERPOLATION)
11198282dbSAlp Dener  x_0 = VecMedian(x_0)
12198282dbSAlp Dener  f_0, g_0 = TaoComputeObjectiveAndGradient(x_0)
13198282dbSAlp Dener  pg_0 = VecBoundGradientProjection(g_0)
14198282dbSAlp Dener  check convergence at pg_0
15198282dbSAlp Dener  niter = 0
16198282dbSAlp Dener  step_accepted = true
17198282dbSAlp Dener 
18198282dbSAlp Dener  while niter <= max_it
19198282dbSAlp Dener     niter += 1
20198282dbSAlp Dener     H_k = TaoComputeHessian(x_k)
21198282dbSAlp Dener     if pc_type == BNK_PC_BFGS
22198282dbSAlp Dener       add correction to BFGS approx
23198282dbSAlp Dener       if scale_type == BNK_SCALE_AHESS
24198282dbSAlp Dener         D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
25198282dbSAlp Dener         scale BFGS with VecReciprocal(D)
26198282dbSAlp Dener       end
27198282dbSAlp Dener     end
28198282dbSAlp Dener 
29198282dbSAlp Dener     if pc_type = BNK_PC_BFGS
30198282dbSAlp Dener       B_k = BFGS
31198282dbSAlp Dener     else
32198282dbSAlp Dener       B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
33198282dbSAlp Dener       B_k = VecReciprocal(B_k)
34198282dbSAlp Dener     end
35198282dbSAlp Dener     w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
36198282dbSAlp Dener     eps = min(eps, norm2(w))
37198282dbSAlp Dener     determine the active and inactive index sets such that
38198282dbSAlp Dener       L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
39198282dbSAlp Dener       U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
40198282dbSAlp Dener       F = {i : l_i = (x_k)_i = u_i}
41198282dbSAlp Dener       A = {L + U + F}
42198282dbSAlp Dener       I = {i : i not in A}
43198282dbSAlp Dener 
44198282dbSAlp Dener     generate the reduced system Hr_k dr_k = -gr_k for variables in I
45198282dbSAlp Dener     if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
46198282dbSAlp Dener       D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
47198282dbSAlp Dener       scale BFGS with VecReciprocal(D)
48198282dbSAlp Dener     end
49198282dbSAlp Dener     solve Hr_k dr_k = -gr_k
50198282dbSAlp Dener     set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
51198282dbSAlp Dener 
52198282dbSAlp Dener     x_{k+1} = VecMedian(x_k + d_k)
53198282dbSAlp Dener     s = x_{k+1} - x_k
54198282dbSAlp Dener     prered = dot(s, 0.5*gr_k - Hr_k*s)
55198282dbSAlp Dener     f_{k+1} = TaoComputeObjective(x_{k+1})
56198282dbSAlp Dener     actred = f_k - f_{k+1}
57198282dbSAlp Dener 
58198282dbSAlp Dener     oldTrust = trust
59198282dbSAlp Dener     step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION)
60198282dbSAlp Dener     if step_accepted
61198282dbSAlp Dener       g_{k+1} = TaoComputeGradient(x_{k+1})
62198282dbSAlp Dener       pg_{k+1} = VecBoundGradientProjection(g_{k+1})
63198282dbSAlp Dener       count the accepted Newton step
64198282dbSAlp Dener     else
65198282dbSAlp Dener       if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
66198282dbSAlp Dener         dr_k = -BFGS*gr_k for variables in I
67198282dbSAlp Dener         if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
68198282dbSAlp Dener           reset the BFGS preconditioner
69198282dbSAlp Dener           calculate scale delta and apply it to BFGS
70198282dbSAlp Dener           dr_k = -BFGS*gr_k for variables in I
71198282dbSAlp Dener           if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
72198282dbSAlp Dener             dr_k = -gr_k for variables in I
73198282dbSAlp Dener           end
74198282dbSAlp Dener         end
75198282dbSAlp Dener       end
76198282dbSAlp Dener 
77198282dbSAlp Dener       x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch()
78198282dbSAlp Dener       if ls_failed
79198282dbSAlp Dener         f_{k+1} = f_k
80198282dbSAlp Dener         x_{k+1} = x_k
81198282dbSAlp Dener         g_{k+1} = g_k
82198282dbSAlp Dener         pg_{k+1} = pg_k
83198282dbSAlp Dener         terminate
84198282dbSAlp Dener       else
85198282dbSAlp Dener         pg_{k+1} = VecBoundGradientProjection(g_{k+1})
86198282dbSAlp Dener         trust = oldTrust
87198282dbSAlp Dener         trust = TaoBNKUpdateTrustRadius(BNK_UPDATE_STEP)
88198282dbSAlp Dener         count the accepted step type (Newton, BFGS, scaled grad or grad)
89198282dbSAlp Dener       end
90198282dbSAlp Dener     end
91198282dbSAlp Dener 
92198282dbSAlp Dener     check convergence at pg_{k+1}
93198282dbSAlp Dener  end
94c14b763aSAlp Dener */
95c14b763aSAlp Dener 
96c14b763aSAlp Dener static PetscErrorCode TaoSolve_BNTL(Tao tao)
97c14b763aSAlp Dener {
98c14b763aSAlp Dener   PetscErrorCode               ierr;
99c14b763aSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
100e465cd6fSAlp Dener   KSPConvergedReason           ksp_reason;
101c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
102c14b763aSAlp Dener 
103*9b6ef848SAlp Dener   PetscReal                    resnorm, oldTrust, prered, actred, stepNorm, steplen;
10462675beeSAlp Dener   PetscBool                    stepAccepted = PETSC_TRUE, shift = PETSC_FALSE;
10528017e9fSAlp Dener   PetscInt                     stepType;
106c14b763aSAlp Dener 
107c14b763aSAlp Dener   PetscFunctionBegin;
10828017e9fSAlp Dener   /* Initialize the preconditioner, KSP solver and trust radius/line search */
109c14b763aSAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
11062675beeSAlp Dener   ierr = TaoBNKInitialize(tao, bnk->init_type);CHKERRQ(ierr);
11128017e9fSAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
112c14b763aSAlp Dener 
113c14b763aSAlp Dener   /* Have not converged; continue with Newton method */
114c14b763aSAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
115c14b763aSAlp Dener     tao->niter++;
116c14b763aSAlp Dener     tao->ksp_its=0;
11762675beeSAlp Dener 
11862675beeSAlp Dener     /* Compute the hessian and update the BFGS preconditioner at the new iterate*/
11962675beeSAlp Dener     ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
120c14b763aSAlp Dener 
1218d5ead36SAlp Dener     /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
12262675beeSAlp Dener     ierr = TaoBNKComputeStep(tao, shift, &ksp_reason);CHKERRQ(ierr);
123c14b763aSAlp Dener 
124c14b763aSAlp Dener     /* Store current solution before it changes */
125c14b763aSAlp Dener     oldTrust = tao->trust;
126c14b763aSAlp Dener     bnk->fold = bnk->f;
127c14b763aSAlp Dener     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
128c14b763aSAlp Dener     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
129c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
130c14b763aSAlp Dener 
131c14b763aSAlp Dener     /* Temporarily accept the step and project it into the bounds */
132c14b763aSAlp Dener     ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
133c14b763aSAlp Dener     ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
134c14b763aSAlp Dener 
135c14b763aSAlp Dener     /* Check if the projection changed the step direction */
136c14b763aSAlp Dener     ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
1378d5ead36SAlp Dener     ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
138c14b763aSAlp Dener     ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr);
139c14b763aSAlp Dener     if (stepNorm != bnk->dnorm) {
1408d5ead36SAlp Dener       /* Projection changed the step, so we have to recompute predicted reduction.
1418d5ead36SAlp Dener          However, we deliberately do not change the step norm and the trust radius
1428d5ead36SAlp Dener          in order for the safeguard to more closely mimic a piece-wise linesearch
1438d5ead36SAlp Dener          along the bounds. */
14428017e9fSAlp Dener       ierr = MatMult(bnk->H_inactive, tao->stepdirection, bnk->Xwork);CHKERRQ(ierr);
145198282dbSAlp Dener       ierr = VecAYPX(bnk->Xwork, -0.5, bnk->G_inactive);CHKERRQ(ierr);
146c14b763aSAlp Dener       ierr = VecDot(bnk->Xwork, tao->stepdirection, &prered);
147c14b763aSAlp Dener     } else {
148c14b763aSAlp Dener       /* Step did not change, so we can just recover the pre-computed prediction */
149c14b763aSAlp Dener       ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
150c14b763aSAlp Dener     }
151c14b763aSAlp Dener     prered = -prered;
152c14b763aSAlp Dener 
153c14b763aSAlp Dener     /* Compute the actual reduction and update the trust radius */
154c14b763aSAlp Dener     ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
155c14b763aSAlp Dener     actred = bnk->fold - bnk->f;
15628017e9fSAlp Dener     ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr);
157c14b763aSAlp Dener 
158c14b763aSAlp Dener     if (stepAccepted) {
159c14b763aSAlp Dener       /* Step is good, evaluate the gradient and the hessian */
1608d5ead36SAlp Dener       steplen = 1.0;
161e465cd6fSAlp Dener       ++bnk->newt;
162c14b763aSAlp Dener       ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
163c14b763aSAlp Dener       ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
164c14b763aSAlp Dener     } else {
165c14b763aSAlp Dener       /* Trust-region rejected the step. Revert the solution. */
166c14b763aSAlp Dener       bnk->f = bnk->fold;
167c14b763aSAlp Dener       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
168c14b763aSAlp Dener 
169c14b763aSAlp Dener       /* Trigger the line search */
170e465cd6fSAlp Dener       ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr);
171c14b763aSAlp Dener       ierr = TaoBNKPerformLineSearch(tao, stepType, &steplen, &ls_reason);CHKERRQ(ierr);
172c14b763aSAlp Dener       if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
173c14b763aSAlp Dener         /* Line search failed, revert solution and terminate */
174c14b763aSAlp Dener         bnk->f = bnk->fold;
175c14b763aSAlp Dener         ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
176c14b763aSAlp Dener         ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
177c14b763aSAlp Dener         ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
178c14b763aSAlp Dener         tao->trust = 0.0;
179c14b763aSAlp Dener         tao->reason = TAO_DIVERGED_LS_FAILURE;
180c14b763aSAlp Dener       } else {
181198282dbSAlp Dener         /* compute the projected gradient */
182198282dbSAlp Dener         ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
183*9b6ef848SAlp Dener         ierr = VecNorm(tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
184*9b6ef848SAlp Dener         if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
185c14b763aSAlp Dener         /* Line search succeeded so we should update the trust radius based on the LS step length */
186770b7498SAlp Dener         tao->trust = oldTrust;
18728017e9fSAlp Dener         ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, BNK_UPDATE_STEP, stepType, &stepAccepted);CHKERRQ(ierr);
18862675beeSAlp Dener         /* count the accepted step type */
18962675beeSAlp Dener         ierr = TaoBNKAddStepCounts(tao, stepType);CHKERRQ(ierr);
190c14b763aSAlp Dener       }
191c14b763aSAlp Dener     }
192c14b763aSAlp Dener 
193c14b763aSAlp Dener     /*  Check for termination */
194*9b6ef848SAlp Dener     ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr);
195*9b6ef848SAlp Dener     ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr);
196*9b6ef848SAlp Dener     ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
197*9b6ef848SAlp Dener     ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr);
198c14b763aSAlp Dener     ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr);
199c14b763aSAlp Dener   }
200c14b763aSAlp Dener   PetscFunctionReturn(0);
201c14b763aSAlp Dener }
202c14b763aSAlp Dener 
203df278d8fSAlp Dener /*------------------------------------------------------------*/
204df278d8fSAlp Dener 
205*9b6ef848SAlp Dener PETSC_INTERN PetscErrorCode TaoSetUp_BNTL(Tao tao)
206*9b6ef848SAlp Dener {
207*9b6ef848SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
208*9b6ef848SAlp Dener   PetscErrorCode ierr;
209*9b6ef848SAlp Dener 
210*9b6ef848SAlp Dener   PetscFunctionBegin;
211*9b6ef848SAlp Dener   ierr = TaoSetUp_BNK(tao);CHKERRQ(ierr);
212*9b6ef848SAlp Dener   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)");
213*9b6ef848SAlp Dener   PetscFunctionReturn(0);
214*9b6ef848SAlp Dener }
215*9b6ef848SAlp Dener 
216*9b6ef848SAlp Dener /*------------------------------------------------------------*/
217*9b6ef848SAlp Dener 
218*9b6ef848SAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BNTL(Tao tao)
219c14b763aSAlp Dener {
220c14b763aSAlp Dener   TAO_BNK        *bnk;
221c14b763aSAlp Dener   PetscErrorCode ierr;
222c14b763aSAlp Dener 
223c14b763aSAlp Dener   PetscFunctionBegin;
224c14b763aSAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
225c14b763aSAlp Dener   tao->ops->solve=TaoSolve_BNTL;
226*9b6ef848SAlp Dener   tao->ops->setup=TaoSetUp_BNTL;
227c14b763aSAlp Dener 
228c14b763aSAlp Dener   bnk = (TAO_BNK *)tao->data;
229c14b763aSAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
230c14b763aSAlp Dener   bnk->sval = 0.0; /* disable Hessian shifting */
231c14b763aSAlp Dener   PetscFunctionReturn(0);
232c14b763aSAlp Dener }