xref: /petsc/src/tao/unconstrained/impls/lmvm/lmvm.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
3 
4 #define LMVM_STEP_BFGS     0
5 #define LMVM_STEP_GRAD     1
6 
7 static PetscErrorCode TaoSolve_LMVM(Tao tao)
8 {
9   TAO_LMVM                     *lmP = (TAO_LMVM *)tao->data;
10   PetscReal                    f, fold, gdx, gnorm;
11   PetscReal                    step = 1.0;
12   PetscInt                     stepType = LMVM_STEP_GRAD, nupdates;
13   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
14 
15   PetscFunctionBegin;
16 
17   if (tao->XL || tao->XU || tao->ops->computebounds) {
18     PetscCall(PetscInfo(tao,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n"));
19   }
20 
21   /*  Check convergence criteria */
22   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
23   PetscCall(TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm));
24 
25   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
26 
27   tao->reason = TAO_CONTINUE_ITERATING;
28   PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
29   PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step));
30   PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
31   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
32 
33   /*  Set counter for gradient/reset steps */
34   if (!lmP->recycle) {
35     lmP->bfgs = 0;
36     lmP->grad = 0;
37     PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
38   }
39 
40   /*  Have not converged; continue with Newton method */
41   while (tao->reason == TAO_CONTINUE_ITERATING) {
42     /* Call general purpose update function */
43     if (tao->ops->update) PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update));
44 
45     /*  Compute direction */
46     if (lmP->H0) {
47       PetscCall(MatLMVMSetJ0(lmP->M, lmP->H0));
48       stepType = LMVM_STEP_BFGS;
49     }
50     PetscCall(MatLMVMUpdate(lmP->M,tao->solution,tao->gradient));
51     PetscCall(MatSolve(lmP->M, tao->gradient, lmP->D));
52     PetscCall(MatLMVMGetUpdateCount(lmP->M, &nupdates));
53     if (nupdates > 0) stepType = LMVM_STEP_BFGS;
54 
55     /*  Check for success (descent direction) */
56     PetscCall(VecDot(lmP->D, tao->gradient, &gdx));
57     if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
58       /* Step is not descent or direction produced not a number
59          We can assert bfgsUpdates > 1 in this case because
60          the first solve produces the scaled gradient direction,
61          which is guaranteed to be descent
62 
63          Use steepest descent direction (scaled)
64       */
65 
66       PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
67       PetscCall(MatLMVMClearJ0(lmP->M));
68       PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
69       PetscCall(MatSolve(lmP->M,tao->gradient, lmP->D));
70 
71       /* On a reset, the direction cannot be not a number; it is a
72          scaled gradient step.  No need to check for this condition. */
73       stepType = LMVM_STEP_GRAD;
74     }
75     PetscCall(VecScale(lmP->D, -1.0));
76 
77     /*  Perform the linesearch */
78     fold = f;
79     PetscCall(VecCopy(tao->solution, lmP->Xold));
80     PetscCall(VecCopy(tao->gradient, lmP->Gold));
81 
82     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status));
83     PetscCall(TaoAddLineSearchCounts(tao));
84 
85     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_STEP_GRAD)) {
86       /*  Reset factors and use scaled gradient step */
87       f = fold;
88       PetscCall(VecCopy(lmP->Xold, tao->solution));
89       PetscCall(VecCopy(lmP->Gold, tao->gradient));
90 
91       /*  Failed to obtain acceptable iterate with BFGS step */
92       /*  Attempt to use the scaled gradient direction */
93 
94       PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
95       PetscCall(MatLMVMClearJ0(lmP->M));
96       PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
97       PetscCall(MatSolve(lmP->M, tao->solution, tao->gradient));
98 
99       /* On a reset, the direction cannot be not a number; it is a
100           scaled gradient step.  No need to check for this condition. */
101       stepType = LMVM_STEP_GRAD;
102       PetscCall(VecScale(lmP->D, -1.0));
103 
104       /*  Perform the linesearch */
105       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status));
106       PetscCall(TaoAddLineSearchCounts(tao));
107     }
108 
109     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
110       /*  Failed to find an improving point */
111       f = fold;
112       PetscCall(VecCopy(lmP->Xold, tao->solution));
113       PetscCall(VecCopy(lmP->Gold, tao->gradient));
114       step = 0.0;
115       tao->reason = TAO_DIVERGED_LS_FAILURE;
116     } else {
117       /* LS found valid step, so tally up step type */
118       switch (stepType) {
119       case LMVM_STEP_BFGS:
120         ++lmP->bfgs;
121         break;
122       case LMVM_STEP_GRAD:
123         ++lmP->grad;
124         break;
125       default:
126         break;
127       }
128       /*  Compute new gradient norm */
129       PetscCall(TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm));
130     }
131 
132     /* Check convergence */
133     tao->niter++;
134     PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
135     PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step));
136     PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
137   }
138   PetscFunctionReturn(0);
139 }
140 
141 static PetscErrorCode TaoSetUp_LMVM(Tao tao)
142 {
143   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
144   PetscInt       n,N;
145   PetscBool      is_spd;
146 
147   PetscFunctionBegin;
148   /* Existence of tao->solution checked in TaoSetUp() */
149   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution,&tao->gradient));
150   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution,&tao->stepdirection));
151   if (!lmP->D) PetscCall(VecDuplicate(tao->solution,&lmP->D));
152   if (!lmP->Xold) PetscCall(VecDuplicate(tao->solution,&lmP->Xold));
153   if (!lmP->Gold) PetscCall(VecDuplicate(tao->solution,&lmP->Gold));
154 
155   /*  Create matrix for the limited memory approximation */
156   PetscCall(VecGetLocalSize(tao->solution,&n));
157   PetscCall(VecGetSize(tao->solution,&N));
158   PetscCall(MatSetSizes(lmP->M, n, n, N, N));
159   PetscCall(MatLMVMAllocate(lmP->M,tao->solution,tao->gradient));
160   PetscCall(MatGetOption(lmP->M, MAT_SPD, &is_spd));
161   PetscCheck(is_spd,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix is not symmetric positive-definite.");
162 
163   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
164   if (lmP->H0) PetscCall(MatLMVMSetJ0(lmP->M, lmP->H0));
165 
166   PetscFunctionReturn(0);
167 }
168 
169 /* ---------------------------------------------------------- */
170 static PetscErrorCode TaoDestroy_LMVM(Tao tao)
171 {
172   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
173 
174   PetscFunctionBegin;
175   if (tao->setupcalled) {
176     PetscCall(VecDestroy(&lmP->Xold));
177     PetscCall(VecDestroy(&lmP->Gold));
178     PetscCall(VecDestroy(&lmP->D));
179   }
180   PetscCall(MatDestroy(&lmP->M));
181   if (lmP->H0) PetscCall(PetscObjectDereference((PetscObject)lmP->H0));
182   PetscCall(PetscFree(tao->data));
183 
184   PetscFunctionReturn(0);
185 }
186 
187 /*------------------------------------------------------------*/
188 static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptionItems *PetscOptionsObject,Tao tao)
189 {
190   TAO_LMVM       *lm = (TAO_LMVM *)tao->data;
191 
192   PetscFunctionBegin;
193   PetscOptionsHeadBegin(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");
194   PetscCall(PetscOptionsBool("-tao_lmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",lm->recycle,&lm->recycle,NULL));
195   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
196   PetscCall(MatSetFromOptions(lm->M));
197   PetscOptionsHeadEnd();
198   PetscFunctionReturn(0);
199 }
200 
201 /*------------------------------------------------------------*/
202 static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer)
203 {
204   TAO_LMVM       *lm = (TAO_LMVM *)tao->data;
205   PetscBool      isascii;
206   PetscInt       recycled_its;
207 
208   PetscFunctionBegin;
209   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
210   if (isascii) {
211     PetscCall(PetscViewerASCIIPrintf(viewer, "  Gradient steps: %" PetscInt_FMT "\n", lm->grad));
212     if (lm->recycle) {
213       PetscCall(PetscViewerASCIIPrintf(viewer, "  Recycle: on\n"));
214       recycled_its = lm->bfgs + lm->grad;
215       PetscCall(PetscViewerASCIIPrintf(viewer, "  Total recycled iterations: %" PetscInt_FMT "\n", recycled_its));
216     }
217   }
218   PetscFunctionReturn(0);
219 }
220 
221 /* ---------------------------------------------------------- */
222 
223 /*MC
224   TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton
225   optimization solver for unconstrained minimization. It solves
226   the Newton step
227           Hkdk = - gk
228 
229   using an approximation Bk in place of Hk, where Bk is composed using
230   the BFGS update formula. A More-Thuente line search is then used
231   to computed the steplength in the dk direction
232 
233   Options Database Keys:
234 +   -tao_lmvm_recycle - enable recycling LMVM updates between TaoSolve() calls
235 -   -tao_lmvm_no_scale - (developer) disables diagonal Broyden scaling on the LMVM approximation
236 
237   Level: beginner
238 M*/
239 
240 PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao)
241 {
242   TAO_LMVM       *lmP;
243   const char     *morethuente_type = TAOLINESEARCHMT;
244 
245   PetscFunctionBegin;
246   tao->ops->setup = TaoSetUp_LMVM;
247   tao->ops->solve = TaoSolve_LMVM;
248   tao->ops->view = TaoView_LMVM;
249   tao->ops->setfromoptions = TaoSetFromOptions_LMVM;
250   tao->ops->destroy = TaoDestroy_LMVM;
251 
252   PetscCall(PetscNewLog(tao,&lmP));
253   lmP->D = NULL;
254   lmP->M = NULL;
255   lmP->Xold = NULL;
256   lmP->Gold = NULL;
257   lmP->H0   = NULL;
258   lmP->recycle = PETSC_FALSE;
259 
260   tao->data = (void*)lmP;
261   /* Override default settings (unless already changed) */
262   if (!tao->max_it_changed) tao->max_it = 2000;
263   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
264 
265   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch));
266   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
267   PetscCall(TaoLineSearchSetType(tao->linesearch,morethuente_type));
268   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao));
269   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix));
270 
271   PetscCall(KSPInitializePackage());
272   PetscCall(MatCreate(((PetscObject)tao)->comm, &lmP->M));
273   PetscCall(PetscObjectIncrementTabLevel((PetscObject)lmP->M, (PetscObject)tao, 1));
274   PetscCall(MatSetType(lmP->M, MATLMVMBFGS));
275   PetscCall(MatSetOptionsPrefix(lmP->M, "tao_lmvm_"));
276   PetscFunctionReturn(0);
277 }
278