xref: /petsc/src/tao/bound/impls/bnk/bntl.c (revision 03047865b8d8757cf1cf9cda45785c1537b01dc1)
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  x_0 = VecMedian(x_0)
13198282dbSAlp Dener  f_0, g_0 = TaoComputeObjectiveAndGradient(x_0)
14c4b75bccSAlp Dener  pg_0 = project(g_0)
15198282dbSAlp Dener  check convergence at pg_0
16c4b75bccSAlp Dener  needH = TaoBNKInitialize(default:BNK_INIT_INTERPOLATION)
17198282dbSAlp Dener  niter = 0
18198282dbSAlp Dener  step_accepted = true
19198282dbSAlp Dener 
20198282dbSAlp Dener  while niter <= max_it
21198282dbSAlp Dener     niter += 1
22c4b75bccSAlp Dener 
23c4b75bccSAlp Dener     if needH
24c4b75bccSAlp Dener       If max_cg_steps > 0
25c4b75bccSAlp Dener         x_k, g_k, pg_k = TaoSolve(BNCG)
26c4b75bccSAlp Dener       end
27c4b75bccSAlp Dener 
28198282dbSAlp Dener       H_k = TaoComputeHessian(x_k)
29198282dbSAlp Dener       if pc_type == BNK_PC_BFGS
30198282dbSAlp Dener         add correction to BFGS approx
31198282dbSAlp Dener         if scale_type == BNK_SCALE_AHESS
32198282dbSAlp Dener           D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
33198282dbSAlp Dener           scale BFGS with VecReciprocal(D)
34198282dbSAlp Dener         end
35198282dbSAlp Dener       end
36c4b75bccSAlp Dener       needH = False
37c4b75bccSAlp Dener     end
38198282dbSAlp Dener 
39198282dbSAlp Dener     if pc_type = BNK_PC_BFGS
40198282dbSAlp Dener       B_k = BFGS
41198282dbSAlp Dener     else
42198282dbSAlp Dener       B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
43198282dbSAlp Dener       B_k = VecReciprocal(B_k)
44198282dbSAlp Dener     end
45198282dbSAlp Dener     w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
46198282dbSAlp Dener     eps = min(eps, norm2(w))
47198282dbSAlp Dener     determine the active and inactive index sets such that
48198282dbSAlp Dener       L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
49198282dbSAlp Dener       U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
50198282dbSAlp Dener       F = {i : l_i = (x_k)_i = u_i}
51198282dbSAlp Dener       A = {L + U + F}
52c4b75bccSAlp Dener       IA = {i : i not in A}
53198282dbSAlp Dener 
54c4b75bccSAlp Dener     generate the reduced system Hr_k dr_k = -gr_k for variables in IA
55198282dbSAlp Dener     if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
56198282dbSAlp Dener       D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
57198282dbSAlp Dener       scale BFGS with VecReciprocal(D)
58198282dbSAlp Dener     end
59198282dbSAlp Dener     solve Hr_k dr_k = -gr_k
60198282dbSAlp Dener     set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
61198282dbSAlp Dener 
62198282dbSAlp Dener     x_{k+1} = VecMedian(x_k + d_k)
63198282dbSAlp Dener     s = x_{k+1} - x_k
64198282dbSAlp Dener     prered = dot(s, 0.5*gr_k - Hr_k*s)
65198282dbSAlp Dener     f_{k+1} = TaoComputeObjective(x_{k+1})
66198282dbSAlp Dener     actred = f_k - f_{k+1}
67198282dbSAlp Dener 
68198282dbSAlp Dener     oldTrust = trust
69198282dbSAlp Dener     step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION)
70198282dbSAlp Dener     if step_accepted
71198282dbSAlp Dener       g_{k+1} = TaoComputeGradient(x_{k+1})
72c4b75bccSAlp Dener       pg_{k+1} = project(g_{k+1})
73198282dbSAlp Dener       count the accepted Newton step
74198282dbSAlp Dener     else
75198282dbSAlp Dener       if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
76198282dbSAlp Dener         dr_k = -BFGS*gr_k for variables in I
77198282dbSAlp Dener         if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
78198282dbSAlp Dener           reset the BFGS preconditioner
79198282dbSAlp Dener           calculate scale delta and apply it to BFGS
80198282dbSAlp Dener           dr_k = -BFGS*gr_k for variables in I
81198282dbSAlp Dener           if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
82198282dbSAlp Dener             dr_k = -gr_k for variables in I
83198282dbSAlp Dener           end
84198282dbSAlp Dener         end
85198282dbSAlp Dener       end
86198282dbSAlp Dener 
87198282dbSAlp Dener       x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch()
88198282dbSAlp Dener       if ls_failed
89198282dbSAlp Dener         f_{k+1} = f_k
90198282dbSAlp Dener         x_{k+1} = x_k
91198282dbSAlp Dener         g_{k+1} = g_k
92198282dbSAlp Dener         pg_{k+1} = pg_k
93198282dbSAlp Dener         terminate
94198282dbSAlp Dener       else
95c4b75bccSAlp Dener         pg_{k+1} = project(g_{k+1})
96198282dbSAlp Dener         trust = oldTrust
97198282dbSAlp Dener         trust = TaoBNKUpdateTrustRadius(BNK_UPDATE_STEP)
98198282dbSAlp Dener         count the accepted step type (Newton, BFGS, scaled grad or grad)
99198282dbSAlp Dener       end
100198282dbSAlp Dener     end
101198282dbSAlp Dener 
102198282dbSAlp Dener     check convergence at pg_{k+1}
103198282dbSAlp Dener  end
104c14b763aSAlp Dener */
105c14b763aSAlp Dener 
TaoSolve_BNTL(Tao tao)106d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSolve_BNTL(Tao tao)
107d71ae5a4SJacob Faibussowitsch {
108c14b763aSAlp Dener   TAO_BNK                     *bnk = (TAO_BNK *)tao->data;
109e465cd6fSAlp Dener   KSPConvergedReason           ksp_reason;
110c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
111c14b763aSAlp Dener 
11289da521bSAlp Dener   PetscReal oldTrust, prered, actred, steplen, resnorm;
113937a31a1SAlp Dener   PetscBool cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE;
114c4b75bccSAlp Dener   PetscInt  stepType, nDiff;
115c14b763aSAlp Dener 
116c14b763aSAlp Dener   PetscFunctionBegin;
11728017e9fSAlp Dener   /* Initialize the preconditioner, KSP solver and trust radius/line search */
118c14b763aSAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1199566063dSJacob Faibussowitsch   PetscCall(TaoBNKInitialize(tao, bnk->init_type, &needH));
1203ba16761SJacob Faibussowitsch   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
121c14b763aSAlp Dener 
122c14b763aSAlp Dener   /* Have not converged; continue with Newton method */
123c14b763aSAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
124e1e80dc8SAlp Dener     /* Call general purpose update function */
125e1e80dc8SAlp Dener     if (tao->ops->update) {
126dbbe0bcdSBarry Smith       PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
127270bebe6SStefano Zampini       PetscCall(TaoComputeObjective(tao, tao->solution, &bnk->f));
128e1e80dc8SAlp Dener     }
12962675beeSAlp Dener 
13089da521bSAlp Dener     if (needH && bnk->inactive_idx) {
131e031d6f5SAlp Dener       /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
1329566063dSJacob Faibussowitsch       PetscCall(TaoBNKTakeCGSteps(tao, &cgTerminate));
133e031d6f5SAlp Dener       if (cgTerminate) {
134e031d6f5SAlp Dener         tao->reason = bnk->bncg->reason;
1353ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
136e031d6f5SAlp Dener       }
13708752603SAlp Dener       /* Compute the hessian and update the BFGS preconditioner at the new iterate */
1389566063dSJacob Faibussowitsch       PetscCall((*bnk->computehessian)(tao));
139937a31a1SAlp Dener       needH = PETSC_FALSE;
140937a31a1SAlp Dener     }
141c14b763aSAlp Dener 
1428d5ead36SAlp Dener     /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
1439566063dSJacob Faibussowitsch     PetscCall((*bnk->computestep)(tao, shift, &ksp_reason, &stepType));
144c14b763aSAlp Dener 
145c14b763aSAlp Dener     /* Store current solution before it changes */
146c14b763aSAlp Dener     oldTrust  = tao->trust;
147c14b763aSAlp Dener     bnk->fold = bnk->f;
1489566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, bnk->Xold));
1499566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, bnk->Gold));
1509566063dSJacob Faibussowitsch     PetscCall(VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old));
151c14b763aSAlp Dener 
152c14b763aSAlp Dener     /* Temporarily accept the step and project it into the bounds */
1539566063dSJacob Faibussowitsch     PetscCall(VecAXPY(tao->solution, 1.0, tao->stepdirection));
1549566063dSJacob Faibussowitsch     PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
155c14b763aSAlp Dener 
156c14b763aSAlp Dener     /* Check if the projection changed the step direction */
157c4b75bccSAlp Dener     if (nDiff > 0) {
158c4b75bccSAlp Dener       /* Projection changed the step, so we have to recompute the step and
159c4b75bccSAlp Dener          the predicted reduction. Leave the trust radius unchanged. */
1609566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->solution, tao->stepdirection));
1619566063dSJacob Faibussowitsch       PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold));
1629566063dSJacob Faibussowitsch       PetscCall(TaoBNKRecomputePred(tao, tao->stepdirection, &prered));
163c14b763aSAlp Dener     } else {
164c14b763aSAlp Dener       /* Step did not change, so we can just recover the pre-computed prediction */
1659566063dSJacob Faibussowitsch       PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
166c14b763aSAlp Dener     }
167c14b763aSAlp Dener     prered = -prered;
168c14b763aSAlp Dener 
169c14b763aSAlp Dener     /* Compute the actual reduction and update the trust radius */
1709566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(tao, tao->solution, &bnk->f));
171*76c63389SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(bnk->f), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
172c14b763aSAlp Dener     actred = bnk->fold - bnk->f;
1739566063dSJacob Faibussowitsch     PetscCall(TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted));
174c14b763aSAlp Dener 
175c14b763aSAlp Dener     if (stepAccepted) {
176c14b763aSAlp Dener       /* Step is good, evaluate the gradient and the hessian */
1778d5ead36SAlp Dener       steplen = 1.0;
178937a31a1SAlp Dener       needH   = PETSC_TRUE;
179e465cd6fSAlp Dener       ++bnk->newt;
1809566063dSJacob Faibussowitsch       PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient));
1819566063dSJacob Faibussowitsch       PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
1829566063dSJacob Faibussowitsch       PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
183976ed0a4SStefano Zampini       if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
1849566063dSJacob Faibussowitsch       PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
185c14b763aSAlp Dener     } else {
186c14b763aSAlp Dener       /* Trust-region rejected the step. Revert the solution. */
187c14b763aSAlp Dener       bnk->f = bnk->fold;
1889566063dSJacob Faibussowitsch       PetscCall(VecCopy(bnk->Xold, tao->solution));
189c14b763aSAlp Dener       /* Trigger the line search */
1909566063dSJacob Faibussowitsch       PetscCall(TaoBNKSafeguardStep(tao, ksp_reason, &stepType));
1919566063dSJacob Faibussowitsch       PetscCall(TaoBNKPerformLineSearch(tao, &stepType, &steplen, &ls_reason));
192c14b763aSAlp Dener       if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
193c14b763aSAlp Dener         /* Line search failed, revert solution and terminate */
194c0f10754SAlp Dener         stepAccepted = PETSC_FALSE;
195937a31a1SAlp Dener         needH        = PETSC_FALSE;
196c14b763aSAlp Dener         bnk->f       = bnk->fold;
1979566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->Xold, tao->solution));
1989566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->Gold, tao->gradient));
1999566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));
200c14b763aSAlp Dener         tao->trust  = 0.0;
201c14b763aSAlp Dener         tao->reason = TAO_DIVERGED_LS_FAILURE;
202c14b763aSAlp Dener       } else {
203937a31a1SAlp Dener         /* new iterate so we need to recompute the Hessian */
204937a31a1SAlp Dener         needH = PETSC_TRUE;
205198282dbSAlp Dener         /* compute the projected gradient */
2069566063dSJacob Faibussowitsch         PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
2079566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
208976ed0a4SStefano Zampini         if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
2099566063dSJacob Faibussowitsch         PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
210c14b763aSAlp Dener         /* Line search succeeded so we should update the trust radius based on the LS step length */
211770b7498SAlp Dener         tao->trust = oldTrust;
2129566063dSJacob Faibussowitsch         PetscCall(TaoBNKUpdateTrustRadius(tao, prered, actred, BNK_UPDATE_STEP, stepType, &stepAccepted));
21362675beeSAlp Dener         /* count the accepted step type */
2149566063dSJacob Faibussowitsch         PetscCall(TaoBNKAddStepCounts(tao, stepType));
215c14b763aSAlp Dener       }
216c14b763aSAlp Dener     }
217c14b763aSAlp Dener 
218c14b763aSAlp Dener     /*  Check for termination */
2199566063dSJacob Faibussowitsch     PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
2209566063dSJacob Faibussowitsch     PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
221*76c63389SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
2220f0abf79SStefano Zampini     ++tao->niter;
2239566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
2249566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen));
225dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
226c14b763aSAlp Dener   }
2273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
228c14b763aSAlp Dener }
229c14b763aSAlp Dener 
TaoSetUp_BNTL(Tao tao)230d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_BNTL(Tao tao)
231d71ae5a4SJacob Faibussowitsch {
2322e6e4ca1SStefano Zampini   KSP       ksp;
2330cd8b6e2SStefano Zampini   PetscBool valid;
2345eb5f4d6SAlp Dener 
2355eb5f4d6SAlp Dener   PetscFunctionBegin;
2369566063dSJacob Faibussowitsch   PetscCall(TaoSetUp_BNK(tao));
2379566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao, &ksp));
2380cd8b6e2SStefano Zampini   PetscCall(PetscObjectHasFunction((PetscObject)ksp, "KSPCGSetRadius_C", &valid));
2393c859ba3SBarry Smith   PetscCheck(valid, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Not for KSP type %s. Must use a trust-region CG method for KSP (e.g. KSPNASH, KSPSTCG, KSPGLTR)", ((PetscObject)ksp)->type_name);
2403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2415eb5f4d6SAlp Dener }
2425eb5f4d6SAlp Dener 
TaoSetFromOptions_BNTL(Tao tao,PetscOptionItems PetscOptionsObject)243ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_BNTL(Tao tao, PetscOptionItems PetscOptionsObject)
244d71ae5a4SJacob Faibussowitsch {
2459b6ef848SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
2469b6ef848SAlp Dener 
2479b6ef848SAlp Dener   PetscFunctionBegin;
248dbbe0bcdSBarry Smith   PetscCall(TaoSetFromOptions_BNK(tao, PetscOptionsObject));
249e0ed867bSAlp Dener   if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
2503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2519b6ef848SAlp Dener }
2529b6ef848SAlp Dener 
2533850be85SAlp Dener /*MC
2543850be85SAlp Dener   TAOBNTL - Bounded Newton Trust Region method with line-search fall-back for nonlinear
2553850be85SAlp Dener             minimization with bound constraints.
2569b6ef848SAlp Dener 
2573850be85SAlp Dener   Options Database Keys:
2583850be85SAlp Dener   + -tao_bnk_max_cg_its   - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
2593850be85SAlp Dener   . -tao_bnk_init_type   - trust radius initialization method ("constant", "direction", "interpolation")
2603850be85SAlp Dener   . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
2613850be85SAlp Dener   - -tao_bnk_as_type     - active-set estimation method ("none", "bertsekas")
2623850be85SAlp Dener 
2633850be85SAlp Dener   Level: beginner
264606f75f6SBarry Smith 
265606f75f6SBarry Smith   Developer Note:
266606f75f6SBarry Smith   One should control the maximum number of cg iterations through the standard pc_max_it option not with a special
267606f75f6SBarry Smith   ad hoc option
268606f75f6SBarry Smith 
2693850be85SAlp Dener M*/
TaoCreate_BNTL(Tao tao)270d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BNTL(Tao tao)
271d71ae5a4SJacob Faibussowitsch {
272c14b763aSAlp Dener   TAO_BNK *bnk;
273c14b763aSAlp Dener 
274c14b763aSAlp Dener   PetscFunctionBegin;
2759566063dSJacob Faibussowitsch   PetscCall(TaoCreate_BNK(tao));
276c14b763aSAlp Dener   tao->ops->solve          = TaoSolve_BNTL;
2775eb5f4d6SAlp Dener   tao->ops->setup          = TaoSetUp_BNTL;
278e0ed867bSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNTL;
279c14b763aSAlp Dener 
280c14b763aSAlp Dener   bnk              = (TAO_BNK *)tao->data;
281c14b763aSAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
2823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
283c14b763aSAlp Dener }
284