xref: /petsc/src/tao/bound/impls/bnk/bnls.c (revision 0875260307e97403c513ee2ad2ca062f01ff5a8a)
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  ------------------------------------------------------------
9eb910715SAlp Dener 
10198282dbSAlp Dener  x_0 = VecMedian(x_0)
11198282dbSAlp Dener  f_0, g_0 = TaoComputeObjectiveAndGradient(x_0)
12198282dbSAlp Dener  pg_0 = VecBoundGradientProjection(g_0)
13198282dbSAlp Dener  check convergence at pg_0
14198282dbSAlp Dener  trust = max_radius
15198282dbSAlp Dener  niter = 0
16198282dbSAlp Dener 
17198282dbSAlp Dener  while niter < max_it
18198282dbSAlp Dener     niter += 1
19198282dbSAlp Dener     H_k = TaoComputeHessian(x_k)
20198282dbSAlp Dener     if pc_type == BNK_PC_BFGS
21198282dbSAlp Dener       add correction to BFGS approx
22198282dbSAlp Dener       if scale_type == BNK_SCALE_AHESS
23198282dbSAlp Dener          D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
24198282dbSAlp Dener          scale BFGS with VecReciprocal(D)
25198282dbSAlp Dener       end
26198282dbSAlp Dener     end
27198282dbSAlp Dener 
28198282dbSAlp Dener     if pc_type = BNK_PC_BFGS
29198282dbSAlp Dener       B_k = BFGS
30198282dbSAlp Dener     else
31198282dbSAlp Dener       B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
32198282dbSAlp Dener       B_k = VecReciprocal(B_k)
33198282dbSAlp Dener     end
34198282dbSAlp Dener     w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
35198282dbSAlp Dener     eps = min(eps, norm2(w))
36198282dbSAlp Dener     determine the active and inactive index sets such that
37198282dbSAlp Dener       L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
38198282dbSAlp Dener       U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
39198282dbSAlp Dener       F = {i : l_i = (x_k)_i = u_i}
40198282dbSAlp Dener       A = {L + U + F}
41198282dbSAlp Dener       I = {i : i not in A}
42198282dbSAlp Dener 
43198282dbSAlp Dener     generate the reduced system Hr_k dr_k = -gr_k for variables in I
44198282dbSAlp Dener     if p > 0
45198282dbSAlp Dener       Hr_k += p*I
46198282dbSAlp Dener     end
47198282dbSAlp Dener     if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
48198282dbSAlp Dener       D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
49198282dbSAlp Dener       scale BFGS with VecReciprocal(D)
50198282dbSAlp Dener     end
51198282dbSAlp Dener     trust = max_radius
52198282dbSAlp Dener     solve Hr_k dr_k = -gr_k
53198282dbSAlp Dener     set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
54198282dbSAlp Dener 
55198282dbSAlp Dener     if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
56198282dbSAlp Dener       dr_k = -BFGS*gr_k for variables in I
57198282dbSAlp Dener       if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
58198282dbSAlp Dener         reset the BFGS preconditioner
59198282dbSAlp Dener         calculate scale delta and apply it to BFGS
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           dr_k = -gr_k for variables in I
63198282dbSAlp Dener         end
64198282dbSAlp Dener       end
65198282dbSAlp Dener     end
66198282dbSAlp Dener 
67198282dbSAlp Dener     x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch()
68198282dbSAlp Dener     if ls_failed
69198282dbSAlp Dener       f_{k+1} = f_k
70198282dbSAlp Dener       x_{k+1} = x_k
71198282dbSAlp Dener       g_{k+1} = g_k
72198282dbSAlp Dener       pg_{k+1} = pg_k
73198282dbSAlp Dener       terminate
74198282dbSAlp Dener     else
75198282dbSAlp Dener       pg_{k+1} = VecBoundGradientProjection(g_{k+1})
76198282dbSAlp Dener       count the accepted step type (Newton, BFGS, scaled grad or grad)
77198282dbSAlp Dener     end
78198282dbSAlp Dener 
79198282dbSAlp Dener     check convergence at pg_{k+1}
80198282dbSAlp Dener  end
81eb910715SAlp Dener */
82eb910715SAlp Dener 
83eb910715SAlp Dener static PetscErrorCode TaoSolve_BNLS(Tao tao)
84eb910715SAlp Dener {
85eb910715SAlp Dener   PetscErrorCode               ierr;
86eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
87e465cd6fSAlp Dener   KSPConvergedReason           ksp_reason;
88eb910715SAlp Dener   TaoLineSearchConvergedReason ls_reason;
89eb910715SAlp Dener 
909b6ef848SAlp Dener   PetscReal                    resnorm, steplen = 1.0;
91c0f10754SAlp Dener   PetscBool                    cgTerminate, stepAccepted = PETSC_TRUE, shift = PETSC_TRUE;
92eb910715SAlp Dener   PetscInt                     stepType;
93eb910715SAlp Dener 
94eb910715SAlp Dener   PetscFunctionBegin;
9528017e9fSAlp Dener   /* Initialize the preconditioner, KSP solver and trust radius/line search */
96eb910715SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
97c0f10754SAlp Dener   ierr = TaoBNKInitialize(tao, bnk->init_type, &stepAccepted);CHKERRQ(ierr);
9828017e9fSAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
99eb910715SAlp Dener 
100eb910715SAlp Dener   /* Have not converged; continue with Newton method */
101eb910715SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
102eb910715SAlp Dener     ++tao->niter;
103eb910715SAlp Dener     tao->ksp_its=0;
104eb910715SAlp Dener 
105c0f10754SAlp Dener     /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
106c0f10754SAlp Dener     ierr = TaoBNKTakeCGSteps(tao, &cgTerminate);CHKERRQ(ierr);
107c0f10754SAlp Dener     if (cgTerminate) {
108c0f10754SAlp Dener       tao->reason = bnk->bncg->reason;
109c0f10754SAlp Dener       PetscFunctionReturn(0);
110c0f10754SAlp Dener     }
111c0f10754SAlp Dener 
112*08752603SAlp Dener     /* Compute the hessian and update the BFGS preconditioner at the new iterate */
113*08752603SAlp Dener     if (stepAccepted) {ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);}
114fed79b8eSAlp Dener 
1158d5ead36SAlp Dener     /* Use the common BNK kernel to compute the safeguarded Newton step (for inactive variables only) */
11662675beeSAlp Dener     ierr = TaoBNKComputeStep(tao, shift, &ksp_reason);CHKERRQ(ierr);
117e465cd6fSAlp Dener     ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr);
118eb910715SAlp Dener 
119080d2917SAlp Dener     /* Store current solution before it changes */
120080d2917SAlp Dener     bnk->fold = bnk->f;
121eb910715SAlp Dener     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
122eb910715SAlp Dener     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
12309164190SAlp Dener     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
124eb910715SAlp Dener 
125c14b763aSAlp Dener     /* Trigger the line search */
126c14b763aSAlp Dener     ierr = TaoBNKPerformLineSearch(tao, stepType, &steplen, &ls_reason);CHKERRQ(ierr);
127eb910715SAlp Dener 
128eb910715SAlp Dener     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
129eb910715SAlp Dener       /* Failed to find an improving point */
130c0f10754SAlp Dener       stepAccepted = PETSC_FALSE;
131080d2917SAlp Dener       bnk->f = bnk->fold;
132eb910715SAlp Dener       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
133eb910715SAlp Dener       ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
13409164190SAlp Dener       ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
135c14b763aSAlp Dener       steplen = 0.0;
136eb910715SAlp Dener       tao->reason = TAO_DIVERGED_LS_FAILURE;
137e465cd6fSAlp Dener     } else {
138198282dbSAlp Dener       /* compute the projected gradient */
139198282dbSAlp Dener       ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
1409b6ef848SAlp Dener       ierr = VecNorm(tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
1419b6ef848SAlp Dener       if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF, 1, "User provided compute function generated Not-a-Number");
1429b6ef848SAlp Dener       /* update the trust radius based on the step length */
1439b6ef848SAlp Dener       ierr = TaoBNKUpdateTrustRadius(tao, 0.0, 0.0, BNK_UPDATE_STEP, stepType, &stepAccepted);CHKERRQ(ierr);
14462675beeSAlp Dener       /* count the accepted step type */
14562675beeSAlp Dener       ierr = TaoBNKAddStepCounts(tao, stepType);CHKERRQ(ierr);
146eb910715SAlp Dener     }
147eb910715SAlp Dener 
148eb910715SAlp Dener     /*  Check for termination */
1499b6ef848SAlp Dener     ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr);
1509b6ef848SAlp Dener     ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr);
1519b6ef848SAlp Dener     ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
1529b6ef848SAlp Dener     ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr);
153eb910715SAlp Dener     ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr);
154eb910715SAlp Dener   }
155eb910715SAlp Dener   PetscFunctionReturn(0);
156eb910715SAlp Dener }
157eb910715SAlp Dener 
158df278d8fSAlp Dener /*------------------------------------------------------------*/
159df278d8fSAlp Dener 
1609b6ef848SAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BNLS(Tao tao)
161eb910715SAlp Dener {
162fed79b8eSAlp Dener   TAO_BNK        *bnk;
163eb910715SAlp Dener   PetscErrorCode ierr;
164eb910715SAlp Dener 
165eb910715SAlp Dener   PetscFunctionBegin;
166eb910715SAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
167eb910715SAlp Dener   tao->ops->solve = TaoSolve_BNLS;
168fed79b8eSAlp Dener 
169fed79b8eSAlp Dener   bnk = (TAO_BNK *)tao->data;
170e031d6f5SAlp Dener   bnk->init_type = BNK_INIT_DIRECTION;
17166ed3702SAlp Dener   bnk->update_type = BNK_UPDATE_STEP; /* trust region updates based on line search step length */
172eb910715SAlp Dener   PetscFunctionReturn(0);
173eb910715SAlp Dener }