xref: /petsc/src/tao/bound/impls/bnk/bnls.c (revision 198282db2bfc4cc7307b0b13cb98d7a27b116ed7)
1eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
2eb910715SAlp Dener #include <petscksp.h>
3eb910715SAlp Dener 
4eb910715SAlp Dener /*
5*198282dbSAlp Dener  Implements Newton's Method with a line search approach for
6*198282dbSAlp Dener  solving bound constrained minimization problems.
7eb910715SAlp Dener 
8*198282dbSAlp Dener  ------------------------------------------------------------
9eb910715SAlp Dener 
10*198282dbSAlp Dener  x_0 = VecMedian(x_0)
11*198282dbSAlp Dener  f_0, g_0 = TaoComputeObjectiveAndGradient(x_0)
12*198282dbSAlp Dener  pg_0 = VecBoundGradientProjection(g_0)
13*198282dbSAlp Dener  check convergence at pg_0
14*198282dbSAlp Dener  trust = max_radius
15*198282dbSAlp Dener  niter = 0
16*198282dbSAlp Dener 
17*198282dbSAlp Dener  while niter < max_it
18*198282dbSAlp Dener     niter += 1
19*198282dbSAlp Dener     H_k = TaoComputeHessian(x_k)
20*198282dbSAlp Dener     if pc_type == BNK_PC_BFGS
21*198282dbSAlp Dener       add correction to BFGS approx
22*198282dbSAlp Dener       if scale_type == BNK_SCALE_AHESS
23*198282dbSAlp Dener          D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
24*198282dbSAlp Dener          scale BFGS with VecReciprocal(D)
25*198282dbSAlp Dener       end
26*198282dbSAlp Dener     end
27*198282dbSAlp Dener 
28*198282dbSAlp Dener     if pc_type = BNK_PC_BFGS
29*198282dbSAlp Dener       B_k = BFGS
30*198282dbSAlp Dener     else
31*198282dbSAlp Dener       B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
32*198282dbSAlp Dener       B_k = VecReciprocal(B_k)
33*198282dbSAlp Dener     end
34*198282dbSAlp Dener     w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
35*198282dbSAlp Dener     eps = min(eps, norm2(w))
36*198282dbSAlp Dener     determine the active and inactive index sets such that
37*198282dbSAlp Dener       L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
38*198282dbSAlp Dener       U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
39*198282dbSAlp Dener       F = {i : l_i = (x_k)_i = u_i}
40*198282dbSAlp Dener       A = {L + U + F}
41*198282dbSAlp Dener       I = {i : i not in A}
42*198282dbSAlp Dener 
43*198282dbSAlp Dener     generate the reduced system Hr_k dr_k = -gr_k for variables in I
44*198282dbSAlp Dener     if p > 0
45*198282dbSAlp Dener       Hr_k += p*I
46*198282dbSAlp Dener     end
47*198282dbSAlp Dener     if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
48*198282dbSAlp Dener       D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
49*198282dbSAlp Dener       scale BFGS with VecReciprocal(D)
50*198282dbSAlp Dener     end
51*198282dbSAlp Dener     trust = max_radius
52*198282dbSAlp Dener     solve Hr_k dr_k = -gr_k
53*198282dbSAlp Dener     set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
54*198282dbSAlp Dener 
55*198282dbSAlp Dener     if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
56*198282dbSAlp Dener       dr_k = -BFGS*gr_k for variables in I
57*198282dbSAlp Dener       if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
58*198282dbSAlp Dener         reset the BFGS preconditioner
59*198282dbSAlp Dener         calculate scale delta and apply it to BFGS
60*198282dbSAlp Dener         dr_k = -BFGS*gr_k for variables in I
61*198282dbSAlp Dener         if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
62*198282dbSAlp Dener           dr_k = -gr_k for variables in I
63*198282dbSAlp Dener         end
64*198282dbSAlp Dener       end
65*198282dbSAlp Dener     end
66*198282dbSAlp Dener 
67*198282dbSAlp Dener     x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch()
68*198282dbSAlp Dener     if ls_failed
69*198282dbSAlp Dener       f_{k+1} = f_k
70*198282dbSAlp Dener       x_{k+1} = x_k
71*198282dbSAlp Dener       g_{k+1} = g_k
72*198282dbSAlp Dener       pg_{k+1} = pg_k
73*198282dbSAlp Dener       terminate
74*198282dbSAlp Dener     else
75*198282dbSAlp Dener       pg_{k+1} = VecBoundGradientProjection(g_{k+1})
76*198282dbSAlp Dener       count the accepted step type (Newton, BFGS, scaled grad or grad)
77*198282dbSAlp Dener     end
78*198282dbSAlp Dener 
79*198282dbSAlp Dener     check convergence at pg_{k+1}
80*198282dbSAlp 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 
90c14b763aSAlp Dener   PetscReal                    steplen = 1.0;
9162675beeSAlp Dener   PetscBool                    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;
9762675beeSAlp Dener   ierr = TaoBNKInitialize(tao, BNK_INIT_CONSTANT);CHKERRQ(ierr);
98*198282dbSAlp Dener   tao->trust = bnk->max_radius;
9928017e9fSAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
100eb910715SAlp Dener 
101eb910715SAlp Dener   /* Have not converged; continue with Newton method */
102eb910715SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
103eb910715SAlp Dener     ++tao->niter;
104eb910715SAlp Dener     tao->ksp_its=0;
105eb910715SAlp Dener 
10662675beeSAlp Dener     /* Compute the hessian and update the BFGS preconditioner at the new iterate*/
10762675beeSAlp Dener     ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
108fed79b8eSAlp Dener 
1098d5ead36SAlp Dener     /* Use the common BNK kernel to compute the safeguarded Newton step (for inactive variables only) */
11062675beeSAlp Dener     ierr = TaoBNKComputeStep(tao, shift, &ksp_reason);CHKERRQ(ierr);
111e465cd6fSAlp Dener     ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr);
112eb910715SAlp Dener 
113080d2917SAlp Dener     /* Store current solution before it changes */
114080d2917SAlp Dener     bnk->fold = bnk->f;
115eb910715SAlp Dener     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
116eb910715SAlp Dener     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
11709164190SAlp Dener     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
118eb910715SAlp Dener 
119c14b763aSAlp Dener     /* Trigger the line search */
120c14b763aSAlp Dener     ierr = TaoBNKPerformLineSearch(tao, stepType, &steplen, &ls_reason);CHKERRQ(ierr);
121eb910715SAlp Dener 
122eb910715SAlp Dener     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
123eb910715SAlp Dener       /* Failed to find an improving point */
124080d2917SAlp Dener       bnk->f = bnk->fold;
125eb910715SAlp Dener       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
126eb910715SAlp Dener       ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
12709164190SAlp Dener       ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
128c14b763aSAlp Dener       steplen = 0.0;
129eb910715SAlp Dener       tao->reason = TAO_DIVERGED_LS_FAILURE;
130e465cd6fSAlp Dener     } else {
131*198282dbSAlp Dener       /* compute the projected gradient */
132*198282dbSAlp Dener       ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
13362675beeSAlp Dener       /* count the accepted step type */
13462675beeSAlp Dener       ierr = TaoBNKAddStepCounts(tao, stepType);CHKERRQ(ierr);
135eb910715SAlp Dener     }
136eb910715SAlp Dener 
137eb910715SAlp Dener     /*  Check for termination */
138eb910715SAlp Dener     ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
139eb910715SAlp Dener     if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF, 1, "User provided compute function generated Not-a-Number");
140eb910715SAlp Dener     ierr = TaoLogConvergenceHistory(tao, bnk->f, bnk->gnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
141c14b763aSAlp Dener     ierr = TaoMonitor(tao, tao->niter, bnk->f, bnk->gnorm, 0.0, steplen);CHKERRQ(ierr);
142eb910715SAlp Dener     ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr);
143eb910715SAlp Dener   }
144eb910715SAlp Dener   PetscFunctionReturn(0);
145eb910715SAlp Dener }
146eb910715SAlp Dener 
147df278d8fSAlp Dener /*------------------------------------------------------------*/
148df278d8fSAlp Dener 
149eb910715SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao)
150eb910715SAlp Dener {
151fed79b8eSAlp Dener   TAO_BNK        *bnk;
152eb910715SAlp Dener   PetscErrorCode ierr;
153eb910715SAlp Dener 
154eb910715SAlp Dener   PetscFunctionBegin;
155eb910715SAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
156eb910715SAlp Dener   tao->ops->solve = TaoSolve_BNLS;
157fed79b8eSAlp Dener 
158fed79b8eSAlp Dener   bnk = (TAO_BNK *)tao->data;
15928017e9fSAlp Dener   bnk->init_type = BNK_INIT_CONSTANT;
16066ed3702SAlp Dener   bnk->update_type = BNK_UPDATE_STEP; /* trust region updates based on line search step length */
161eb910715SAlp Dener   PetscFunctionReturn(0);
162eb910715SAlp Dener }