xref: /petsc/src/tao/bound/impls/bnk/bntl.c (revision f5766c097db85df4ee62e2acb750cb76ef9efd11)
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) {
127c4b75bccSAlp Dener     ++tao->niter;
12862675beeSAlp Dener 
12989da521bSAlp Dener     if (needH && bnk->inactive_idx) {
130e031d6f5SAlp Dener       /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
131e031d6f5SAlp Dener       ierr = TaoBNKTakeCGSteps(tao, &cgTerminate);CHKERRQ(ierr);
132e031d6f5SAlp Dener       if (cgTerminate) {
133e031d6f5SAlp Dener         tao->reason = bnk->bncg->reason;
134e031d6f5SAlp Dener         PetscFunctionReturn(0);
135e031d6f5SAlp Dener       }
13608752603SAlp Dener       /* Compute the hessian and update the BFGS preconditioner at the new iterate */
137e0ed867bSAlp Dener       ierr = bnk->computehessian(tao);CHKERRQ(ierr);
138937a31a1SAlp Dener       needH = PETSC_FALSE;
139937a31a1SAlp Dener     }
140c14b763aSAlp Dener 
1418d5ead36SAlp Dener     /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
1426b591159SAlp Dener     ierr = bnk->computestep(tao, shift, &ksp_reason, &stepType);CHKERRQ(ierr);
143c14b763aSAlp Dener 
144c14b763aSAlp Dener     /* Store current solution before it changes */
145c14b763aSAlp Dener     oldTrust = tao->trust;
146c14b763aSAlp Dener     bnk->fold = bnk->f;
147c14b763aSAlp Dener     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
148c14b763aSAlp Dener     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
149c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
150c14b763aSAlp Dener 
151c14b763aSAlp Dener     /* Temporarily accept the step and project it into the bounds */
152c14b763aSAlp Dener     ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
1533b063059SAlp Dener     ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
154c14b763aSAlp Dener 
155c14b763aSAlp Dener     /* Check if the projection changed the step direction */
156c4b75bccSAlp Dener     if (nDiff > 0) {
157c4b75bccSAlp Dener       /* Projection changed the step, so we have to recompute the step and
158c4b75bccSAlp Dener          the predicted reduction. Leave the trust radius unchanged. */
159c14b763aSAlp Dener       ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
1608d5ead36SAlp Dener       ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
1615e9b73cbSAlp Dener       ierr = TaoBNKRecomputePred(tao, tao->stepdirection, &prered);CHKERRQ(ierr);
162c14b763aSAlp Dener     } else {
163c14b763aSAlp Dener       /* Step did not change, so we can just recover the pre-computed prediction */
164c14b763aSAlp Dener       ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
165c14b763aSAlp Dener     }
166c14b763aSAlp Dener     prered = -prered;
167c14b763aSAlp Dener 
168c14b763aSAlp Dener     /* Compute the actual reduction and update the trust radius */
169c14b763aSAlp Dener     ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
170b4a30f08SAlp Dener     if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
171c14b763aSAlp Dener     actred = bnk->fold - bnk->f;
17228017e9fSAlp Dener     ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr);
173c14b763aSAlp Dener 
174c14b763aSAlp Dener     if (stepAccepted) {
175c14b763aSAlp Dener       /* Step is good, evaluate the gradient and the hessian */
1768d5ead36SAlp Dener       steplen = 1.0;
177937a31a1SAlp Dener       needH = PETSC_TRUE;
178e465cd6fSAlp Dener       ++bnk->newt;
179c14b763aSAlp Dener       ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
18061be54a6SAlp Dener       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
18161be54a6SAlp Dener       ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
18261be54a6SAlp Dener       ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
183*f5766c09SAlp Dener       ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
184c14b763aSAlp Dener     } else {
185c14b763aSAlp Dener       /* Trust-region rejected the step. Revert the solution. */
186c14b763aSAlp Dener       bnk->f = bnk->fold;
187c14b763aSAlp Dener       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
188c14b763aSAlp Dener       /* Trigger the line search */
189e465cd6fSAlp Dener       ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr);
190937a31a1SAlp Dener       ierr = TaoBNKPerformLineSearch(tao, &stepType, &steplen, &ls_reason);CHKERRQ(ierr);
191c14b763aSAlp Dener       if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
192c14b763aSAlp Dener         /* Line search failed, revert solution and terminate */
193c0f10754SAlp Dener         stepAccepted = PETSC_FALSE;
194937a31a1SAlp Dener         needH = PETSC_FALSE;
195c14b763aSAlp Dener         bnk->f = bnk->fold;
196c14b763aSAlp Dener         ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
197c14b763aSAlp Dener         ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
198c14b763aSAlp Dener         ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
199c14b763aSAlp Dener         tao->trust = 0.0;
200c14b763aSAlp Dener         tao->reason = TAO_DIVERGED_LS_FAILURE;
201c14b763aSAlp Dener       } else {
202937a31a1SAlp Dener         /* new iterate so we need to recompute the Hessian */
203937a31a1SAlp Dener         needH = PETSC_TRUE;
204198282dbSAlp Dener         /* compute the projected gradient */
20561be54a6SAlp Dener         ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
20661be54a6SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
20761be54a6SAlp Dener         ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
208*f5766c09SAlp Dener         ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
209c14b763aSAlp Dener         /* Line search succeeded so we should update the trust radius based on the LS step length */
210770b7498SAlp Dener         tao->trust = oldTrust;
21128017e9fSAlp Dener         ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, BNK_UPDATE_STEP, stepType, &stepAccepted);CHKERRQ(ierr);
21262675beeSAlp Dener         /* count the accepted step type */
21362675beeSAlp Dener         ierr = TaoBNKAddStepCounts(tao, stepType);CHKERRQ(ierr);
214c14b763aSAlp Dener       }
215c14b763aSAlp Dener     }
216c14b763aSAlp Dener 
217c14b763aSAlp Dener     /*  Check for termination */
2180b7db9bbSAlp Dener     ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
2190b7db9bbSAlp Dener     ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
220b4a30f08SAlp Dener     if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
2219b6ef848SAlp Dener     ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
2229b6ef848SAlp Dener     ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr);
223c14b763aSAlp Dener     ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr);
224c14b763aSAlp Dener   }
225c14b763aSAlp Dener   PetscFunctionReturn(0);
226c14b763aSAlp Dener }
227c14b763aSAlp Dener 
228df278d8fSAlp Dener /*------------------------------------------------------------*/
229df278d8fSAlp Dener 
230e0ed867bSAlp Dener static PetscErrorCode TaoSetFromOptions_BNTL(PetscOptionItems *PetscOptionsObject,Tao tao)
2319b6ef848SAlp Dener {
2329b6ef848SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
2339b6ef848SAlp Dener   PetscErrorCode ierr;
2349b6ef848SAlp Dener 
2359b6ef848SAlp Dener   PetscFunctionBegin;
236e0ed867bSAlp Dener   ierr = TaoSetFromOptions_BNK(PetscOptionsObject, tao);CHKERRQ(ierr);
237e0ed867bSAlp Dener   if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
2389b6ef848SAlp 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)");
2399b6ef848SAlp Dener   PetscFunctionReturn(0);
2409b6ef848SAlp Dener }
2419b6ef848SAlp Dener 
2429b6ef848SAlp Dener /*------------------------------------------------------------*/
2439b6ef848SAlp Dener 
244e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNTL(Tao tao)
245c14b763aSAlp Dener {
246c14b763aSAlp Dener   TAO_BNK        *bnk;
247c14b763aSAlp Dener   PetscErrorCode ierr;
248c14b763aSAlp Dener 
249c14b763aSAlp Dener   PetscFunctionBegin;
250c14b763aSAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
251c14b763aSAlp Dener   tao->ops->solve=TaoSolve_BNTL;
252e0ed867bSAlp Dener   tao->ops->setfromoptions=TaoSetFromOptions_BNTL;
253c14b763aSAlp Dener 
254c14b763aSAlp Dener   bnk = (TAO_BNK *)tao->data;
255c14b763aSAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
256c14b763aSAlp Dener   PetscFunctionReturn(0);
257c14b763aSAlp Dener }
258