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