xref: /petsc/src/tao/bound/impls/bnk/bntr.c (revision c4b75bccfe5bcbbfb1eb60f6996a18ddf23ce09b)
1fed79b8eSAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
2fed79b8eSAlp Dener #include <petscksp.h>
3fed79b8eSAlp Dener 
4fed79b8eSAlp Dener /*
5fed79b8eSAlp Dener  Implements Newton's Method with a trust region approach for solving
6fed79b8eSAlp Dener  bound constrained minimization problems.
7fed79b8eSAlp Dener 
8198282dbSAlp Dener  ------------------------------------------------------------
9198282dbSAlp Dener 
10198282dbSAlp Dener  x_0 = VecMedian(x_0)
11198282dbSAlp Dener  f_0, g_0= TaoComputeObjectiveAndGradient(x_0)
12*c4b75bccSAlp Dener  pg_0 = project(g_0)
13198282dbSAlp Dener  check convergence at pg_0
14*c4b75bccSAlp Dener  needH = TaoBNKInitialize(default:BNK_INIT_INTERPOLATION)
15198282dbSAlp Dener  niter = 0
16*c4b75bccSAlp Dener  step_accepted = false
17198282dbSAlp Dener 
18198282dbSAlp Dener  while niter <= max_it
19198282dbSAlp Dener     niter += 1
20*c4b75bccSAlp Dener 
21*c4b75bccSAlp Dener     if needH
22*c4b75bccSAlp Dener       If max_cg_steps > 0
23*c4b75bccSAlp Dener         x_k, g_k, pg_k = TaoSolve(BNCG)
24*c4b75bccSAlp Dener       end
25*c4b75bccSAlp Dener 
26198282dbSAlp Dener       H_k = TaoComputeHessian(x_k)
27198282dbSAlp Dener       if pc_type == BNK_PC_BFGS
28198282dbSAlp Dener         add correction to BFGS approx
29198282dbSAlp Dener         if scale_type == BNK_SCALE_AHESS
30198282dbSAlp Dener           D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
31198282dbSAlp Dener           scale BFGS with VecReciprocal(D)
32198282dbSAlp Dener         end
33198282dbSAlp Dener       end
34*c4b75bccSAlp Dener       needH = False
35198282dbSAlp Dener     end
36198282dbSAlp Dener 
37198282dbSAlp Dener     if pc_type = BNK_PC_BFGS
38198282dbSAlp Dener       B_k = BFGS
39198282dbSAlp Dener     else
40198282dbSAlp Dener       B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
41198282dbSAlp Dener       B_k = VecReciprocal(B_k)
42198282dbSAlp Dener     end
43198282dbSAlp Dener     w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
44198282dbSAlp Dener     eps = min(eps, norm2(w))
45198282dbSAlp Dener     determine the active and inactive index sets such that
46198282dbSAlp Dener       L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
47198282dbSAlp Dener       U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
48198282dbSAlp Dener       F = {i : l_i = (x_k)_i = u_i}
49198282dbSAlp Dener       A = {L + U + F}
50*c4b75bccSAlp Dener       IA = {i : i not in A}
51198282dbSAlp Dener 
52*c4b75bccSAlp Dener     generate the reduced system Hr_k dr_k = -gr_k for variables in IA
53198282dbSAlp Dener     if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
54198282dbSAlp Dener       D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
55198282dbSAlp Dener       scale BFGS with VecReciprocal(D)
56198282dbSAlp Dener     end
57*c4b75bccSAlp Dener 
58*c4b75bccSAlp Dener     while !stepAccepted
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})
72*c4b75bccSAlp Dener         pg_{k+1} = project(g_{k+1})
73198282dbSAlp Dener         count the accepted Newton step
74*c4b75bccSAlp Dener         needH = True
75198282dbSAlp Dener       else
76198282dbSAlp Dener         f_{k+1} = f_k
77198282dbSAlp Dener         x_{k+1} = x_k
78198282dbSAlp Dener         g_{k+1} = g_k
79198282dbSAlp Dener         pg_{k+1} = pg_k
80198282dbSAlp Dener         if trust == oldTrust
81198282dbSAlp Dener           terminate because we cannot shrink the radius any further
82198282dbSAlp Dener         end
83198282dbSAlp Dener       end
84198282dbSAlp Dener 
85198282dbSAlp Dener       check convergence at pg_{k+1}
86198282dbSAlp Dener     end
87*c4b75bccSAlp Dener 
88*c4b75bccSAlp Dener  end
89fed79b8eSAlp Dener */
90fed79b8eSAlp Dener 
91fed79b8eSAlp Dener static PetscErrorCode TaoSolve_BNTR(Tao tao)
92fed79b8eSAlp Dener {
93fed79b8eSAlp Dener   PetscErrorCode               ierr;
94fed79b8eSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
95e465cd6fSAlp Dener   KSPConvergedReason           ksp_reason;
96fed79b8eSAlp Dener 
97*c4b75bccSAlp Dener   PetscReal                    resnorm, oldTrust, prered, actred, steplen;
98937a31a1SAlp Dener   PetscBool                    cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE;
99*c4b75bccSAlp Dener   PetscInt                     stepType = BNK_NEWTON, nDiff;
100fed79b8eSAlp Dener 
101fed79b8eSAlp Dener   PetscFunctionBegin;
10228017e9fSAlp Dener   /* Initialize the preconditioner, KSP solver and trust radius/line search */
103fed79b8eSAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
104937a31a1SAlp Dener   ierr = TaoBNKInitialize(tao, bnk->init_type, &needH);CHKERRQ(ierr);
10528017e9fSAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
106fed79b8eSAlp Dener 
107fed79b8eSAlp Dener   /* Have not converged; continue with Newton method */
108fed79b8eSAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
109*c4b75bccSAlp Dener     ++tao->niter;
110e031d6f5SAlp Dener 
111937a31a1SAlp Dener     if (needH) {
112e031d6f5SAlp Dener       /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
113e031d6f5SAlp Dener       ierr = TaoBNKTakeCGSteps(tao, &cgTerminate);CHKERRQ(ierr);
114e031d6f5SAlp Dener       if (cgTerminate) {
115e031d6f5SAlp Dener         tao->reason = bnk->bncg->reason;
116e031d6f5SAlp Dener         PetscFunctionReturn(0);
117fed79b8eSAlp Dener       }
118937a31a1SAlp Dener       /* Compute the hessian and update the BFGS preconditioner at the new iterate */
119937a31a1SAlp Dener       ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
120937a31a1SAlp Dener       needH = PETSC_FALSE;
121937a31a1SAlp Dener     }
122fed79b8eSAlp Dener 
123fed79b8eSAlp Dener     /* Store current solution before it changes */
124fed79b8eSAlp Dener     oldTrust = tao->trust;
125fed79b8eSAlp Dener     bnk->fold = bnk->f;
126fed79b8eSAlp Dener     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
127fed79b8eSAlp Dener     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
128fed79b8eSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
129fed79b8eSAlp Dener 
130937a31a1SAlp Dener     /* Enter into trust region loops */
131937a31a1SAlp Dener     stepAccepted = PETSC_FALSE;
132937a31a1SAlp Dener     while (!stepAccepted && tao->reason == TAO_CONTINUE_ITERATING) {
133937a31a1SAlp Dener       tao->ksp_its=0;
134937a31a1SAlp Dener 
135937a31a1SAlp Dener       /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
136937a31a1SAlp Dener       ierr = TaoBNKComputeStep(tao, shift, &ksp_reason);CHKERRQ(ierr);
137937a31a1SAlp Dener 
138b1c2d0e3SAlp Dener       /* Temporarily accept the step and project it into the bounds */
139fed79b8eSAlp Dener       ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
140*c4b75bccSAlp Dener       ierr = TaoBoundSolution(tao->XL, tao->XU, tao->solution, &nDiff);CHKERRQ(ierr);
141b1c2d0e3SAlp Dener 
142b1c2d0e3SAlp Dener       /* Check if the projection changed the step direction */
143*c4b75bccSAlp Dener       if (nDiff > 0) {
144*c4b75bccSAlp Dener         /* Projection changed the step, so we have to recompute the step and
145*c4b75bccSAlp Dener            the predicted reduction. Leave the trust radius unchanged. */
146b1c2d0e3SAlp Dener         ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
1478d5ead36SAlp Dener         ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
1485e9b73cbSAlp Dener         ierr = TaoBNKRecomputePred(tao, tao->stepdirection, &prered);CHKERRQ(ierr);
149b1c2d0e3SAlp Dener       } else {
150b1c2d0e3SAlp Dener         /* Step did not change, so we can just recover the pre-computed prediction */
151b1c2d0e3SAlp Dener         ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
152b1c2d0e3SAlp Dener       }
153b1c2d0e3SAlp Dener       prered = -prered;
154b1c2d0e3SAlp Dener 
155b1c2d0e3SAlp Dener       /* Compute the actual reduction and update the trust radius */
156fed79b8eSAlp Dener       ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
157b1c2d0e3SAlp Dener       actred = bnk->fold - bnk->f;
158e761ccfdSAlp Dener       oldTrust = tao->trust;
15928017e9fSAlp Dener       ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr);
160fed79b8eSAlp Dener 
161fed79b8eSAlp Dener       if (stepAccepted) {
162937a31a1SAlp Dener         /* Step is good, evaluate the gradient and flip the need-Hessian switch */
1638d5ead36SAlp Dener         steplen = 1.0;
164937a31a1SAlp Dener         needH = PETSC_TRUE;
165e465cd6fSAlp Dener         ++bnk->newt;
166fed79b8eSAlp Dener         ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
16761be54a6SAlp Dener         ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
16861be54a6SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
16961be54a6SAlp Dener         ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
1709b6ef848SAlp Dener         ierr = VecNorm(tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
1719b6ef848SAlp Dener         if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
172fed79b8eSAlp Dener       } else {
173fed79b8eSAlp Dener         /* Step is bad, revert old solution and re-solve with new radius*/
1748d5ead36SAlp Dener         steplen = 0.0;
175937a31a1SAlp Dener         needH = PETSC_FALSE;
176fed79b8eSAlp Dener         bnk->f = bnk->fold;
177fed79b8eSAlp Dener         ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
178fed79b8eSAlp Dener         ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
179fed79b8eSAlp Dener         ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
18073e4db90SAlp Dener         if (oldTrust == tao->trust) {
18173e4db90SAlp Dener           /* Can't change the radius anymore so just terminate */
182fed79b8eSAlp Dener           tao->reason = TAO_DIVERGED_TR_REDUCTION;
183fed79b8eSAlp Dener         }
184fed79b8eSAlp Dener       }
185fed79b8eSAlp Dener 
186fed79b8eSAlp Dener       /*  Check for termination */
1879b6ef848SAlp Dener       ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr);
1889b6ef848SAlp Dener       ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr);
1899b6ef848SAlp Dener       ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
1909b6ef848SAlp Dener       ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr);
191fed79b8eSAlp Dener       ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr);
192fed79b8eSAlp Dener     }
193937a31a1SAlp Dener   }
194fed79b8eSAlp Dener   PetscFunctionReturn(0);
195fed79b8eSAlp Dener }
196fed79b8eSAlp Dener 
197df278d8fSAlp Dener /*------------------------------------------------------------*/
198df278d8fSAlp Dener 
1999b6ef848SAlp Dener PETSC_INTERN PetscErrorCode TaoSetUp_BNTR(Tao tao)
2009b6ef848SAlp Dener {
2019b6ef848SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
2029b6ef848SAlp Dener   PetscErrorCode ierr;
2039b6ef848SAlp Dener 
2049b6ef848SAlp Dener   PetscFunctionBegin;
2059b6ef848SAlp Dener   ierr = TaoSetUp_BNK(tao);CHKERRQ(ierr);
2069b6ef848SAlp 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)");
2079b6ef848SAlp Dener   PetscFunctionReturn(0);
2089b6ef848SAlp Dener }
2099b6ef848SAlp Dener 
2109b6ef848SAlp Dener /*------------------------------------------------------------*/
2119b6ef848SAlp Dener 
2129b6ef848SAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BNTR(Tao tao)
213fed79b8eSAlp Dener {
214fed79b8eSAlp Dener   TAO_BNK        *bnk;
215fed79b8eSAlp Dener   PetscErrorCode ierr;
216fed79b8eSAlp Dener 
217fed79b8eSAlp Dener   PetscFunctionBegin;
218fed79b8eSAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
219fed79b8eSAlp Dener   tao->ops->solve=TaoSolve_BNTR;
2209b6ef848SAlp Dener   tao->ops->setup=TaoSetUp_BNTR;
221fed79b8eSAlp Dener 
222fed79b8eSAlp Dener   bnk = (TAO_BNK *)tao->data;
22366ed3702SAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
224fed79b8eSAlp Dener   PetscFunctionReturn(0);
225fed79b8eSAlp Dener }