xref: /petsc/src/tao/bound/impls/bnk/bnls.c (revision 03047865b8d8757cf1cf9cda45785c1537b01dc1)
1eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
2eb910715SAlp Dener #include <petscksp.h>
3eb910715SAlp Dener 
4eb910715SAlp Dener /*
5198282dbSAlp Dener  Implements Newton's Method with a line search approach for
6198282dbSAlp Dener  solving bound constrained minimization problems.
7eb910715SAlp Dener 
8198282dbSAlp Dener  x_0 = VecMedian(x_0)
9198282dbSAlp Dener  f_0, g_0 = TaoComputeObjectiveAndGradient(x_0)
10c4b75bccSAlp Dener  pg_0 = project(g_0)
11198282dbSAlp Dener  check convergence at pg_0
12c4b75bccSAlp Dener  needH = TaoBNKInitialize(default:BNK_INIT_DIRECTION)
13198282dbSAlp Dener  niter = 0
14c4b75bccSAlp Dener  step_accepted = true
15198282dbSAlp Dener 
16198282dbSAlp Dener  while niter < max_it
17c4b75bccSAlp Dener     if needH
18c4b75bccSAlp Dener       If max_cg_steps > 0
19c4b75bccSAlp Dener         x_k, g_k, pg_k = TaoSolve(BNCG)
20c4b75bccSAlp Dener       end
21c4b75bccSAlp Dener 
22198282dbSAlp Dener       H_k = TaoComputeHessian(x_k)
23198282dbSAlp Dener       if pc_type == BNK_PC_BFGS
24198282dbSAlp Dener         add correction to BFGS approx
25198282dbSAlp Dener         if scale_type == BNK_SCALE_AHESS
26198282dbSAlp Dener           D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
27198282dbSAlp Dener           scale BFGS with VecReciprocal(D)
28198282dbSAlp Dener         end
29198282dbSAlp Dener       end
30c4b75bccSAlp Dener       needH = False
31c4b75bccSAlp Dener     end
32198282dbSAlp Dener 
33198282dbSAlp Dener     if pc_type = BNK_PC_BFGS
34198282dbSAlp Dener       B_k = BFGS
35198282dbSAlp Dener     else
36198282dbSAlp Dener       B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
37198282dbSAlp Dener       B_k = VecReciprocal(B_k)
38198282dbSAlp Dener     end
39198282dbSAlp Dener     w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
40198282dbSAlp Dener     eps = min(eps, norm2(w))
41198282dbSAlp Dener     determine the active and inactive index sets such that
42198282dbSAlp Dener       L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
43198282dbSAlp Dener       U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
44198282dbSAlp Dener       F = {i : l_i = (x_k)_i = u_i}
45198282dbSAlp Dener       A = {L + U + F}
46c4b75bccSAlp Dener       IA = {i : i not in A}
47198282dbSAlp Dener 
48c4b75bccSAlp Dener     generate the reduced system Hr_k dr_k = -gr_k for variables in IA
49198282dbSAlp Dener     if p > 0
50c4b75bccSAlp Dener       Hr_k += p*
51198282dbSAlp Dener     end
52198282dbSAlp Dener     if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
53198282dbSAlp Dener       D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
54198282dbSAlp Dener       scale BFGS with VecReciprocal(D)
55198282dbSAlp Dener     end
56198282dbSAlp Dener     solve Hr_k dr_k = -gr_k
57198282dbSAlp Dener     set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
58198282dbSAlp Dener 
59198282dbSAlp Dener     if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
60198282dbSAlp Dener       dr_k = -BFGS*gr_k for variables in I
61198282dbSAlp Dener       if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
62198282dbSAlp Dener         reset the BFGS preconditioner
63198282dbSAlp Dener         calculate scale delta and apply it to BFGS
64198282dbSAlp Dener         dr_k = -BFGS*gr_k for variables in I
65198282dbSAlp Dener         if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
66198282dbSAlp Dener           dr_k = -gr_k for variables in I
67198282dbSAlp Dener         end
68198282dbSAlp Dener       end
69198282dbSAlp Dener     end
70198282dbSAlp Dener 
71198282dbSAlp Dener     x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch()
72198282dbSAlp Dener     if ls_failed
73198282dbSAlp Dener       f_{k+1} = f_k
74198282dbSAlp Dener       x_{k+1} = x_k
75198282dbSAlp Dener       g_{k+1} = g_k
76198282dbSAlp Dener       pg_{k+1} = pg_k
77198282dbSAlp Dener       terminate
78198282dbSAlp Dener     else
79c4b75bccSAlp Dener       pg_{k+1} = project(g_{k+1})
80198282dbSAlp Dener       count the accepted step type (Newton, BFGS, scaled grad or grad)
81198282dbSAlp Dener     end
82198282dbSAlp Dener 
830279bc1bSStefano Zampini     niter += 1
84198282dbSAlp Dener     check convergence at pg_{k+1}
85198282dbSAlp Dener  end
86eb910715SAlp Dener */
87eb910715SAlp Dener 
TaoSolve_BNLS(Tao tao)88d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSolve_BNLS(Tao tao)
89d71ae5a4SJacob Faibussowitsch {
90eb910715SAlp Dener   TAO_BNK                     *bnk = (TAO_BNK *)tao->data;
91e465cd6fSAlp Dener   KSPConvergedReason           ksp_reason;
92eb910715SAlp Dener   TaoLineSearchConvergedReason ls_reason;
9389da521bSAlp Dener   PetscReal                    steplen = 1.0, resnorm;
94937a31a1SAlp Dener   PetscBool                    cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_TRUE;
95eb910715SAlp Dener   PetscInt                     stepType;
96eb910715SAlp Dener 
97eb910715SAlp Dener   PetscFunctionBegin;
9828017e9fSAlp Dener   /* Initialize the preconditioner, KSP solver and trust radius/line search */
99eb910715SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1009566063dSJacob Faibussowitsch   PetscCall(TaoBNKInitialize(tao, bnk->init_type, &needH));
1013ba16761SJacob Faibussowitsch   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
102eb910715SAlp Dener 
103eb910715SAlp Dener   /* Have not converged; continue with Newton method */
104eb910715SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
105e1e80dc8SAlp Dener     /* Call general purpose update function */
106e1e80dc8SAlp Dener     if (tao->ops->update) {
107dbbe0bcdSBarry Smith       PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
108270bebe6SStefano Zampini       PetscCall(TaoComputeObjective(tao, tao->solution, &bnk->f));
109e1e80dc8SAlp Dener     }
110eb910715SAlp Dener 
11189da521bSAlp Dener     if (needH && bnk->inactive_idx) {
112c0f10754SAlp Dener       /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
1139566063dSJacob Faibussowitsch       PetscCall(TaoBNKTakeCGSteps(tao, &cgTerminate));
114c0f10754SAlp Dener       if (cgTerminate) {
115c0f10754SAlp Dener         tao->reason = bnk->bncg->reason;
1163ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
117c0f10754SAlp Dener       }
11808752603SAlp Dener       /* Compute the hessian and update the BFGS preconditioner at the new iterate */
1199566063dSJacob Faibussowitsch       PetscCall((*bnk->computehessian)(tao));
120937a31a1SAlp Dener       needH = PETSC_FALSE;
121937a31a1SAlp Dener     }
122fed79b8eSAlp Dener 
1238d5ead36SAlp Dener     /* Use the common BNK kernel to compute the safeguarded Newton step (for inactive variables only) */
1249566063dSJacob Faibussowitsch     PetscCall((*bnk->computestep)(tao, shift, &ksp_reason, &stepType));
1259566063dSJacob Faibussowitsch     PetscCall(TaoBNKSafeguardStep(tao, ksp_reason, &stepType));
126eb910715SAlp Dener 
127080d2917SAlp Dener     /* Store current solution before it changes */
128080d2917SAlp Dener     bnk->fold = bnk->f;
1299566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, bnk->Xold));
1309566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, bnk->Gold));
1319566063dSJacob Faibussowitsch     PetscCall(VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old));
132eb910715SAlp Dener 
133c14b763aSAlp Dener     /* Trigger the line search */
1349566063dSJacob Faibussowitsch     PetscCall(TaoBNKPerformLineSearch(tao, &stepType, &steplen, &ls_reason));
135eb910715SAlp Dener 
136eb910715SAlp Dener     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
137eb910715SAlp Dener       /* Failed to find an improving point */
138937a31a1SAlp Dener       needH  = PETSC_FALSE;
139080d2917SAlp Dener       bnk->f = bnk->fold;
1409566063dSJacob Faibussowitsch       PetscCall(VecCopy(bnk->Xold, tao->solution));
1419566063dSJacob Faibussowitsch       PetscCall(VecCopy(bnk->Gold, tao->gradient));
1429566063dSJacob Faibussowitsch       PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));
143c14b763aSAlp Dener       steplen     = 0.0;
144eb910715SAlp Dener       tao->reason = TAO_DIVERGED_LS_FAILURE;
145e465cd6fSAlp Dener     } else {
146937a31a1SAlp Dener       /* new iterate so we need to recompute the Hessian */
147937a31a1SAlp Dener       needH = PETSC_TRUE;
148198282dbSAlp Dener       /* compute the projected gradient */
1499566063dSJacob Faibussowitsch       PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
1509566063dSJacob Faibussowitsch       PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
151976ed0a4SStefano Zampini       if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
1529566063dSJacob Faibussowitsch       PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
1539b6ef848SAlp Dener       /* update the trust radius based on the step length */
1549566063dSJacob Faibussowitsch       PetscCall(TaoBNKUpdateTrustRadius(tao, 0.0, 0.0, BNK_UPDATE_STEP, stepType, &stepAccepted));
15562675beeSAlp Dener       /* count the accepted step type */
1569566063dSJacob Faibussowitsch       PetscCall(TaoBNKAddStepCounts(tao, stepType));
157937a31a1SAlp Dener       /* active BNCG recycling for next iteration */
1589566063dSJacob Faibussowitsch       PetscCall(TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE));
159eb910715SAlp Dener     }
160eb910715SAlp Dener 
161eb910715SAlp Dener     /*  Check for termination */
1629566063dSJacob Faibussowitsch     PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
1639566063dSJacob Faibussowitsch     PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
164*76c63389SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
1650279bc1bSStefano Zampini     ++tao->niter;
1669566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
1679566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen));
168dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
169eb910715SAlp Dener   }
1703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
171eb910715SAlp Dener }
172eb910715SAlp Dener 
1733850be85SAlp Dener /*MC
1743850be85SAlp Dener   TAOBNLS - Bounded Newton Line Search for nonlinear minimization with bound constraints.
175df278d8fSAlp Dener 
1763850be85SAlp Dener   Options Database Keys:
1773850be85SAlp Dener + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
1783850be85SAlp Dener . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
1793850be85SAlp Dener . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
1803850be85SAlp Dener - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
1813850be85SAlp Dener 
1823850be85SAlp Dener   Level: beginner
1833850be85SAlp Dener M*/
TaoCreate_BNLS(Tao tao)184d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao)
185d71ae5a4SJacob Faibussowitsch {
186fed79b8eSAlp Dener   TAO_BNK *bnk;
187eb910715SAlp Dener 
188eb910715SAlp Dener   PetscFunctionBegin;
1899566063dSJacob Faibussowitsch   PetscCall(TaoCreate_BNK(tao));
190eb910715SAlp Dener   tao->ops->solve = TaoSolve_BNLS;
191fed79b8eSAlp Dener 
192fed79b8eSAlp Dener   bnk              = (TAO_BNK *)tao->data;
193e031d6f5SAlp Dener   bnk->init_type   = BNK_INIT_DIRECTION;
19466ed3702SAlp Dener   bnk->update_type = BNK_UPDATE_STEP; /* trust region updates based on line search step length */
1953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
196eb910715SAlp Dener }
197