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