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