xref: /petsc/src/tao/bound/impls/bnk/bntl.c (revision 8fcddce65efd55a8fe3f87d4c08c15577ce4cbef)
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 
8c4b75bccSAlp Dener  In this variant, the trust region failures trigger a line search with
9c4b75bccSAlp Dener  the existing Newton step instead of re-solving the step with a
10c4b75bccSAlp Dener  different radius.
11c4b75bccSAlp Dener 
12198282dbSAlp Dener  ------------------------------------------------------------
13198282dbSAlp Dener 
14198282dbSAlp Dener  x_0 = VecMedian(x_0)
15198282dbSAlp Dener  f_0, g_0 = TaoComputeObjectiveAndGradient(x_0)
16c4b75bccSAlp Dener  pg_0 = project(g_0)
17198282dbSAlp Dener  check convergence at pg_0
18c4b75bccSAlp Dener  needH = TaoBNKInitialize(default:BNK_INIT_INTERPOLATION)
19198282dbSAlp Dener  niter = 0
20198282dbSAlp Dener  step_accepted = true
21198282dbSAlp Dener 
22198282dbSAlp Dener  while niter <= max_it
23198282dbSAlp Dener     niter += 1
24c4b75bccSAlp Dener 
25c4b75bccSAlp Dener     if needH
26c4b75bccSAlp Dener       If max_cg_steps > 0
27c4b75bccSAlp Dener         x_k, g_k, pg_k = TaoSolve(BNCG)
28c4b75bccSAlp Dener       end
29c4b75bccSAlp Dener 
30198282dbSAlp Dener       H_k = TaoComputeHessian(x_k)
31198282dbSAlp Dener       if pc_type == BNK_PC_BFGS
32198282dbSAlp Dener         add correction to BFGS approx
33198282dbSAlp Dener         if scale_type == BNK_SCALE_AHESS
34198282dbSAlp Dener           D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
35198282dbSAlp Dener           scale BFGS with VecReciprocal(D)
36198282dbSAlp Dener         end
37198282dbSAlp Dener       end
38c4b75bccSAlp Dener       needH = False
39c4b75bccSAlp Dener     end
40198282dbSAlp Dener 
41198282dbSAlp Dener     if pc_type = BNK_PC_BFGS
42198282dbSAlp Dener       B_k = BFGS
43198282dbSAlp Dener     else
44198282dbSAlp Dener       B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
45198282dbSAlp Dener       B_k = VecReciprocal(B_k)
46198282dbSAlp Dener     end
47198282dbSAlp Dener     w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
48198282dbSAlp Dener     eps = min(eps, norm2(w))
49198282dbSAlp Dener     determine the active and inactive index sets such that
50198282dbSAlp Dener       L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
51198282dbSAlp Dener       U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
52198282dbSAlp Dener       F = {i : l_i = (x_k)_i = u_i}
53198282dbSAlp Dener       A = {L + U + F}
54c4b75bccSAlp Dener       IA = {i : i not in A}
55198282dbSAlp Dener 
56c4b75bccSAlp Dener     generate the reduced system Hr_k dr_k = -gr_k for variables in IA
57198282dbSAlp Dener     if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
58198282dbSAlp Dener       D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
59198282dbSAlp Dener       scale BFGS with VecReciprocal(D)
60198282dbSAlp Dener     end
61198282dbSAlp Dener     solve Hr_k dr_k = -gr_k
62198282dbSAlp Dener     set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
63198282dbSAlp Dener 
64198282dbSAlp Dener     x_{k+1} = VecMedian(x_k + d_k)
65198282dbSAlp Dener     s = x_{k+1} - x_k
66198282dbSAlp Dener     prered = dot(s, 0.5*gr_k - Hr_k*s)
67198282dbSAlp Dener     f_{k+1} = TaoComputeObjective(x_{k+1})
68198282dbSAlp Dener     actred = f_k - f_{k+1}
69198282dbSAlp Dener 
70198282dbSAlp Dener     oldTrust = trust
71198282dbSAlp Dener     step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION)
72198282dbSAlp Dener     if step_accepted
73198282dbSAlp Dener       g_{k+1} = TaoComputeGradient(x_{k+1})
74c4b75bccSAlp Dener       pg_{k+1} = project(g_{k+1})
75198282dbSAlp Dener       count the accepted Newton step
76198282dbSAlp Dener     else
77198282dbSAlp Dener       if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
78198282dbSAlp Dener         dr_k = -BFGS*gr_k for variables in I
79198282dbSAlp Dener         if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
80198282dbSAlp Dener           reset the BFGS preconditioner
81198282dbSAlp Dener           calculate scale delta and apply it to BFGS
82198282dbSAlp Dener           dr_k = -BFGS*gr_k for variables in I
83198282dbSAlp Dener           if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
84198282dbSAlp Dener             dr_k = -gr_k for variables in I
85198282dbSAlp Dener           end
86198282dbSAlp Dener         end
87198282dbSAlp Dener       end
88198282dbSAlp Dener 
89198282dbSAlp Dener       x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch()
90198282dbSAlp Dener       if ls_failed
91198282dbSAlp Dener         f_{k+1} = f_k
92198282dbSAlp Dener         x_{k+1} = x_k
93198282dbSAlp Dener         g_{k+1} = g_k
94198282dbSAlp Dener         pg_{k+1} = pg_k
95198282dbSAlp Dener         terminate
96198282dbSAlp Dener       else
97c4b75bccSAlp Dener         pg_{k+1} = project(g_{k+1})
98198282dbSAlp Dener         trust = oldTrust
99198282dbSAlp Dener         trust = TaoBNKUpdateTrustRadius(BNK_UPDATE_STEP)
100198282dbSAlp Dener         count the accepted step type (Newton, BFGS, scaled grad or grad)
101198282dbSAlp Dener       end
102198282dbSAlp Dener     end
103198282dbSAlp Dener 
104198282dbSAlp Dener     check convergence at pg_{k+1}
105198282dbSAlp Dener  end
106c14b763aSAlp Dener */
107c14b763aSAlp Dener 
108e0ed867bSAlp Dener PetscErrorCode TaoSolve_BNTL(Tao tao)
109c14b763aSAlp Dener {
110c14b763aSAlp Dener   PetscErrorCode               ierr;
111c14b763aSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
112e465cd6fSAlp Dener   KSPConvergedReason           ksp_reason;
113c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
114c14b763aSAlp Dener 
11589da521bSAlp Dener   PetscReal                    oldTrust, prered, actred, steplen, resnorm;
116937a31a1SAlp Dener   PetscBool                    cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE;
117c4b75bccSAlp Dener   PetscInt                     stepType, nDiff;
118c14b763aSAlp Dener 
119c14b763aSAlp Dener   PetscFunctionBegin;
12028017e9fSAlp Dener   /* Initialize the preconditioner, KSP solver and trust radius/line search */
121c14b763aSAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
122937a31a1SAlp Dener   ierr = TaoBNKInitialize(tao, bnk->init_type, &needH);CHKERRQ(ierr);
12328017e9fSAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
124c14b763aSAlp Dener 
125c14b763aSAlp Dener   /* Have not converged; continue with Newton method */
126c14b763aSAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
127e1e80dc8SAlp Dener     /* Call general purpose update function */
128e1e80dc8SAlp Dener     if (tao->ops->update) {
129*8fcddce6SStefano Zampini       ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr);
130e1e80dc8SAlp Dener     }
131c4b75bccSAlp Dener     ++tao->niter;
13262675beeSAlp Dener 
13389da521bSAlp Dener     if (needH && bnk->inactive_idx) {
134e031d6f5SAlp Dener       /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
135e031d6f5SAlp Dener       ierr = TaoBNKTakeCGSteps(tao, &cgTerminate);CHKERRQ(ierr);
136e031d6f5SAlp Dener       if (cgTerminate) {
137e031d6f5SAlp Dener         tao->reason = bnk->bncg->reason;
138e031d6f5SAlp Dener         PetscFunctionReturn(0);
139e031d6f5SAlp Dener       }
14008752603SAlp Dener       /* Compute the hessian and update the BFGS preconditioner at the new iterate */
141f7bf01afSAlp Dener       ierr = (*bnk->computehessian)(tao);CHKERRQ(ierr);
142937a31a1SAlp Dener       needH = PETSC_FALSE;
143937a31a1SAlp Dener     }
144c14b763aSAlp Dener 
1458d5ead36SAlp Dener     /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
146f7bf01afSAlp Dener     ierr = (*bnk->computestep)(tao, shift, &ksp_reason, &stepType);CHKERRQ(ierr);
147c14b763aSAlp Dener 
148c14b763aSAlp Dener     /* Store current solution before it changes */
149c14b763aSAlp Dener     oldTrust = tao->trust;
150c14b763aSAlp Dener     bnk->fold = bnk->f;
151c14b763aSAlp Dener     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
152c14b763aSAlp Dener     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
153c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
154c14b763aSAlp Dener 
155c14b763aSAlp Dener     /* Temporarily accept the step and project it into the bounds */
156c14b763aSAlp Dener     ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
1573b063059SAlp Dener     ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
158c14b763aSAlp Dener 
159c14b763aSAlp Dener     /* Check if the projection changed the step direction */
160c4b75bccSAlp Dener     if (nDiff > 0) {
161c4b75bccSAlp Dener       /* Projection changed the step, so we have to recompute the step and
162c4b75bccSAlp Dener          the predicted reduction. Leave the trust radius unchanged. */
163c14b763aSAlp Dener       ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
1648d5ead36SAlp Dener       ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
1655e9b73cbSAlp Dener       ierr = TaoBNKRecomputePred(tao, tao->stepdirection, &prered);CHKERRQ(ierr);
166c14b763aSAlp Dener     } else {
167c14b763aSAlp Dener       /* Step did not change, so we can just recover the pre-computed prediction */
168c14b763aSAlp Dener       ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
169c14b763aSAlp Dener     }
170c14b763aSAlp Dener     prered = -prered;
171c14b763aSAlp Dener 
172c14b763aSAlp Dener     /* Compute the actual reduction and update the trust radius */
173c14b763aSAlp Dener     ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
174b4a30f08SAlp Dener     if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
175c14b763aSAlp Dener     actred = bnk->fold - bnk->f;
17628017e9fSAlp Dener     ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr);
177c14b763aSAlp Dener 
178c14b763aSAlp Dener     if (stepAccepted) {
179c14b763aSAlp Dener       /* Step is good, evaluate the gradient and the hessian */
1808d5ead36SAlp Dener       steplen = 1.0;
181937a31a1SAlp Dener       needH = PETSC_TRUE;
182e465cd6fSAlp Dener       ++bnk->newt;
183c14b763aSAlp Dener       ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
18461be54a6SAlp Dener       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
18561be54a6SAlp Dener       ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
18661be54a6SAlp Dener       ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
187f5766c09SAlp Dener       ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
188c14b763aSAlp Dener     } else {
189c14b763aSAlp Dener       /* Trust-region rejected the step. Revert the solution. */
190c14b763aSAlp Dener       bnk->f = bnk->fold;
191c14b763aSAlp Dener       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
192c14b763aSAlp Dener       /* Trigger the line search */
193e465cd6fSAlp Dener       ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr);
194937a31a1SAlp Dener       ierr = TaoBNKPerformLineSearch(tao, &stepType, &steplen, &ls_reason);CHKERRQ(ierr);
195c14b763aSAlp Dener       if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
196c14b763aSAlp Dener         /* Line search failed, revert solution and terminate */
197c0f10754SAlp Dener         stepAccepted = PETSC_FALSE;
198937a31a1SAlp Dener         needH = PETSC_FALSE;
199c14b763aSAlp Dener         bnk->f = bnk->fold;
200c14b763aSAlp Dener         ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
201c14b763aSAlp Dener         ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
202c14b763aSAlp Dener         ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
203c14b763aSAlp Dener         tao->trust = 0.0;
204c14b763aSAlp Dener         tao->reason = TAO_DIVERGED_LS_FAILURE;
205c14b763aSAlp Dener       } else {
206937a31a1SAlp Dener         /* new iterate so we need to recompute the Hessian */
207937a31a1SAlp Dener         needH = PETSC_TRUE;
208198282dbSAlp Dener         /* compute the projected gradient */
20961be54a6SAlp Dener         ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
21061be54a6SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
21161be54a6SAlp Dener         ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
212f5766c09SAlp Dener         ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
213c14b763aSAlp Dener         /* Line search succeeded so we should update the trust radius based on the LS step length */
214770b7498SAlp Dener         tao->trust = oldTrust;
21528017e9fSAlp Dener         ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, BNK_UPDATE_STEP, stepType, &stepAccepted);CHKERRQ(ierr);
21662675beeSAlp Dener         /* count the accepted step type */
21762675beeSAlp Dener         ierr = TaoBNKAddStepCounts(tao, stepType);CHKERRQ(ierr);
218c14b763aSAlp Dener       }
219c14b763aSAlp Dener     }
220c14b763aSAlp Dener 
221c14b763aSAlp Dener     /*  Check for termination */
2220b7db9bbSAlp Dener     ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
2230b7db9bbSAlp Dener     ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
224b4a30f08SAlp Dener     if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
2259b6ef848SAlp Dener     ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
2269b6ef848SAlp Dener     ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr);
227c14b763aSAlp Dener     ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr);
228c14b763aSAlp Dener   }
229c14b763aSAlp Dener   PetscFunctionReturn(0);
230c14b763aSAlp Dener }
231c14b763aSAlp Dener 
232df278d8fSAlp Dener /*------------------------------------------------------------*/
233e0ed867bSAlp Dener static PetscErrorCode TaoSetFromOptions_BNTL(PetscOptionItems *PetscOptionsObject,Tao tao)
2349b6ef848SAlp Dener {
2359b6ef848SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
2369b6ef848SAlp Dener   PetscErrorCode ierr;
2379b6ef848SAlp Dener 
2389b6ef848SAlp Dener   PetscFunctionBegin;
239e0ed867bSAlp Dener   ierr = TaoSetFromOptions_BNK(PetscOptionsObject, tao);CHKERRQ(ierr);
240e0ed867bSAlp Dener   if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
2419b6ef848SAlp 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)");
2429b6ef848SAlp Dener   PetscFunctionReturn(0);
2439b6ef848SAlp Dener }
2449b6ef848SAlp Dener 
2459b6ef848SAlp Dener /*------------------------------------------------------------*/
2463850be85SAlp Dener /*MC
2473850be85SAlp Dener   TAOBNTL - Bounded Newton Trust Region method with line-search fall-back for nonlinear
2483850be85SAlp Dener             minimization with bound constraints.
2499b6ef848SAlp Dener 
2503850be85SAlp Dener   Options Database Keys:
2513850be85SAlp Dener   + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
2523850be85SAlp Dener   . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
2533850be85SAlp Dener   . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
2543850be85SAlp Dener   - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
2553850be85SAlp Dener 
2563850be85SAlp Dener   Level: beginner
2573850be85SAlp Dener M*/
258e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNTL(Tao tao)
259c14b763aSAlp Dener {
260c14b763aSAlp Dener   TAO_BNK        *bnk;
261c14b763aSAlp Dener   PetscErrorCode ierr;
262c14b763aSAlp Dener 
263c14b763aSAlp Dener   PetscFunctionBegin;
264c14b763aSAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
265c14b763aSAlp Dener   tao->ops->solve=TaoSolve_BNTL;
266e0ed867bSAlp Dener   tao->ops->setfromoptions=TaoSetFromOptions_BNTL;
267c14b763aSAlp Dener 
268c14b763aSAlp Dener   bnk = (TAO_BNK *)tao->data;
269c14b763aSAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
270c14b763aSAlp Dener   PetscFunctionReturn(0);
271c14b763aSAlp Dener }
272