xref: /petsc/src/tao/bound/impls/bnk/bntl.c (revision 976ed0a4a7e2c9a504a034b07aaf489e2c2d55c5)
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 
108d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSolve_BNTL(Tao tao)
109d71ae5a4SJacob Faibussowitsch {
110c14b763aSAlp Dener   TAO_BNK                     *bnk = (TAO_BNK *)tao->data;
111e465cd6fSAlp Dener   KSPConvergedReason           ksp_reason;
112c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
113c14b763aSAlp Dener 
11489da521bSAlp Dener   PetscReal oldTrust, prered, actred, steplen, resnorm;
115937a31a1SAlp Dener   PetscBool cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE;
116c4b75bccSAlp Dener   PetscInt  stepType, nDiff;
117c14b763aSAlp Dener 
118c14b763aSAlp Dener   PetscFunctionBegin;
11928017e9fSAlp Dener   /* Initialize the preconditioner, KSP solver and trust radius/line search */
120c14b763aSAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1219566063dSJacob Faibussowitsch   PetscCall(TaoBNKInitialize(tao, bnk->init_type, &needH));
1223ba16761SJacob Faibussowitsch   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
123c14b763aSAlp Dener 
124c14b763aSAlp Dener   /* Have not converged; continue with Newton method */
125c14b763aSAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
126e1e80dc8SAlp Dener     /* Call general purpose update function */
127e1e80dc8SAlp Dener     if (tao->ops->update) {
128dbbe0bcdSBarry Smith       PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
1297494f0b1SStefano Zampini       PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient));
130e1e80dc8SAlp Dener     }
13162675beeSAlp Dener 
13289da521bSAlp Dener     if (needH && bnk->inactive_idx) {
133e031d6f5SAlp Dener       /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
1349566063dSJacob Faibussowitsch       PetscCall(TaoBNKTakeCGSteps(tao, &cgTerminate));
135e031d6f5SAlp Dener       if (cgTerminate) {
136e031d6f5SAlp Dener         tao->reason = bnk->bncg->reason;
1373ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
138e031d6f5SAlp Dener       }
13908752603SAlp Dener       /* Compute the hessian and update the BFGS preconditioner at the new iterate */
1409566063dSJacob Faibussowitsch       PetscCall((*bnk->computehessian)(tao));
141937a31a1SAlp Dener       needH = PETSC_FALSE;
142937a31a1SAlp Dener     }
143c14b763aSAlp Dener 
1448d5ead36SAlp Dener     /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
1459566063dSJacob Faibussowitsch     PetscCall((*bnk->computestep)(tao, shift, &ksp_reason, &stepType));
146c14b763aSAlp Dener 
147c14b763aSAlp Dener     /* Store current solution before it changes */
148c14b763aSAlp Dener     oldTrust  = tao->trust;
149c14b763aSAlp Dener     bnk->fold = bnk->f;
1509566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, bnk->Xold));
1519566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, bnk->Gold));
1529566063dSJacob Faibussowitsch     PetscCall(VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old));
153c14b763aSAlp Dener 
154c14b763aSAlp Dener     /* Temporarily accept the step and project it into the bounds */
1559566063dSJacob Faibussowitsch     PetscCall(VecAXPY(tao->solution, 1.0, tao->stepdirection));
1569566063dSJacob Faibussowitsch     PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
157c14b763aSAlp Dener 
158c14b763aSAlp Dener     /* Check if the projection changed the step direction */
159c4b75bccSAlp Dener     if (nDiff > 0) {
160c4b75bccSAlp Dener       /* Projection changed the step, so we have to recompute the step and
161c4b75bccSAlp Dener          the predicted reduction. Leave the trust radius unchanged. */
1629566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->solution, tao->stepdirection));
1639566063dSJacob Faibussowitsch       PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold));
1649566063dSJacob Faibussowitsch       PetscCall(TaoBNKRecomputePred(tao, tao->stepdirection, &prered));
165c14b763aSAlp Dener     } else {
166c14b763aSAlp Dener       /* Step did not change, so we can just recover the pre-computed prediction */
1679566063dSJacob Faibussowitsch       PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
168c14b763aSAlp Dener     }
169c14b763aSAlp Dener     prered = -prered;
170c14b763aSAlp Dener 
171c14b763aSAlp Dener     /* Compute the actual reduction and update the trust radius */
1729566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(tao, tao->solution, &bnk->f));
1733c859ba3SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(bnk->f), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
174c14b763aSAlp Dener     actred = bnk->fold - bnk->f;
1759566063dSJacob Faibussowitsch     PetscCall(TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted));
176c14b763aSAlp Dener 
177c14b763aSAlp Dener     if (stepAccepted) {
178c14b763aSAlp Dener       /* Step is good, evaluate the gradient and the hessian */
1798d5ead36SAlp Dener       steplen = 1.0;
180937a31a1SAlp Dener       needH   = PETSC_TRUE;
181e465cd6fSAlp Dener       ++bnk->newt;
1829566063dSJacob Faibussowitsch       PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient));
1839566063dSJacob Faibussowitsch       PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
1849566063dSJacob Faibussowitsch       PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
185*976ed0a4SStefano Zampini       if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
1869566063dSJacob Faibussowitsch       PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
187c14b763aSAlp Dener     } else {
188c14b763aSAlp Dener       /* Trust-region rejected the step. Revert the solution. */
189c14b763aSAlp Dener       bnk->f = bnk->fold;
1909566063dSJacob Faibussowitsch       PetscCall(VecCopy(bnk->Xold, tao->solution));
191c14b763aSAlp Dener       /* Trigger the line search */
1929566063dSJacob Faibussowitsch       PetscCall(TaoBNKSafeguardStep(tao, ksp_reason, &stepType));
1939566063dSJacob Faibussowitsch       PetscCall(TaoBNKPerformLineSearch(tao, &stepType, &steplen, &ls_reason));
194c14b763aSAlp Dener       if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
195c14b763aSAlp Dener         /* Line search failed, revert solution and terminate */
196c0f10754SAlp Dener         stepAccepted = PETSC_FALSE;
197937a31a1SAlp Dener         needH        = PETSC_FALSE;
198c14b763aSAlp Dener         bnk->f       = bnk->fold;
1999566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->Xold, tao->solution));
2009566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->Gold, tao->gradient));
2019566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));
202c14b763aSAlp Dener         tao->trust  = 0.0;
203c14b763aSAlp Dener         tao->reason = TAO_DIVERGED_LS_FAILURE;
204c14b763aSAlp Dener       } else {
205937a31a1SAlp Dener         /* new iterate so we need to recompute the Hessian */
206937a31a1SAlp Dener         needH = PETSC_TRUE;
207198282dbSAlp Dener         /* compute the projected gradient */
2089566063dSJacob Faibussowitsch         PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
2099566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
210*976ed0a4SStefano Zampini         if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
2119566063dSJacob Faibussowitsch         PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
212c14b763aSAlp Dener         /* Line search succeeded so we should update the trust radius based on the LS step length */
213770b7498SAlp Dener         tao->trust = oldTrust;
2149566063dSJacob Faibussowitsch         PetscCall(TaoBNKUpdateTrustRadius(tao, prered, actred, BNK_UPDATE_STEP, stepType, &stepAccepted));
21562675beeSAlp Dener         /* count the accepted step type */
2169566063dSJacob Faibussowitsch         PetscCall(TaoBNKAddStepCounts(tao, stepType));
217c14b763aSAlp Dener       }
218c14b763aSAlp Dener     }
219c14b763aSAlp Dener 
220c14b763aSAlp Dener     /*  Check for termination */
2219566063dSJacob Faibussowitsch     PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
2229566063dSJacob Faibussowitsch     PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
2233c859ba3SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
2240f0abf79SStefano Zampini     ++tao->niter;
2259566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
2269566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen));
227dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
228c14b763aSAlp Dener   }
2293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
230c14b763aSAlp Dener }
231c14b763aSAlp Dener 
232df278d8fSAlp Dener /*------------------------------------------------------------*/
233d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_BNTL(Tao tao)
234d71ae5a4SJacob Faibussowitsch {
2352e6e4ca1SStefano Zampini   KSP               ksp;
2362e6e4ca1SStefano Zampini   PetscVoidFunction valid;
2375eb5f4d6SAlp Dener 
2385eb5f4d6SAlp Dener   PetscFunctionBegin;
2399566063dSJacob Faibussowitsch   PetscCall(TaoSetUp_BNK(tao));
2409566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao, &ksp));
2419566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)ksp, "KSPCGSetRadius_C", &valid));
2423c859ba3SBarry 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);
2433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2445eb5f4d6SAlp Dener }
2455eb5f4d6SAlp Dener 
2465eb5f4d6SAlp Dener /*------------------------------------------------------------*/
247d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetFromOptions_BNTL(Tao tao, PetscOptionItems *PetscOptionsObject)
248d71ae5a4SJacob Faibussowitsch {
2499b6ef848SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
2509b6ef848SAlp Dener 
2519b6ef848SAlp Dener   PetscFunctionBegin;
252dbbe0bcdSBarry Smith   PetscCall(TaoSetFromOptions_BNK(tao, PetscOptionsObject));
253e0ed867bSAlp Dener   if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
2543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2559b6ef848SAlp Dener }
2569b6ef848SAlp Dener 
2579b6ef848SAlp Dener /*------------------------------------------------------------*/
2583850be85SAlp Dener /*MC
2593850be85SAlp Dener   TAOBNTL - Bounded Newton Trust Region method with line-search fall-back for nonlinear
2603850be85SAlp Dener             minimization with bound constraints.
2619b6ef848SAlp Dener 
2623850be85SAlp Dener   Options Database Keys:
2633850be85SAlp Dener   + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
2643850be85SAlp Dener   . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
2653850be85SAlp Dener   . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
2663850be85SAlp Dener   - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
2673850be85SAlp Dener 
2683850be85SAlp Dener   Level: beginner
2693850be85SAlp Dener M*/
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