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