xref: /petsc/src/tao/unconstrained/impls/lmvm/lmvm.c (revision b94d7ded0a05f1bbd5e48daa6f92b28259c75b44)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
2aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
3a7e14dcfSSatish Balay 
4cd929ea3SAlp Dener #define LMVM_STEP_BFGS     0
5cd929ea3SAlp Dener #define LMVM_STEP_GRAD     1
6a7e14dcfSSatish Balay 
7441846f8SBarry Smith static PetscErrorCode TaoSolve_LMVM(Tao tao)
8a7e14dcfSSatish Balay {
9a7e14dcfSSatish Balay   TAO_LMVM                     *lmP = (TAO_LMVM *)tao->data;
10a7e14dcfSSatish Balay   PetscReal                    f, fold, gdx, gnorm;
11a7e14dcfSSatish Balay   PetscReal                    step = 1.0;
12cd929ea3SAlp Dener   PetscInt                     stepType = LMVM_STEP_GRAD, nupdates;
13e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
14a7e14dcfSSatish Balay 
15a7e14dcfSSatish Balay   PetscFunctionBegin;
16a7e14dcfSSatish Balay 
17a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
189566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n"));
19a7e14dcfSSatish Balay   }
20a7e14dcfSSatish Balay 
21a7e14dcfSSatish Balay   /*  Check convergence criteria */
229566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
239566063dSJacob Faibussowitsch   PetscCall(TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm));
24a9603a14SPatrick Farrell 
253c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
26a7e14dcfSSatish Balay 
273ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
289566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
299566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step));
309566063dSJacob Faibussowitsch   PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
313ecd9318SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
32a7e14dcfSSatish Balay 
33a7e14dcfSSatish Balay   /*  Set counter for gradient/reset steps */
34cd929ea3SAlp Dener   if (!lmP->recycle) {
35a7e14dcfSSatish Balay     lmP->bfgs = 0;
36a7e14dcfSSatish Balay     lmP->grad = 0;
379566063dSJacob Faibussowitsch     PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
38de6ffafeSAlp Dener   }
39a7e14dcfSSatish Balay 
40a7e14dcfSSatish Balay   /*  Have not converged; continue with Newton method */
413ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
42e1e80dc8SAlp Dener     /* Call general purpose update function */
431baa6e33SBarry Smith     if (tao->ops->update) PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update));
44e1e80dc8SAlp Dener 
45a7e14dcfSSatish Balay     /*  Compute direction */
46cd929ea3SAlp Dener     if (lmP->H0) {
479566063dSJacob Faibussowitsch       PetscCall(MatLMVMSetJ0(lmP->M, lmP->H0));
48cd929ea3SAlp Dener       stepType = LMVM_STEP_BFGS;
49cd929ea3SAlp Dener     }
509566063dSJacob Faibussowitsch     PetscCall(MatLMVMUpdate(lmP->M,tao->solution,tao->gradient));
519566063dSJacob Faibussowitsch     PetscCall(MatSolve(lmP->M, tao->gradient, lmP->D));
529566063dSJacob Faibussowitsch     PetscCall(MatLMVMGetUpdateCount(lmP->M, &nupdates));
53cd929ea3SAlp Dener     if (nupdates > 0) stepType = LMVM_STEP_BFGS;
54a7e14dcfSSatish Balay 
55a7e14dcfSSatish Balay     /*  Check for success (descent direction) */
569566063dSJacob Faibussowitsch     PetscCall(VecDot(lmP->D, tao->gradient, &gdx));
57a7e14dcfSSatish Balay     if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
58a7e14dcfSSatish Balay       /* Step is not descent or direction produced not a number
59a7e14dcfSSatish Balay          We can assert bfgsUpdates > 1 in this case because
60a7e14dcfSSatish Balay          the first solve produces the scaled gradient direction,
61a7e14dcfSSatish Balay          which is guaranteed to be descent
62a7e14dcfSSatish Balay 
63a7e14dcfSSatish Balay          Use steepest descent direction (scaled)
64a7e14dcfSSatish Balay       */
65a7e14dcfSSatish Balay 
669566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
679566063dSJacob Faibussowitsch       PetscCall(MatLMVMClearJ0(lmP->M));
689566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
699566063dSJacob Faibussowitsch       PetscCall(MatSolve(lmP->M,tao->gradient, lmP->D));
70a7e14dcfSSatish Balay 
71a7e14dcfSSatish Balay       /* On a reset, the direction cannot be not a number; it is a
72a7e14dcfSSatish Balay          scaled gradient step.  No need to check for this condition. */
73cd929ea3SAlp Dener       stepType = LMVM_STEP_GRAD;
74a7e14dcfSSatish Balay     }
759566063dSJacob Faibussowitsch     PetscCall(VecScale(lmP->D, -1.0));
76a7e14dcfSSatish Balay 
77a7e14dcfSSatish Balay     /*  Perform the linesearch */
78a7e14dcfSSatish Balay     fold = f;
799566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, lmP->Xold));
809566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, lmP->Gold));
81a7e14dcfSSatish Balay 
829566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status));
839566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
84a7e14dcfSSatish Balay 
85cd929ea3SAlp Dener     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_STEP_GRAD)) {
86a7e14dcfSSatish Balay       /*  Reset factors and use scaled gradient step */
87a7e14dcfSSatish Balay       f = fold;
889566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Xold, tao->solution));
899566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Gold, tao->gradient));
90a7e14dcfSSatish Balay 
91a7e14dcfSSatish Balay       /*  Failed to obtain acceptable iterate with BFGS step */
92a7e14dcfSSatish Balay       /*  Attempt to use the scaled gradient direction */
93a7e14dcfSSatish Balay 
949566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
959566063dSJacob Faibussowitsch       PetscCall(MatLMVMClearJ0(lmP->M));
969566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
979566063dSJacob Faibussowitsch       PetscCall(MatSolve(lmP->M, tao->solution, tao->gradient));
98a7e14dcfSSatish Balay 
99a7e14dcfSSatish Balay       /* On a reset, the direction cannot be not a number; it is a
100a7e14dcfSSatish Balay           scaled gradient step.  No need to check for this condition. */
101cd929ea3SAlp Dener       stepType = LMVM_STEP_GRAD;
1029566063dSJacob Faibussowitsch       PetscCall(VecScale(lmP->D, -1.0));
103a7e14dcfSSatish Balay 
104a7e14dcfSSatish Balay       /*  Perform the linesearch */
1059566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status));
1069566063dSJacob Faibussowitsch       PetscCall(TaoAddLineSearchCounts(tao));
107a7e14dcfSSatish Balay     }
108a7e14dcfSSatish Balay 
10987f595a5SBarry Smith     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
110a7e14dcfSSatish Balay       /*  Failed to find an improving point */
111a7e14dcfSSatish Balay       f = fold;
1129566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Xold, tao->solution));
1139566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Gold, tao->gradient));
114a7e14dcfSSatish Balay       step = 0.0;
115a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
116cd929ea3SAlp Dener     } else {
117cd929ea3SAlp Dener       /* LS found valid step, so tally up step type */
118cd929ea3SAlp Dener       switch (stepType) {
119cd929ea3SAlp Dener       case LMVM_STEP_BFGS:
120cd929ea3SAlp Dener         ++lmP->bfgs;
121cd929ea3SAlp Dener         break;
122cd929ea3SAlp Dener       case LMVM_STEP_GRAD:
123cd929ea3SAlp Dener         ++lmP->grad;
124cd929ea3SAlp Dener         break;
125cd929ea3SAlp Dener       default:
126cd929ea3SAlp Dener         break;
127cd929ea3SAlp Dener       }
128cd929ea3SAlp Dener       /*  Compute new gradient norm */
1299566063dSJacob Faibussowitsch       PetscCall(TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm));
130a7e14dcfSSatish Balay     }
131a9603a14SPatrick Farrell 
132cd929ea3SAlp Dener     /* Check convergence */
1338931d482SJason Sarich     tao->niter++;
1349566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
1359566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step));
1369566063dSJacob Faibussowitsch     PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
137a7e14dcfSSatish Balay   }
138a7e14dcfSSatish Balay   PetscFunctionReturn(0);
139a7e14dcfSSatish Balay }
14087f595a5SBarry Smith 
141441846f8SBarry Smith static PetscErrorCode TaoSetUp_LMVM(Tao tao)
142a7e14dcfSSatish Balay {
143a7e14dcfSSatish Balay   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
144a7e14dcfSSatish Balay   PetscInt       n,N;
145*b94d7dedSBarry Smith   PetscBool      is_set,is_spd;
146a7e14dcfSSatish Balay 
147a7e14dcfSSatish Balay   PetscFunctionBegin;
148a7e14dcfSSatish Balay   /* Existence of tao->solution checked in TaoSetUp() */
1499566063dSJacob Faibussowitsch   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution,&tao->gradient));
1509566063dSJacob Faibussowitsch   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution,&tao->stepdirection));
1519566063dSJacob Faibussowitsch   if (!lmP->D) PetscCall(VecDuplicate(tao->solution,&lmP->D));
1529566063dSJacob Faibussowitsch   if (!lmP->Xold) PetscCall(VecDuplicate(tao->solution,&lmP->Xold));
1539566063dSJacob Faibussowitsch   if (!lmP->Gold) PetscCall(VecDuplicate(tao->solution,&lmP->Gold));
154a7e14dcfSSatish Balay 
155a7e14dcfSSatish Balay   /*  Create matrix for the limited memory approximation */
1569566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(tao->solution,&n));
1579566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tao->solution,&N));
1589566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(lmP->M, n, n, N, N));
1599566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(lmP->M,tao->solution,tao->gradient));
160*b94d7dedSBarry Smith   PetscCall(MatIsSPDKnown(lmP->M, &is_set,&is_spd));
161*b94d7dedSBarry Smith   PetscCheck(is_set && is_spd,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix is not symmetric positive-definite.");
162a9603a14SPatrick Farrell 
163a9603a14SPatrick Farrell   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
1641baa6e33SBarry Smith   if (lmP->H0) PetscCall(MatLMVMSetJ0(lmP->M, lmP->H0));
165a7e14dcfSSatish Balay   PetscFunctionReturn(0);
166a7e14dcfSSatish Balay }
167a7e14dcfSSatish Balay 
168a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
169441846f8SBarry Smith static PetscErrorCode TaoDestroy_LMVM(Tao tao)
170a7e14dcfSSatish Balay {
171a7e14dcfSSatish Balay   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
172a7e14dcfSSatish Balay 
173a7e14dcfSSatish Balay   PetscFunctionBegin;
174a7e14dcfSSatish Balay   if (tao->setupcalled) {
1759566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->Xold));
1769566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->Gold));
1779566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->D));
178a7e14dcfSSatish Balay   }
1799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lmP->M));
1801baa6e33SBarry Smith   if (lmP->H0) PetscCall(PetscObjectDereference((PetscObject)lmP->H0));
1819566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
182a7e14dcfSSatish Balay   PetscFunctionReturn(0);
183a7e14dcfSSatish Balay }
184a7e14dcfSSatish Balay 
185a7e14dcfSSatish Balay /*------------------------------------------------------------*/
1864416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptionItems *PetscOptionsObject,Tao tao)
187a7e14dcfSSatish Balay {
188cd929ea3SAlp Dener   TAO_LMVM       *lm = (TAO_LMVM *)tao->data;
189a7e14dcfSSatish Balay 
190a7e14dcfSSatish Balay   PetscFunctionBegin;
191d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");
1929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_lmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",lm->recycle,&lm->recycle,NULL));
1939566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
1949566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(lm->M));
195d0609cedSBarry Smith   PetscOptionsHeadEnd();
196a7e14dcfSSatish Balay   PetscFunctionReturn(0);
197a7e14dcfSSatish Balay }
198a7e14dcfSSatish Balay 
199a7e14dcfSSatish Balay /*------------------------------------------------------------*/
200441846f8SBarry Smith static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer)
201a7e14dcfSSatish Balay {
202a7e14dcfSSatish Balay   TAO_LMVM       *lm = (TAO_LMVM *)tao->data;
203cd929ea3SAlp Dener   PetscBool      isascii;
2044d6623b4SAlp Dener   PetscInt       recycled_its;
205a7e14dcfSSatish Balay 
206a7e14dcfSSatish Balay   PetscFunctionBegin;
2079566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
208a7e14dcfSSatish Balay   if (isascii) {
20963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Gradient steps: %" PetscInt_FMT "\n", lm->grad));
210cd929ea3SAlp Dener     if (lm->recycle) {
2119566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Recycle: on\n"));
212cd929ea3SAlp Dener       recycled_its = lm->bfgs + lm->grad;
21363a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Total recycled iterations: %" PetscInt_FMT "\n", recycled_its));
214a0bfee83SAlp Dener     }
215a7e14dcfSSatish Balay   }
216a7e14dcfSSatish Balay   PetscFunctionReturn(0);
217a7e14dcfSSatish Balay }
218a7e14dcfSSatish Balay 
219a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
220a7e14dcfSSatish Balay 
2214aa34175SJason Sarich /*MC
2224aa34175SJason Sarich   TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton
2234aa34175SJason Sarich   optimization solver for unconstrained minimization. It solves
2244aa34175SJason Sarich   the Newton step
2254aa34175SJason Sarich           Hkdk = - gk
2264aa34175SJason Sarich 
2274aa34175SJason Sarich   using an approximation Bk in place of Hk, where Bk is composed using
2284aa34175SJason Sarich   the BFGS update formula. A More-Thuente line search is then used
2294aa34175SJason Sarich   to computed the steplength in the dk direction
2300ad3a497SAlp Dener 
2314aa34175SJason Sarich   Options Database Keys:
232a2b725a8SWilliam Gropp +   -tao_lmvm_recycle - enable recycling LMVM updates between TaoSolve() calls
233a2b725a8SWilliam Gropp -   -tao_lmvm_no_scale - (developer) disables diagonal Broyden scaling on the LMVM approximation
2344aa34175SJason Sarich 
2351eb8069cSJason Sarich   Level: beginner
2364aa34175SJason Sarich M*/
2374aa34175SJason Sarich 
238728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao)
239a7e14dcfSSatish Balay {
240a7e14dcfSSatish Balay   TAO_LMVM       *lmP;
2418caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
242a7e14dcfSSatish Balay 
243a7e14dcfSSatish Balay   PetscFunctionBegin;
244a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_LMVM;
245a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_LMVM;
246a7e14dcfSSatish Balay   tao->ops->view = TaoView_LMVM;
247a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_LMVM;
248a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_LMVM;
249a7e14dcfSSatish Balay 
2509566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao,&lmP));
25183c8fe1dSLisandro Dalcin   lmP->D = NULL;
25283c8fe1dSLisandro Dalcin   lmP->M = NULL;
25383c8fe1dSLisandro Dalcin   lmP->Xold = NULL;
25483c8fe1dSLisandro Dalcin   lmP->Gold = NULL;
255a9603a14SPatrick Farrell   lmP->H0   = NULL;
256cd929ea3SAlp Dener   lmP->recycle = PETSC_FALSE;
257a7e14dcfSSatish Balay 
258a7e14dcfSSatish Balay   tao->data = (void*)lmP;
2596552cf8aSJason Sarich   /* Override default settings (unless already changed) */
2606552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 2000;
2616552cf8aSJason Sarich   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
262a7e14dcfSSatish Balay 
2639566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch));
2649566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
2659566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch,morethuente_type));
2669566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao));
2679566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix));
268cd929ea3SAlp Dener 
2699566063dSJacob Faibussowitsch   PetscCall(KSPInitializePackage());
2709566063dSJacob Faibussowitsch   PetscCall(MatCreate(((PetscObject)tao)->comm, &lmP->M));
2719566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)lmP->M, (PetscObject)tao, 1));
2729566063dSJacob Faibussowitsch   PetscCall(MatSetType(lmP->M, MATLMVMBFGS));
2739566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(lmP->M, "tao_lmvm_"));
274a7e14dcfSSatish Balay   PetscFunctionReturn(0);
275a7e14dcfSSatish Balay }
276