xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision 2b8d69ca7ea5fe9190df62c1dce3bbd66fce84dd)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/matrix/lmvmmat.h>
3 #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
4 #include <../src/tao/bound/impls/blmvm/blmvm.h>
5 
6 /*------------------------------------------------------------*/
7 #undef __FUNCT__
8 #define __FUNCT__ "TaoSolve_BLMVM"
9 static PetscErrorCode TaoSolve_BLMVM(Tao tao)
10 {
11   PetscErrorCode               ierr;
12   TAO_BLMVM                    *blmP = (TAO_BLMVM *)tao->data;
13   TaoConvergedReason           reason = TAO_CONTINUE_ITERATING;
14   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
15   PetscReal                    f, fold, gdx, gnorm;
16   PetscReal                    stepsize = 1.0,delta;
17 
18   PetscFunctionBegin;
19   /*  Project initial point onto bounds */
20   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
21   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
22   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
23 
24 
25   /* Check convergence criteria */
26   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr);
27   ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
28 
29   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
30   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf pr NaN");
31 
32   ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr);
33   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
34 
35   /* Set initial scaling for the function */
36   if (f != 0.0) {
37     delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
38   } else {
39     delta = 2.0 / (gnorm*gnorm);
40   }
41   ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
42 
43   /* Set counter for gradient/reset steps */
44   blmP->grad = 0;
45   blmP->reset = 0;
46 
47   /* Have not converged; continue with Newton method */
48   while (reason == TAO_CONTINUE_ITERATING) {
49     /* Compute direction */
50     ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
51     ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
52     ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
53 
54     /* Check for success (descent direction) */
55     ierr = VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);CHKERRQ(ierr);
56     if (gdx <= 0) {
57       /* Step is not descent or solve was not successful
58          Use steepest descent direction (scaled) */
59       ++blmP->grad;
60 
61       if (f != 0.0) {
62         delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
63       } else {
64         delta = 2.0 / (gnorm*gnorm);
65       }
66       ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
67       ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr);
68       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
69       ierr = MatLMVMSolve(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       if (f != 0.0) {
91         delta = 2.0* PetscAbsScalar(f) / (gnorm*gnorm);
92       } else {
93         delta = 2.0/ (gnorm*gnorm);
94       }
95       ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
96       ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr);
97       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
98       ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
99       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
100 
101       /* This may be incorrect; linesearch has values fo stepmax and stepmin
102          that should be reset. */
103       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
104       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);CHKERRQ(ierr);
105       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
106 
107       if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
108         tao->reason = TAO_DIVERGED_LS_FAILURE;
109         break;
110       }
111     }
112 
113     /* Check for converged */
114     ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
115     ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
116 
117 
118     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
119     tao->niter++;
120     ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr);
121   }
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "TaoSetup_BLMVM"
127 static PetscErrorCode TaoSetup_BLMVM(Tao tao)
128 {
129   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
130   PetscInt       n,N;
131   PetscErrorCode ierr;
132   KSP            H0ksp;
133 
134   PetscFunctionBegin;
135   /* Existence of tao->solution checked in TaoSetup() */
136   ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr);
137   ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr);
138   ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr);
139 
140   if (!tao->stepdirection) {
141     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
142   }
143   if (!tao->gradient) {
144     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
145   }
146   if (!tao->XL) {
147     ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr);
148     ierr = VecSet(tao->XL,PETSC_NINFINITY);CHKERRQ(ierr);
149   }
150   if (!tao->XU) {
151     ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr);
152     ierr = VecSet(tao->XU,PETSC_INFINITY);CHKERRQ(ierr);
153   }
154   /* Create matrix for the limited memory approximation */
155   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
156   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
157   ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&blmP->M);CHKERRQ(ierr);
158   ierr = MatLMVMAllocateVectors(blmP->M,tao->solution);CHKERRQ(ierr);
159 
160   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
161   if (blmP->H0) {
162     const char *prefix;
163     PC H0pc;
164 
165     ierr = MatLMVMSetH0(blmP->M, blmP->H0);CHKERRQ(ierr);
166     ierr = MatLMVMGetH0KSP(blmP->M, &H0ksp);CHKERRQ(ierr);
167 
168     ierr = TaoGetOptionsPrefix(tao, &prefix);CHKERRQ(ierr);
169     ierr = KSPSetOptionsPrefix(H0ksp, prefix);CHKERRQ(ierr);
170     ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0ksp, "tao_h0_");CHKERRQ(ierr);
171     ierr = KSPGetPC(H0ksp, &H0pc);CHKERRQ(ierr);
172     ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0pc,  "tao_h0_");CHKERRQ(ierr);
173 
174     ierr = KSPSetFromOptions(H0ksp);CHKERRQ(ierr);
175     ierr = KSPSetUp(H0ksp);CHKERRQ(ierr);
176   }
177 
178   PetscFunctionReturn(0);
179 }
180 
181 /* ---------------------------------------------------------- */
182 #undef __FUNCT__
183 #define __FUNCT__ "TaoDestroy_BLMVM"
184 static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
185 {
186   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
187   PetscErrorCode ierr;
188 
189   PetscFunctionBegin;
190   if (tao->setupcalled) {
191     ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
192     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
193     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
194     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
195   }
196 
197   if (blmP->H0) {
198     PetscObjectDereference((PetscObject)blmP->H0);
199   }
200 
201   ierr = PetscFree(tao->data);CHKERRQ(ierr);
202   PetscFunctionReturn(0);
203 }
204 
205 /*------------------------------------------------------------*/
206 #undef __FUNCT__
207 #define __FUNCT__ "TaoSetFromOptions_BLMVM"
208 static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
209 {
210   PetscErrorCode ierr;
211 
212   PetscFunctionBegin;
213   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
214   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
215   ierr = PetscOptionsTail();CHKERRQ(ierr);
216   PetscFunctionReturn(0);
217 }
218 
219 
220 /*------------------------------------------------------------*/
221 #undef __FUNCT__
222 #define __FUNCT__ "TaoView_BLMVM"
223 static int TaoView_BLMVM(Tao tao, PetscViewer viewer)
224 {
225   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
226   PetscBool      isascii;
227   PetscErrorCode ierr;
228 
229   PetscFunctionBegin;
230   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
231   if (isascii) {
232     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
233     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
234     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
235   }
236   PetscFunctionReturn(0);
237 }
238 
239 #undef __FUNCT__
240 #define __FUNCT__ "TaoComputeDual_BLMVM"
241 static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
242 {
243   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
244   PetscErrorCode ierr;
245 
246   PetscFunctionBegin;
247   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
248   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
249   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
250   if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
251 
252   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
253   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
254   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
255   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
256 
257   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
258   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
259   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
260   PetscFunctionReturn(0);
261 }
262 
263 /* ---------------------------------------------------------- */
264 /*MC
265   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
266          for nonlinear minimization with bound constraints. It is an extension
267          of TAOLMVM
268 
269   Options Database Keys:
270 +     -tao_lmm_vectors - number of vectors to use for approximation
271 .     -tao_lmm_scale_type - "none","scalar","broyden"
272 .     -tao_lmm_limit_type - "none","average","relative","absolute"
273 .     -tao_lmm_rescale_type - "none","scalar","gl"
274 .     -tao_lmm_limit_mu - mu limiting factor
275 .     -tao_lmm_limit_nu - nu limiting factor
276 .     -tao_lmm_delta_min - minimum delta value
277 .     -tao_lmm_delta_max - maximum delta value
278 .     -tao_lmm_broyden_phi - phi factor for Broyden scaling
279 .     -tao_lmm_scalar_alpha - alpha factor for scalar scaling
280 .     -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal
281 .     -tao_lmm_rescale_beta - beta factor for rescaling diagonal
282 .     -tao_lmm_scalar_history - amount of history for scalar scaling
283 .     -tao_lmm_rescale_history - amount of history for rescaling diagonal
284 -     -tao_lmm_eps - rejection tolerance
285 
286   Level: beginner
287 M*/
288 #undef __FUNCT__
289 #define __FUNCT__ "TaoCreate_BLMVM"
290 PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
291 {
292   TAO_BLMVM      *blmP;
293   const char     *morethuente_type = TAOLINESEARCHMT;
294   PetscErrorCode ierr;
295 
296   PetscFunctionBegin;
297   tao->ops->setup = TaoSetup_BLMVM;
298   tao->ops->solve = TaoSolve_BLMVM;
299   tao->ops->view = TaoView_BLMVM;
300   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
301   tao->ops->destroy = TaoDestroy_BLMVM;
302   tao->ops->computedual = TaoComputeDual_BLMVM;
303 
304   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
305   blmP->H0 = NULL;
306   tao->data = (void*)blmP;
307 
308   /* Override default settings (unless already changed) */
309   if (!tao->max_it_changed) tao->max_it = 2000;
310   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
311 
312   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
313   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
314   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
315   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
316   PetscFunctionReturn(0);
317 }
318 
319 #undef __FUNCT__
320 #define __FUNCT__ "TaoLMVMSetH0"
321 PETSC_EXTERN PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0)
322 {
323   TAO_LMVM       *lmP;
324   TAO_BLMVM      *blmP;
325   const TaoType  type;
326   PetscBool is_lmvm, is_blmvm;
327 
328   PetscErrorCode ierr;
329 
330   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
331   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
332   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
333 
334   if (is_lmvm) {
335     lmP = (TAO_LMVM *)tao->data;
336     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
337     lmP->H0 = H0;
338   } else if (is_blmvm) {
339     blmP = (TAO_BLMVM *)tao->data;
340     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
341     blmP->H0 = H0;
342   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
343 
344   PetscFunctionReturn(0);
345 }
346 
347 #undef __FUNCT__
348 #define __FUNCT__ "TaoLMVMGetH0"
349 PETSC_EXTERN PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0)
350 {
351   TAO_LMVM       *lmP;
352   TAO_BLMVM      *blmP;
353   const TaoType  type;
354   PetscBool      is_lmvm, is_blmvm;
355   Mat            M;
356 
357   PetscErrorCode ierr;
358 
359   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
360   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
361   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
362 
363   if (is_lmvm) {
364     lmP = (TAO_LMVM *)tao->data;
365     M = lmP->M;
366   } else if (is_blmvm) {
367     blmP = (TAO_BLMVM *)tao->data;
368     M = blmP->M;
369   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
370 
371   ierr = MatLMVMGetH0(M, H0);CHKERRQ(ierr);
372   PetscFunctionReturn(0);
373 }
374 
375 #undef __FUNCT__
376 #define __FUNCT__ "TaoLMVMGetH0KSP"
377 PETSC_EXTERN PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp)
378 {
379   TAO_LMVM       *lmP;
380   TAO_BLMVM      *blmP;
381   const TaoType  type;
382   PetscBool      is_lmvm, is_blmvm;
383   Mat            M;
384   PetscErrorCode ierr;
385 
386   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
387   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
388   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
389 
390   if (is_lmvm) {
391     lmP = (TAO_LMVM *)tao->data;
392     M = lmP->M;
393   } else if (is_blmvm) {
394     blmP = (TAO_BLMVM *)tao->data;
395     M = blmP->M;
396   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
397 
398   ierr = MatLMVMGetH0KSP(M, ksp);CHKERRQ(ierr);
399   PetscFunctionReturn(0);
400 }
401