xref: /petsc/src/tao/unconstrained/impls/lmvm/lmvm.c (revision 4f02bc6a7df4e4b03c783003a8e0899ee35fcb05)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/matrix/lmvmmat.h>
3 #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
4 
5 #define LMVM_BFGS                0
6 #define LMVM_SCALED_GRADIENT     1
7 #define LMVM_GRADIENT            2
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "TaoSolve_LMVM"
11 static PetscErrorCode TaoSolve_LMVM(Tao tao)
12 {
13   TAO_LMVM                     *lmP = (TAO_LMVM *)tao->data;
14   PetscReal                    f, fold, gdx, gnorm;
15   PetscReal                    step = 1.0;
16   PetscReal                    delta;
17   PetscErrorCode               ierr;
18   PetscInt                     stepType;
19   TaoConvergedReason           reason = TAO_CONTINUE_ITERATING;
20   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
21 
22   PetscFunctionBegin;
23 
24   if (tao->XL || tao->XU || tao->ops->computebounds) {
25     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n");CHKERRQ(ierr);
26   }
27 
28   /*  Check convergence criteria */
29   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
30   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
31 
32   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
33 
34   ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
35   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
36 
37   /*  Set initial scaling for the function */
38   if (f != 0.0) {
39     delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
40   } else {
41     delta = 2.0 / (gnorm*gnorm);
42   }
43   ierr = MatLMVMSetDelta(lmP->M,delta);CHKERRQ(ierr);
44 
45   /*  Set counter for gradient/reset steps */
46   lmP->bfgs = 0;
47   lmP->sgrad = 0;
48   lmP->grad = 0;
49 
50   /*  Have not converged; continue with Newton method */
51   while (reason == TAO_CONTINUE_ITERATING) {
52     /*  Compute direction */
53     ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr);
54     ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr);
55     ++lmP->bfgs;
56 
57     /*  Check for success (descent direction) */
58     ierr = VecDot(lmP->D, tao->gradient, &gdx);CHKERRQ(ierr);
59     if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
60       /* Step is not descent or direction produced not a number
61          We can assert bfgsUpdates > 1 in this case because
62          the first solve produces the scaled gradient direction,
63          which is guaranteed to be descent
64 
65          Use steepest descent direction (scaled)
66       */
67 
68       ++lmP->grad;
69 
70       if (f != 0.0) {
71         delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
72       } else {
73         delta = 2.0 / (gnorm*gnorm);
74       }
75       ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr);
76       ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr);
77       ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
78       ierr = MatLMVMSolve(lmP->M,tao->gradient, lmP->D);CHKERRQ(ierr);
79 
80       /* On a reset, the direction cannot be not a number; it is a
81          scaled gradient step.  No need to check for this condition. */
82 
83       lmP->bfgs = 1;
84       ++lmP->sgrad;
85       stepType = LMVM_SCALED_GRADIENT;
86     } else {
87       if (1 == lmP->bfgs) {
88         /*  The first BFGS direction is always the scaled gradient */
89         ++lmP->sgrad;
90         stepType = LMVM_SCALED_GRADIENT;
91       } else {
92         ++lmP->bfgs;
93         stepType = LMVM_BFGS;
94       }
95     }
96     ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr);
97 
98     /*  Perform the linesearch */
99     fold = f;
100     ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr);
101     ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr);
102 
103     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status);CHKERRQ(ierr);
104     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
105 
106     while (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_GRADIENT)) {
107       /*  Linesearch failed */
108       /*  Reset factors and use scaled gradient step */
109       f = fold;
110       ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr);
111       ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr);
112 
113       switch(stepType) {
114       case LMVM_BFGS:
115         /*  Failed to obtain acceptable iterate with BFGS step */
116         /*  Attempt to use the scaled gradient direction */
117 
118         if (f != 0.0) {
119           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
120         } else {
121           delta = 2.0 / (gnorm*gnorm);
122         }
123         ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr);
124         ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr);
125         ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
126         ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr);
127 
128         /* On a reset, the direction cannot be not a number; it is a
129            scaled gradient step.  No need to check for this condition. */
130 
131         lmP->bfgs = 1;
132         ++lmP->sgrad;
133         stepType = LMVM_SCALED_GRADIENT;
134         break;
135 
136       case LMVM_SCALED_GRADIENT:
137         /* The scaled gradient step did not produce a new iterate;
138            attempt to use the gradient direction.
139            Need to make sure we are not using a different diagonal scaling */
140         ierr = MatLMVMSetDelta(lmP->M, 1.0);CHKERRQ(ierr);
141         ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr);
142         ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
143         ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr);
144 
145         lmP->bfgs = 1;
146         ++lmP->grad;
147         stepType = LMVM_GRADIENT;
148         break;
149       }
150       ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr);
151 
152       /*  Perform the linesearch */
153       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status);CHKERRQ(ierr);
154       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
155     }
156 
157     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
158       /*  Failed to find an improving point */
159       f = fold;
160       ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr);
161       ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr);
162       step = 0.0;
163       reason = TAO_DIVERGED_LS_FAILURE;
164       tao->reason = TAO_DIVERGED_LS_FAILURE;
165     }
166 
167     /*  Check for termination */
168     ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
169 
170     tao->niter++;
171     ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step,&reason);CHKERRQ(ierr);
172   }
173   PetscFunctionReturn(0);
174 }
175 
176 #undef __FUNCT__
177 #define __FUNCT__ "TaoSetUp_LMVM"
178 static PetscErrorCode TaoSetUp_LMVM(Tao tao)
179 {
180   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
181   PetscInt       n,N;
182   PetscErrorCode ierr;
183   KSP            H0ksp;
184 
185   PetscFunctionBegin;
186   /* Existence of tao->solution checked in TaoSetUp() */
187   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);  }
188   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);  }
189   if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr);  }
190   if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr);  }
191   if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr);  }
192 
193   /*  Create matrix for the limited memory approximation */
194   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
195   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
196   ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&lmP->M);CHKERRQ(ierr);
197   ierr = MatLMVMAllocateVectors(lmP->M,tao->solution);CHKERRQ(ierr);
198 
199   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
200   if (lmP->H0) {
201     const char *prefix;
202     PC H0pc;
203 
204     ierr = MatLMVMSetH0(lmP->M, lmP->H0);CHKERRQ(ierr);
205     ierr = MatLMVMGetH0KSP(lmP->M, &H0ksp);CHKERRQ(ierr);
206 
207     ierr = TaoGetOptionsPrefix(tao, &prefix);CHKERRQ(ierr);
208     ierr = KSPSetOptionsPrefix(H0ksp, prefix);CHKERRQ(ierr);
209     ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0ksp, "tao_h0_");CHKERRQ(ierr);
210     ierr = KSPGetPC(H0ksp, &H0pc);CHKERRQ(ierr);
211     ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0pc,  "tao_h0_");CHKERRQ(ierr);
212 
213     ierr = KSPSetFromOptions(H0ksp);CHKERRQ(ierr);
214     ierr = KSPSetUp(H0ksp);CHKERRQ(ierr);
215   }
216 
217   PetscFunctionReturn(0);
218 }
219 
220 /* ---------------------------------------------------------- */
221 #undef __FUNCT__
222 #define __FUNCT__ "TaoDestroy_LMVM"
223 static PetscErrorCode TaoDestroy_LMVM(Tao tao)
224 {
225   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
226   PetscErrorCode ierr;
227 
228   PetscFunctionBegin;
229   if (tao->setupcalled) {
230     ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr);
231     ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr);
232     ierr = VecDestroy(&lmP->D);CHKERRQ(ierr);
233     ierr = MatDestroy(&lmP->M);CHKERRQ(ierr);
234   }
235 
236   if (lmP->H0) {
237     ierr = PetscObjectDereference((PetscObject)lmP->H0);CHKERRQ(ierr);
238   }
239 
240   ierr = PetscFree(tao->data);CHKERRQ(ierr);
241 
242   PetscFunctionReturn(0);
243 }
244 
245 /*------------------------------------------------------------*/
246 #undef __FUNCT__
247 #define __FUNCT__ "TaoSetFromOptions_LMVM"
248 static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptionItems *PetscOptionsObject,Tao tao)
249 {
250   PetscErrorCode ierr;
251 
252   PetscFunctionBegin;
253   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");CHKERRQ(ierr);
254   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
255   ierr = PetscOptionsTail();CHKERRQ(ierr);
256   PetscFunctionReturn(0);
257 }
258 
259 /*------------------------------------------------------------*/
260 #undef __FUNCT__
261 #define __FUNCT__ "TaoView_LMVM"
262 static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer)
263 {
264   TAO_LMVM       *lm = (TAO_LMVM *)tao->data;
265   PetscBool      isascii;
266   PetscErrorCode ierr;
267 
268   PetscFunctionBegin;
269   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
270   if (isascii) {
271     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
272     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", lm->bfgs);CHKERRQ(ierr);
273     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", lm->sgrad);CHKERRQ(ierr);
274     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lm->grad);CHKERRQ(ierr);
275     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
276   }
277   PetscFunctionReturn(0);
278 }
279 
280 /* ---------------------------------------------------------- */
281 
282 /*MC
283      TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton
284      optimization solver for unconstrained minimization. It solves
285      the Newton step
286               Hkdk = - gk
287 
288      using an approximation Bk in place of Hk, where Bk is composed using
289      the BFGS update formula. A More-Thuente line search is then used
290      to computed the steplength in the dk direction
291   Options Database Keys:
292 +     -tao_lmm_vectors - number of vectors to use for approximation
293 .     -tao_lmm_scale_type - "none","scalar","broyden"
294 .     -tao_lmm_limit_type - "none","average","relative","absolute"
295 .     -tao_lmm_rescale_type - "none","scalar","gl"
296 .     -tao_lmm_limit_mu - mu limiting factor
297 .     -tao_lmm_limit_nu - nu limiting factor
298 .     -tao_lmm_delta_min - minimum delta value
299 .     -tao_lmm_delta_max - maximum delta value
300 .     -tao_lmm_broyden_phi - phi factor for Broyden scaling
301 .     -tao_lmm_scalar_alpha - alpha factor for scalar scaling
302 .     -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal
303 .     -tao_lmm_rescale_beta - beta factor for rescaling diagonal
304 .     -tao_lmm_scalar_history - amount of history for scalar scaling
305 .     -tao_lmm_rescale_history - amount of history for rescaling diagonal
306 -     -tao_lmm_eps - rejection tolerance
307 
308   Level: beginner
309 M*/
310 
311 #undef __FUNCT__
312 #define __FUNCT__ "TaoCreate_LMVM"
313 PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao)
314 {
315   TAO_LMVM       *lmP;
316   const char     *morethuente_type = TAOLINESEARCHMT;
317   PetscErrorCode ierr;
318 
319   PetscFunctionBegin;
320   tao->ops->setup = TaoSetUp_LMVM;
321   tao->ops->solve = TaoSolve_LMVM;
322   tao->ops->view = TaoView_LMVM;
323   tao->ops->setfromoptions = TaoSetFromOptions_LMVM;
324   tao->ops->destroy = TaoDestroy_LMVM;
325 
326   ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr);
327   lmP->D = 0;
328   lmP->M = 0;
329   lmP->Xold = 0;
330   lmP->Gold = 0;
331   lmP->H0   = NULL;
332 
333   tao->data = (void*)lmP;
334   /* Override default settings (unless already changed) */
335   if (!tao->max_it_changed) tao->max_it = 2000;
336   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
337 
338   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
339   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
340   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
341   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
342   PetscFunctionReturn(0);
343 }
344