xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision fdf6c4e30aafdbc795e4f76379caa977fd5cdf5a)
1 #include <petsctaolinesearch.h>      /*I "petsctaolinesearch.h" I*/
2 #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
3 #include <../src/tao/bound/impls/blmvm/blmvm.h>
4 
5 /*------------------------------------------------------------*/
6 static PetscErrorCode TaoSolve_BLMVM(Tao tao)
7 {
8   TAO_BLMVM                    *blmP = (TAO_BLMVM *)tao->data;
9   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
10   PetscReal                    f, fold, gdx, gnorm, gnorm2;
11   PetscReal                    stepsize = 1.0,delta;
12 
13   PetscFunctionBegin;
14   /*  Project initial point onto bounds */
15   PetscCall(TaoComputeVariableBounds(tao));
16   PetscCall(VecMedian(tao->XL,tao->solution,tao->XU,tao->solution));
17   PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU));
18 
19   /* Check convergence criteria */
20   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient));
21   PetscCall(VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient));
22 
23   PetscCall(TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm));
24   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
25 
26   tao->reason = TAO_CONTINUE_ITERATING;
27   PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
28   PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize));
29   PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
30   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
31 
32   /* Set counter for gradient/reset steps */
33   if (!blmP->recycle) {
34     blmP->grad = 0;
35     blmP->reset = 0;
36     PetscCall(MatLMVMReset(blmP->M, PETSC_FALSE));
37   }
38 
39   /* Have not converged; continue with Newton method */
40   while (tao->reason == TAO_CONTINUE_ITERATING) {
41     /* Call general purpose update function */
42     if (tao->ops->update) {
43       PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update));
44       PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
45     }
46     /* Compute direction */
47     gnorm2 = gnorm*gnorm;
48     if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON;
49     if (f == 0.0) {
50       delta = 2.0 / gnorm2;
51     } else {
52       delta = 2.0 * PetscAbsScalar(f) / gnorm2;
53     }
54     PetscCall(MatLMVMSymBroydenSetDelta(blmP->M, delta));
55     PetscCall(MatLMVMUpdate(blmP->M, tao->solution, tao->gradient));
56     PetscCall(MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection));
57     PetscCall(VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient));
58 
59     /* Check for success (descent direction) */
60     PetscCall(VecDot(blmP->unprojected_gradient, tao->gradient, &gdx));
61     if (gdx <= 0) {
62       /* Step is not descent or solve was not successful
63          Use steepest descent direction (scaled) */
64       ++blmP->grad;
65 
66       PetscCall(MatLMVMReset(blmP->M, PETSC_FALSE));
67       PetscCall(MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient));
68       PetscCall(MatSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection));
69     }
70     PetscCall(VecScale(tao->stepdirection,-1.0));
71 
72     /* Perform the linesearch */
73     fold = f;
74     PetscCall(VecCopy(tao->solution, blmP->Xold));
75     PetscCall(VecCopy(blmP->unprojected_gradient, blmP->Gold));
76     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0));
77     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status));
78     PetscCall(TaoAddLineSearchCounts(tao));
79 
80     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
81       /* Linesearch failed
82          Reset factors and use scaled (projected) gradient step */
83       ++blmP->reset;
84 
85       f = fold;
86       PetscCall(VecCopy(blmP->Xold, tao->solution));
87       PetscCall(VecCopy(blmP->Gold, blmP->unprojected_gradient));
88 
89       PetscCall(MatLMVMReset(blmP->M, PETSC_FALSE));
90       PetscCall(MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient));
91       PetscCall(MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection));
92       PetscCall(VecScale(tao->stepdirection, -1.0));
93 
94       /* This may be incorrect; linesearch has values for stepmax and stepmin
95          that should be reset. */
96       PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0));
97       PetscCall(TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status));
98       PetscCall(TaoAddLineSearchCounts(tao));
99 
100       if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
101         tao->reason = TAO_DIVERGED_LS_FAILURE;
102         break;
103       }
104     }
105 
106     /* Check for converged */
107     PetscCall(VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient));
108     PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
109     PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Not-a-Number");
110     tao->niter++;
111     PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
112     PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize));
113     PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
114   }
115   PetscFunctionReturn(0);
116 }
117 
118 static PetscErrorCode TaoSetup_BLMVM(Tao tao)
119 {
120   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
121 
122   PetscFunctionBegin;
123   /* Existence of tao->solution checked in TaoSetup() */
124   PetscCall(VecDuplicate(tao->solution,&blmP->Xold));
125   PetscCall(VecDuplicate(tao->solution,&blmP->Gold));
126   PetscCall(VecDuplicate(tao->solution, &blmP->unprojected_gradient));
127   if (!tao->stepdirection) {
128     PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
129   }
130   if (!tao->gradient) {
131     PetscCall(VecDuplicate(tao->solution,&tao->gradient));
132   }
133   /* Allocate matrix for the limited memory approximation */
134   PetscCall(MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient));
135 
136   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
137   if (blmP->H0) PetscCall(MatLMVMSetJ0(blmP->M, blmP->H0));
138   PetscFunctionReturn(0);
139 }
140 
141 /* ---------------------------------------------------------- */
142 static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
143 {
144   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
145 
146   PetscFunctionBegin;
147   if (tao->setupcalled) {
148     PetscCall(VecDestroy(&blmP->unprojected_gradient));
149     PetscCall(VecDestroy(&blmP->Xold));
150     PetscCall(VecDestroy(&blmP->Gold));
151   }
152   PetscCall(MatDestroy(&blmP->M));
153   if (blmP->H0) {
154     PetscObjectDereference((PetscObject)blmP->H0);
155   }
156   PetscCall(PetscFree(tao->data));
157   PetscFunctionReturn(0);
158 }
159 
160 /*------------------------------------------------------------*/
161 static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
162 {
163   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
164   PetscBool      is_spd,is_set;
165 
166   PetscFunctionBegin;
167   PetscOptionsHeadBegin(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");
168   PetscCall(PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL));
169   PetscOptionsHeadEnd();
170   PetscCall(MatSetOptionsPrefix(blmP->M, ((PetscObject)tao)->prefix));
171   PetscCall(MatAppendOptionsPrefix(blmP->M, "tao_blmvm_"));
172   PetscCall(MatSetFromOptions(blmP->M));
173   PetscCall(MatIsSPDKnown(blmP->M, &is_set, &is_spd));
174   PetscCheck(is_set && is_spd,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
175   PetscFunctionReturn(0);
176 }
177 
178 /*------------------------------------------------------------*/
179 static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer)
180 {
181   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
182   PetscBool      isascii;
183 
184   PetscFunctionBegin;
185   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
186   if (isascii) {
187     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", lmP->grad));
188     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
189     PetscCall(MatView(lmP->M, viewer));
190     PetscCall(PetscViewerPopFormat(viewer));
191   }
192   PetscFunctionReturn(0);
193 }
194 
195 static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
196 {
197   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
198 
199   PetscFunctionBegin;
200   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
201   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
202   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
203   PetscCheck(tao->gradient && blm->unprojected_gradient,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.");
204 
205   PetscCall(VecCopy(tao->gradient,DXL));
206   PetscCall(VecAXPY(DXL,-1.0,blm->unprojected_gradient));
207   PetscCall(VecSet(DXU,0.0));
208   PetscCall(VecPointwiseMax(DXL,DXL,DXU));
209 
210   PetscCall(VecCopy(blm->unprojected_gradient,DXU));
211   PetscCall(VecAXPY(DXU,-1.0,tao->gradient));
212   PetscCall(VecAXPY(DXU,1.0,DXL));
213   PetscFunctionReturn(0);
214 }
215 
216 /* ---------------------------------------------------------- */
217 /*MC
218   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
219          for nonlinear minimization with bound constraints. It is an extension
220          of TAOLMVM
221 
222   Options Database Keys:
223 .     -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls
224 
225   Level: beginner
226 M*/
227 PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
228 {
229   TAO_BLMVM      *blmP;
230   const char     *morethuente_type = TAOLINESEARCHMT;
231 
232   PetscFunctionBegin;
233   tao->ops->setup = TaoSetup_BLMVM;
234   tao->ops->solve = TaoSolve_BLMVM;
235   tao->ops->view = TaoView_BLMVM;
236   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
237   tao->ops->destroy = TaoDestroy_BLMVM;
238   tao->ops->computedual = TaoComputeDual_BLMVM;
239 
240   PetscCall(PetscNewLog(tao,&blmP));
241   blmP->H0 = NULL;
242   blmP->recycle = PETSC_FALSE;
243   tao->data = (void*)blmP;
244 
245   /* Override default settings (unless already changed) */
246   if (!tao->max_it_changed) tao->max_it = 2000;
247   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
248 
249   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
250   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
251   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
252   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao));
253 
254   PetscCall(KSPInitializePackage());
255   PetscCall(MatCreate(((PetscObject)tao)->comm, &blmP->M));
256   PetscCall(MatSetType(blmP->M, MATLMVMBFGS));
257   PetscCall(PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1));
258   PetscFunctionReturn(0);
259 }
260 
261 /*@
262   TaoLMVMRecycle - Enable/disable recycling of the QN history between subsequent TaoSolve calls.
263 
264   Input Parameters:
265 +  tao  - the Tao solver context
266 -  flg - Boolean flag for recycling (PETSC_TRUE or PETSC_FALSE)
267 
268   Level: intermediate
269 @*/
270 PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg)
271 {
272   TAO_LMVM       *lmP;
273   TAO_BLMVM      *blmP;
274   PetscBool      is_lmvm, is_blmvm;
275 
276   PetscFunctionBegin;
277   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm));
278   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm));
279   if (is_lmvm) {
280     lmP = (TAO_LMVM *)tao->data;
281     lmP->recycle = flg;
282   } else if (is_blmvm) {
283     blmP = (TAO_BLMVM *)tao->data;
284     blmP->recycle = flg;
285   }
286   PetscFunctionReturn(0);
287 }
288 
289 /*@
290   TaoLMVMSetH0 - Set the initial Hessian for the QN approximation
291 
292   Input Parameters:
293 +  tao  - the Tao solver context
294 -  H0 - Mat object for the initial Hessian
295 
296   Level: advanced
297 
298 .seealso: `TaoLMVMGetH0()`, `TaoLMVMGetH0KSP()`
299 @*/
300 PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0)
301 {
302   TAO_LMVM       *lmP;
303   TAO_BLMVM      *blmP;
304   PetscBool      is_lmvm, is_blmvm;
305 
306   PetscFunctionBegin;
307   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm));
308   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm));
309   if (is_lmvm) {
310     lmP = (TAO_LMVM *)tao->data;
311     PetscCall(PetscObjectReference((PetscObject)H0));
312     lmP->H0 = H0;
313   } else if (is_blmvm) {
314     blmP = (TAO_BLMVM *)tao->data;
315     PetscCall(PetscObjectReference((PetscObject)H0));
316     blmP->H0 = H0;
317   }
318   PetscFunctionReturn(0);
319 }
320 
321 /*@
322   TaoLMVMGetH0 - Get the matrix object for the QN initial Hessian
323 
324   Input Parameters:
325 .  tao  - the Tao solver context
326 
327   Output Parameters:
328 .  H0 - Mat object for the initial Hessian
329 
330   Level: advanced
331 
332 .seealso: `TaoLMVMSetH0()`, `TaoLMVMGetH0KSP()`
333 @*/
334 PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0)
335 {
336   TAO_LMVM       *lmP;
337   TAO_BLMVM      *blmP;
338   PetscBool      is_lmvm, is_blmvm;
339   Mat            M;
340 
341   PetscFunctionBegin;
342   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm));
343   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm));
344   if (is_lmvm) {
345     lmP = (TAO_LMVM *)tao->data;
346     M = lmP->M;
347   } else if (is_blmvm) {
348     blmP = (TAO_BLMVM *)tao->data;
349     M = blmP->M;
350   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM.");
351   PetscCall(MatLMVMGetJ0(M, H0));
352   PetscFunctionReturn(0);
353 }
354 
355 /*@
356   TaoLMVMGetH0KSP - Get the iterative solver for applying the inverse of the QN initial Hessian
357 
358   Input Parameters:
359 .  tao  - the Tao solver context
360 
361   Output Parameters:
362 .  ksp - KSP solver context for the initial Hessian
363 
364   Level: advanced
365 
366 .seealso: `TaoLMVMGetH0()`, `TaoLMVMGetH0KSP()`
367 @*/
368 PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp)
369 {
370   TAO_LMVM       *lmP;
371   TAO_BLMVM      *blmP;
372   PetscBool      is_lmvm, is_blmvm;
373   Mat            M;
374 
375   PetscFunctionBegin;
376   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm));
377   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm));
378   if (is_lmvm) {
379     lmP = (TAO_LMVM *)tao->data;
380     M = lmP->M;
381   } else if (is_blmvm) {
382     blmP = (TAO_BLMVM *)tao->data;
383     M = blmP->M;
384   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM.");
385   PetscCall(MatLMVMGetJ0KSP(M, ksp));
386   PetscFunctionReturn(0);
387 }
388