xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision 40cbb1a031ea8f2be4fe2b92dc842b003ad37be3)
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) {
138     PetscCall(MatLMVMSetJ0(blmP->M, blmP->H0));
139   }
140   PetscFunctionReturn(0);
141 }
142 
143 /* ---------------------------------------------------------- */
144 static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
145 {
146   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
147 
148   PetscFunctionBegin;
149   if (tao->setupcalled) {
150     PetscCall(VecDestroy(&blmP->unprojected_gradient));
151     PetscCall(VecDestroy(&blmP->Xold));
152     PetscCall(VecDestroy(&blmP->Gold));
153   }
154   PetscCall(MatDestroy(&blmP->M));
155   if (blmP->H0) {
156     PetscObjectDereference((PetscObject)blmP->H0);
157   }
158   PetscCall(PetscFree(tao->data));
159   PetscFunctionReturn(0);
160 }
161 
162 /*------------------------------------------------------------*/
163 static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
164 {
165   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
166   PetscBool      is_spd;
167 
168   PetscFunctionBegin;
169   PetscOptionsHeadBegin(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");
170   PetscCall(PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL));
171   PetscOptionsHeadEnd();
172   PetscCall(MatSetOptionsPrefix(blmP->M, ((PetscObject)tao)->prefix));
173   PetscCall(MatAppendOptionsPrefix(blmP->M, "tao_blmvm_"));
174   PetscCall(MatSetFromOptions(blmP->M));
175   PetscCall(MatGetOption(blmP->M, MAT_SPD, &is_spd));
176   PetscCheck(is_spd,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
177   PetscFunctionReturn(0);
178 }
179 
180 /*------------------------------------------------------------*/
181 static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer)
182 {
183   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
184   PetscBool      isascii;
185 
186   PetscFunctionBegin;
187   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
188   if (isascii) {
189     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", lmP->grad));
190     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
191     PetscCall(MatView(lmP->M, viewer));
192     PetscCall(PetscViewerPopFormat(viewer));
193   }
194   PetscFunctionReturn(0);
195 }
196 
197 static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
198 {
199   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
200 
201   PetscFunctionBegin;
202   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
203   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
204   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
205   PetscCheck(tao->gradient && blm->unprojected_gradient,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.");
206 
207   PetscCall(VecCopy(tao->gradient,DXL));
208   PetscCall(VecAXPY(DXL,-1.0,blm->unprojected_gradient));
209   PetscCall(VecSet(DXU,0.0));
210   PetscCall(VecPointwiseMax(DXL,DXL,DXU));
211 
212   PetscCall(VecCopy(blm->unprojected_gradient,DXU));
213   PetscCall(VecAXPY(DXU,-1.0,tao->gradient));
214   PetscCall(VecAXPY(DXU,1.0,DXL));
215   PetscFunctionReturn(0);
216 }
217 
218 /* ---------------------------------------------------------- */
219 /*MC
220   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
221          for nonlinear minimization with bound constraints. It is an extension
222          of TAOLMVM
223 
224   Options Database Keys:
225 .     -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls
226 
227   Level: beginner
228 M*/
229 PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
230 {
231   TAO_BLMVM      *blmP;
232   const char     *morethuente_type = TAOLINESEARCHMT;
233 
234   PetscFunctionBegin;
235   tao->ops->setup = TaoSetup_BLMVM;
236   tao->ops->solve = TaoSolve_BLMVM;
237   tao->ops->view = TaoView_BLMVM;
238   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
239   tao->ops->destroy = TaoDestroy_BLMVM;
240   tao->ops->computedual = TaoComputeDual_BLMVM;
241 
242   PetscCall(PetscNewLog(tao,&blmP));
243   blmP->H0 = NULL;
244   blmP->recycle = PETSC_FALSE;
245   tao->data = (void*)blmP;
246 
247   /* Override default settings (unless already changed) */
248   if (!tao->max_it_changed) tao->max_it = 2000;
249   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
250 
251   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
252   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
253   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
254   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao));
255 
256   PetscCall(KSPInitializePackage());
257   PetscCall(MatCreate(((PetscObject)tao)->comm, &blmP->M));
258   PetscCall(MatSetType(blmP->M, MATLMVMBFGS));
259   PetscCall(PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1));
260   PetscFunctionReturn(0);
261 }
262 
263 /*@
264   TaoLMVMRecycle - Enable/disable recycling of the QN history between subsequent TaoSolve calls.
265 
266   Input Parameters:
267 +  tao  - the Tao solver context
268 -  flg - Boolean flag for recycling (PETSC_TRUE or PETSC_FALSE)
269 
270   Level: intermediate
271 @*/
272 PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg)
273 {
274   TAO_LMVM       *lmP;
275   TAO_BLMVM      *blmP;
276   PetscBool      is_lmvm, is_blmvm;
277 
278   PetscFunctionBegin;
279   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm));
280   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm));
281   if (is_lmvm) {
282     lmP = (TAO_LMVM *)tao->data;
283     lmP->recycle = flg;
284   } else if (is_blmvm) {
285     blmP = (TAO_BLMVM *)tao->data;
286     blmP->recycle = flg;
287   }
288   PetscFunctionReturn(0);
289 }
290 
291 /*@
292   TaoLMVMSetH0 - Set the initial Hessian for the QN approximation
293 
294   Input Parameters:
295 +  tao  - the Tao solver context
296 -  H0 - Mat object for the initial Hessian
297 
298   Level: advanced
299 
300 .seealso: `TaoLMVMGetH0()`, `TaoLMVMGetH0KSP()`
301 @*/
302 PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0)
303 {
304   TAO_LMVM       *lmP;
305   TAO_BLMVM      *blmP;
306   PetscBool      is_lmvm, is_blmvm;
307 
308   PetscFunctionBegin;
309   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm));
310   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm));
311   if (is_lmvm) {
312     lmP = (TAO_LMVM *)tao->data;
313     PetscCall(PetscObjectReference((PetscObject)H0));
314     lmP->H0 = H0;
315   } else if (is_blmvm) {
316     blmP = (TAO_BLMVM *)tao->data;
317     PetscCall(PetscObjectReference((PetscObject)H0));
318     blmP->H0 = H0;
319   }
320   PetscFunctionReturn(0);
321 }
322 
323 /*@
324   TaoLMVMGetH0 - Get the matrix object for the QN initial Hessian
325 
326   Input Parameters:
327 .  tao  - the Tao solver context
328 
329   Output Parameters:
330 .  H0 - Mat object for the initial Hessian
331 
332   Level: advanced
333 
334 .seealso: `TaoLMVMSetH0()`, `TaoLMVMGetH0KSP()`
335 @*/
336 PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0)
337 {
338   TAO_LMVM       *lmP;
339   TAO_BLMVM      *blmP;
340   PetscBool      is_lmvm, is_blmvm;
341   Mat            M;
342 
343   PetscFunctionBegin;
344   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm));
345   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm));
346   if (is_lmvm) {
347     lmP = (TAO_LMVM *)tao->data;
348     M = lmP->M;
349   } else if (is_blmvm) {
350     blmP = (TAO_BLMVM *)tao->data;
351     M = blmP->M;
352   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM.");
353   PetscCall(MatLMVMGetJ0(M, H0));
354   PetscFunctionReturn(0);
355 }
356 
357 /*@
358   TaoLMVMGetH0KSP - Get the iterative solver for applying the inverse of the QN initial Hessian
359 
360   Input Parameters:
361 .  tao  - the Tao solver context
362 
363   Output Parameters:
364 .  ksp - KSP solver context for the initial Hessian
365 
366   Level: advanced
367 
368 .seealso: `TaoLMVMGetH0()`, `TaoLMVMGetH0KSP()`
369 @*/
370 PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp)
371 {
372   TAO_LMVM       *lmP;
373   TAO_BLMVM      *blmP;
374   PetscBool      is_lmvm, is_blmvm;
375   Mat            M;
376 
377   PetscFunctionBegin;
378   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm));
379   PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm));
380   if (is_lmvm) {
381     lmP = (TAO_LMVM *)tao->data;
382     M = lmP->M;
383   } else if (is_blmvm) {
384     blmP = (TAO_BLMVM *)tao->data;
385     M = blmP->M;
386   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM.");
387   PetscCall(MatLMVMGetJ0KSP(M, ksp));
388   PetscFunctionReturn(0);
389 }
390