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