xref: /petsc/src/tao/bound/impls/bnk/bnls.c (revision b1c2d0e33e2826e0558267820284dcdf0eb67a2f)
1 #include <../src/tao/bound/impls/bnk/bnk.h>
2 #include <petscksp.h>
3 
4 /*
5  Implements Newton's Method with a line search approach for solving
6  bound constrained minimization problems. A projected More'-Thuente line
7  search is used to guarantee that the bfgs preconditioner remains positive
8  definite.
9 
10  The method can shift the Hessian matrix. The shifting procedure is
11  adapted from the PATH algorithm for solving complementarity
12  problems.
13 
14  The linear system solve should be done with a conjugate gradient
15  method, although any method can be used.
16 */
17 
18 static PetscErrorCode TaoSolve_BNLS(Tao tao)
19 {
20   PetscErrorCode               ierr;
21   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
22   TaoLineSearchConvergedReason ls_reason;
23 
24   PetscReal                    f_full, gdx, prered, actred;
25   PetscReal                    step = 1.0;
26   PetscReal                    delta;
27   PetscReal                    e_min;
28 
29   PetscBool                    trustAccept;
30   PetscInt                     stepType;
31   PetscInt                     bfgsUpdates = 0;
32 
33   PetscFunctionBegin;
34   /*   Project the current point onto the feasible set */
35   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
36   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
37 
38   /* Project the initial point onto the feasible region */
39   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
40 
41   /* Check convergence criteria */
42   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
43   ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
44   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
45   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
46 
47   tao->reason = TAO_CONTINUE_ITERATING;
48   ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
49   ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr);
50   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
51   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
52 
53   /* Initialize the preconditioner and trust radius */
54   ierr = TaoBNKInitialize(tao);CHKERRQ(ierr);
55 
56   /* Have not converged; continue with Newton method */
57   while (tao->reason == TAO_CONTINUE_ITERATING) {
58     ++tao->niter;
59     tao->ksp_its=0;
60 
61     /* Compute the Hessian */
62     ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
63 
64     /* Update the BFGS preconditioner */
65     if (BNK_PC_BFGS == bnk->pc_type) {
66       if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) {
67         /* Obtain diagonal for the bfgs preconditioner  */
68         ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
69         ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
70         ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
71         ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
72       }
73       /* Update the limited memory preconditioner and get existing # of updates */
74       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
75       ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
76     }
77 
78     /* Use the common BNK kernel to compute the safeguarded Newton step */
79     ierr = TaoBNKComputeStep(tao, PETSC_TRUE, &stepType);CHKERRQ(ierr);
80 
81     /* Store current solution before it changes */
82     bnk->fold = bnk->f;
83     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
84     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
85     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
86 
87     /* Perform the linesearch */
88     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr);
89     ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
90     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
91 
92     while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != BNK_GRADIENT) {
93       /* Linesearch failed, revert solution */
94       bnk->f = bnk->fold;
95       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
96       ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
97       ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
98 
99       switch(stepType) {
100       case BNK_NEWTON:
101         /* Failed to obtain acceptable iterate with Newton 1step
102            Update the perturbation for next time */
103         if (bnk->pert <= 0.0) {
104           /* Initialize the perturbation */
105           bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
106           if (bnk->is_gltr) {
107             ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
108             bnk->pert = PetscMax(bnk->pert, -e_min);
109           }
110         } else {
111           /* Increase the perturbation */
112           bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
113         }
114 
115         if (BNK_PC_BFGS != bnk->pc_type) {
116           /* We don't have the bfgs matrix around and being updated
117              Must use gradient direction in this case */
118           ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
119           ++bnk->grad;
120           stepType = BNK_GRADIENT;
121         } else {
122           /* Attempt to use the BFGS direction */
123           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
124           ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
125           /* Check for success (descent direction) */
126           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
127           if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
128             /* BFGS direction is not descent or direction produced not a number
129                We can assert bfgsUpdates > 1 in this case
130                Use steepest descent direction (scaled) */
131 
132             if (bnk->f != 0.0) {
133               delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
134             } else {
135               delta = 2.0 / (bnk->gnorm*bnk->gnorm);
136             }
137             ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
138             ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
139             ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
140             ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
141             ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
142 
143             bfgsUpdates = 1;
144             ++bnk->sgrad;
145             stepType = BNK_SCALED_GRADIENT;
146           } else {
147             ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
148             if (1 == bfgsUpdates) {
149               /* The first BFGS direction is always the scaled gradient */
150               ++bnk->sgrad;
151               stepType = BNK_SCALED_GRADIENT;
152             } else {
153               ++bnk->bfgs;
154               stepType = BNK_BFGS;
155             }
156           }
157         }
158         break;
159 
160       case BNK_BFGS:
161         /* Can only enter if pc_type == BNK_PC_BFGS
162            Failed to obtain acceptable iterate with BFGS step
163            Attempt to use the scaled gradient direction */
164 
165         if (bnk->f != 0.0) {
166           delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
167         } else {
168           delta = 2.0 / (bnk->gnorm*bnk->gnorm);
169         }
170         ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
171         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
172         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
173         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
174         ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
175 
176         bfgsUpdates = 1;
177         ++bnk->sgrad;
178         stepType = BNK_SCALED_GRADIENT;
179         break;
180 
181       case BNK_SCALED_GRADIENT:
182         /* Can only enter if pc_type == BNK_PC_BFGS
183            The scaled gradient step did not produce a new iterate;
184            attemp to use the gradient direction.
185            Need to make sure we are not using a different diagonal scaling */
186 
187         ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr);
188         ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr);
189         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
190         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
191         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
192         ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
193 
194         bfgsUpdates = 1;
195         ++bnk->grad;
196         stepType = BNK_GRADIENT;
197         break;
198       }
199       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
200 
201       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr);
202       ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
203       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
204     }
205 
206     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
207       /* Failed to find an improving point */
208       bnk->f = bnk->fold;
209       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
210       ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
211       ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
212       step = 0.0;
213       tao->reason = TAO_DIVERGED_LS_FAILURE;
214       break;
215     }
216 
217     /* update trust radius */
218     ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr);
219     ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
220     prered = -prered;
221     actred = bnk->fold - f_full;
222     ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, stepType, &trustAccept);CHKERRQ(ierr);
223 
224     /*  Check for termination */
225     ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
226     if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
227     ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
228     ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr);
229     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
230   }
231   PetscFunctionReturn(0);
232 }
233 
234 PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao)
235 {
236   TAO_BNK        *bnk;
237   PetscErrorCode ierr;
238 
239   PetscFunctionBegin;
240   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
241   tao->ops->solve = TaoSolve_BNLS;
242 
243   bnk = (TAO_BNK *)tao->data;
244   bnk->update_type = BNK_UPDATE_STEP; /* trust region updates based on line search step length */
245   PetscFunctionReturn(0);
246 }