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