xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision 63f3c55c12ae2f190c391f6df6c540efe018a44a)
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     /* 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 = MatSymBrdnSetDelta(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 
111 
112     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
113     tao->niter++;
114     ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
115     ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);CHKERRQ(ierr);
116     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
117   }
118   PetscFunctionReturn(0);
119 }
120 
121 static PetscErrorCode TaoSetup_BLMVM(Tao tao)
122 {
123   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
124   PetscErrorCode ierr;
125 
126   PetscFunctionBegin;
127   /* Existence of tao->solution checked in TaoSetup() */
128   ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr);
129   ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr);
130   ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr);
131 
132   if (!tao->stepdirection) {
133     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
134   }
135   if (!tao->gradient) {
136     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
137   }
138   if (!tao->XL) {
139     ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr);
140     ierr = VecSet(tao->XL,PETSC_NINFINITY);CHKERRQ(ierr);
141   }
142   if (!tao->XU) {
143     ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr);
144     ierr = VecSet(tao->XU,PETSC_INFINITY);CHKERRQ(ierr);
145   }
146   /* Allocate matrix for the limited memory approximation */
147   ierr = MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient);CHKERRQ(ierr);
148 
149   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
150   if (blmP->H0) {
151     ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr);
152   }
153   PetscFunctionReturn(0);
154 }
155 
156 /* ---------------------------------------------------------- */
157 static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
158 {
159   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
160   PetscErrorCode ierr;
161 
162   PetscFunctionBegin;
163   if (tao->setupcalled) {
164     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
165     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
166     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
167   }
168   ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
169   if (blmP->H0) {
170     PetscObjectDereference((PetscObject)blmP->H0);
171   }
172   ierr = PetscFree(tao->data);CHKERRQ(ierr);
173   PetscFunctionReturn(0);
174 }
175 
176 /*------------------------------------------------------------*/
177 static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
178 {
179   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
180   PetscErrorCode ierr;
181   PetscBool      is_spd;
182 
183   PetscFunctionBegin;
184   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
185   ierr = PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL);CHKERRQ(ierr);
186   ierr = PetscOptionsTail();CHKERRQ(ierr);
187   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
188   ierr = MatSetFromOptions(blmP->M);CHKERRQ(ierr);
189   ierr = MatGetOption(blmP->M, MAT_SPD, &is_spd);CHKERRQ(ierr);
190   if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
191   PetscFunctionReturn(0);
192 }
193 
194 
195 /*------------------------------------------------------------*/
196 static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer)
197 {
198   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
199   PetscBool      isascii;
200   PetscErrorCode ierr;
201 
202   PetscFunctionBegin;
203   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
204   if (isascii) {
205     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
206     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
207     ierr = MatView(lmP->M, viewer);CHKERRQ(ierr);
208     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
209   }
210   PetscFunctionReturn(0);
211 }
212 
213 static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
214 {
215   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
216   PetscErrorCode ierr;
217 
218   PetscFunctionBegin;
219   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
220   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
221   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
222   if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
223 
224   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
225   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
226   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
227   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
228 
229   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
230   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
231   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
232   PetscFunctionReturn(0);
233 }
234 
235 /* ---------------------------------------------------------- */
236 /*MC
237   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
238          for nonlinear minimization with bound constraints. It is an extension
239          of TAOLMVM
240 
241   Options Database Keys:
242 .     -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls
243 
244   Level: beginner
245 M*/
246 PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
247 {
248   TAO_BLMVM      *blmP;
249   const char     *morethuente_type = TAOLINESEARCHMT;
250   PetscErrorCode ierr;
251 
252   PetscFunctionBegin;
253   tao->ops->setup = TaoSetup_BLMVM;
254   tao->ops->solve = TaoSolve_BLMVM;
255   tao->ops->view = TaoView_BLMVM;
256   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
257   tao->ops->destroy = TaoDestroy_BLMVM;
258   tao->ops->computedual = TaoComputeDual_BLMVM;
259 
260   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
261   blmP->H0 = NULL;
262   blmP->recycle = PETSC_FALSE;
263   tao->data = (void*)blmP;
264 
265   /* Override default settings (unless already changed) */
266   if (!tao->max_it_changed) tao->max_it = 2000;
267   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
268 
269   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
270   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
271   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
272   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
273   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
274 
275   ierr = KSPInitializePackage();CHKERRQ(ierr);
276   ierr = MatCreate(((PetscObject)tao)->comm, &blmP->M);CHKERRQ(ierr);
277   ierr = MatSetType(blmP->M, MATLMVMBFGS);CHKERRQ(ierr);
278   ierr = PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1);CHKERRQ(ierr);
279   ierr = MatSetOptionsPrefix(blmP->M, "tao_blmvm_");CHKERRQ(ierr);
280   PetscFunctionReturn(0);
281 }
282 
283 PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg)
284 {
285   TAO_LMVM       *lmP;
286   TAO_BLMVM      *blmP;
287   TaoType        type;
288   PetscBool      is_lmvm, is_blmvm;
289   PetscErrorCode ierr;
290 
291   PetscFunctionBegin;
292   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
293   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
294   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
295 
296   if (is_lmvm) {
297     lmP = (TAO_LMVM *)tao->data;
298     lmP->recycle = flg;
299   } else if (is_blmvm) {
300     blmP = (TAO_BLMVM *)tao->data;
301     blmP->recycle = flg;
302   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
303   PetscFunctionReturn(0);
304 }
305 
306 PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0)
307 {
308   TAO_LMVM       *lmP;
309   TAO_BLMVM      *blmP;
310   TaoType        type;
311   PetscBool      is_lmvm, is_blmvm;
312   PetscErrorCode ierr;
313 
314   PetscFunctionBegin;
315   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
316   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
317   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
318 
319   if (is_lmvm) {
320     lmP = (TAO_LMVM *)tao->data;
321     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
322     lmP->H0 = H0;
323   } else if (is_blmvm) {
324     blmP = (TAO_BLMVM *)tao->data;
325     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
326     blmP->H0 = H0;
327   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
328   PetscFunctionReturn(0);
329 }
330 
331 PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0)
332 {
333   TAO_LMVM       *lmP;
334   TAO_BLMVM      *blmP;
335   TaoType        type;
336   PetscBool      is_lmvm, is_blmvm;
337   Mat            M;
338 
339   PetscErrorCode ierr;
340 
341   PetscFunctionBegin;
342   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
343   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
344   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
345 
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_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
353   ierr = MatLMVMGetJ0(M, H0);CHKERRQ(ierr);
354   PetscFunctionReturn(0);
355 }
356 
357 PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp)
358 {
359   TAO_LMVM       *lmP;
360   TAO_BLMVM      *blmP;
361   TaoType        type;
362   PetscBool      is_lmvm, is_blmvm;
363   Mat            M;
364   PetscErrorCode ierr;
365 
366   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
367   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
368   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
369 
370   if (is_lmvm) {
371     lmP = (TAO_LMVM *)tao->data;
372     M = lmP->M;
373   } else if (is_blmvm) {
374     blmP = (TAO_BLMVM *)tao->data;
375     M = blmP->M;
376   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
377   ierr = MatLMVMGetJ0KSP(M, ksp);CHKERRQ(ierr);
378   PetscFunctionReturn(0);
379 }
380