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