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