xref: /petsc/src/tao/unconstrained/impls/lmvm/lmvm.c (revision bbb72809d0d086fcc837899733ecccb12a5516ec)
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 
7d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_LMVM(Tao tao)
8d71ae5a4SJacob Faibussowitsch {
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;
1648a46eb9SPierre Jolivet   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"));
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay   /*  Check convergence criteria */
199566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
209566063dSJacob Faibussowitsch   PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
21a9603a14SPatrick Farrell 
223c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
23a7e14dcfSSatish Balay 
243ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
259566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
269566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
27dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
283ba16761SJacob Faibussowitsch   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
29a7e14dcfSSatish Balay 
30a7e14dcfSSatish Balay   /*  Set counter for gradient/reset steps */
31cd929ea3SAlp Dener   if (!lmP->recycle) {
32a7e14dcfSSatish Balay     lmP->bfgs = 0;
33a7e14dcfSSatish Balay     lmP->grad = 0;
349566063dSJacob Faibussowitsch     PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
35de6ffafeSAlp Dener   }
36a7e14dcfSSatish Balay 
37a7e14dcfSSatish Balay   /*  Have not converged; continue with Newton method */
383ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
39e1e80dc8SAlp Dener     /* Call general purpose update function */
40dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
41e1e80dc8SAlp Dener 
42a7e14dcfSSatish Balay     /*  Compute direction */
43cd929ea3SAlp Dener     if (lmP->H0) {
449566063dSJacob Faibussowitsch       PetscCall(MatLMVMSetJ0(lmP->M, lmP->H0));
45cd929ea3SAlp Dener       stepType = LMVM_STEP_BFGS;
46cd929ea3SAlp Dener     }
479566063dSJacob Faibussowitsch     PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
489566063dSJacob Faibussowitsch     PetscCall(MatSolve(lmP->M, tao->gradient, lmP->D));
499566063dSJacob Faibussowitsch     PetscCall(MatLMVMGetUpdateCount(lmP->M, &nupdates));
50cd929ea3SAlp Dener     if (nupdates > 0) stepType = LMVM_STEP_BFGS;
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay     /*  Check for success (descent direction) */
53*bbb72809SHansol Suh     PetscCall(VecDotRealPart(lmP->D, tao->gradient, &gdx));
54a7e14dcfSSatish Balay     if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
55a7e14dcfSSatish Balay       /* Step is not descent or direction produced not a number
56a7e14dcfSSatish Balay          We can assert bfgsUpdates > 1 in this case because
57a7e14dcfSSatish Balay          the first solve produces the scaled gradient direction,
58a7e14dcfSSatish Balay          which is guaranteed to be descent
59a7e14dcfSSatish Balay 
60a7e14dcfSSatish Balay          Use steepest descent direction (scaled)
61a7e14dcfSSatish Balay       */
62a7e14dcfSSatish Balay 
639566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
649566063dSJacob Faibussowitsch       PetscCall(MatLMVMClearJ0(lmP->M));
659566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
669566063dSJacob Faibussowitsch       PetscCall(MatSolve(lmP->M, tao->gradient, lmP->D));
67a7e14dcfSSatish Balay 
68a7e14dcfSSatish Balay       /* On a reset, the direction cannot be not a number; it is a
69a7e14dcfSSatish Balay          scaled gradient step.  No need to check for this condition. */
70cd929ea3SAlp Dener       stepType = LMVM_STEP_GRAD;
71a7e14dcfSSatish Balay     }
729566063dSJacob Faibussowitsch     PetscCall(VecScale(lmP->D, -1.0));
73a7e14dcfSSatish Balay 
74a7e14dcfSSatish Balay     /*  Perform the linesearch */
75a7e14dcfSSatish Balay     fold = f;
769566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, lmP->Xold));
779566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, lmP->Gold));
78a7e14dcfSSatish Balay 
799566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status));
809566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
81a7e14dcfSSatish Balay 
82cd929ea3SAlp Dener     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_STEP_GRAD)) {
83a7e14dcfSSatish Balay       /*  Reset factors and use scaled gradient step */
84a7e14dcfSSatish Balay       f = fold;
859566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Xold, tao->solution));
869566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Gold, tao->gradient));
87a7e14dcfSSatish Balay 
88a7e14dcfSSatish Balay       /*  Failed to obtain acceptable iterate with BFGS step */
89a7e14dcfSSatish Balay       /*  Attempt to use the scaled gradient direction */
90a7e14dcfSSatish Balay 
919566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
929566063dSJacob Faibussowitsch       PetscCall(MatLMVMClearJ0(lmP->M));
939566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
949566063dSJacob Faibussowitsch       PetscCall(MatSolve(lmP->M, tao->solution, tao->gradient));
95a7e14dcfSSatish Balay 
96a7e14dcfSSatish Balay       /* On a reset, the direction cannot be not a number; it is a
97a7e14dcfSSatish Balay           scaled gradient step.  No need to check for this condition. */
98cd929ea3SAlp Dener       stepType = LMVM_STEP_GRAD;
999566063dSJacob Faibussowitsch       PetscCall(VecScale(lmP->D, -1.0));
100a7e14dcfSSatish Balay 
101a7e14dcfSSatish Balay       /*  Perform the linesearch */
1029566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status));
1039566063dSJacob Faibussowitsch       PetscCall(TaoAddLineSearchCounts(tao));
104a7e14dcfSSatish Balay     }
105a7e14dcfSSatish Balay 
10687f595a5SBarry Smith     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
107a7e14dcfSSatish Balay       /*  Failed to find an improving point */
108a7e14dcfSSatish Balay       f = fold;
1099566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Xold, tao->solution));
1109566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Gold, tao->gradient));
111a7e14dcfSSatish Balay       step        = 0.0;
112a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
113cd929ea3SAlp Dener     } else {
114cd929ea3SAlp Dener       /* LS found valid step, so tally up step type */
115cd929ea3SAlp Dener       switch (stepType) {
116d71ae5a4SJacob Faibussowitsch       case LMVM_STEP_BFGS:
117d71ae5a4SJacob Faibussowitsch         ++lmP->bfgs;
118d71ae5a4SJacob Faibussowitsch         break;
119d71ae5a4SJacob Faibussowitsch       case LMVM_STEP_GRAD:
120d71ae5a4SJacob Faibussowitsch         ++lmP->grad;
121d71ae5a4SJacob Faibussowitsch         break;
122d71ae5a4SJacob Faibussowitsch       default:
123d71ae5a4SJacob Faibussowitsch         break;
124cd929ea3SAlp Dener       }
125cd929ea3SAlp Dener       /*  Compute new gradient norm */
1269566063dSJacob Faibussowitsch       PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
127a7e14dcfSSatish Balay     }
128a9603a14SPatrick Farrell 
129cd929ea3SAlp Dener     /* Check convergence */
1308931d482SJason Sarich     tao->niter++;
1319566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
1329566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
133dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
134a7e14dcfSSatish Balay   }
1353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
136a7e14dcfSSatish Balay }
13787f595a5SBarry Smith 
138d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_LMVM(Tao tao)
139d71ae5a4SJacob Faibussowitsch {
140a7e14dcfSSatish Balay   TAO_LMVM *lmP = (TAO_LMVM *)tao->data;
141a7e14dcfSSatish Balay   PetscInt  n, N;
142b94d7dedSBarry Smith   PetscBool is_set, is_spd;
143a7e14dcfSSatish Balay 
144a7e14dcfSSatish Balay   PetscFunctionBegin;
145a7e14dcfSSatish Balay   /* Existence of tao->solution checked in TaoSetUp() */
1469566063dSJacob Faibussowitsch   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
1479566063dSJacob Faibussowitsch   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
1489566063dSJacob Faibussowitsch   if (!lmP->D) PetscCall(VecDuplicate(tao->solution, &lmP->D));
1499566063dSJacob Faibussowitsch   if (!lmP->Xold) PetscCall(VecDuplicate(tao->solution, &lmP->Xold));
1509566063dSJacob Faibussowitsch   if (!lmP->Gold) PetscCall(VecDuplicate(tao->solution, &lmP->Gold));
151a7e14dcfSSatish Balay 
152a7e14dcfSSatish Balay   /*  Create matrix for the limited memory approximation */
1539566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(tao->solution, &n));
1549566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tao->solution, &N));
1559566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(lmP->M, n, n, N, N));
1569566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(lmP->M, tao->solution, tao->gradient));
157b94d7dedSBarry Smith   PetscCall(MatIsSPDKnown(lmP->M, &is_set, &is_spd));
158b94d7dedSBarry Smith   PetscCheck(is_set && is_spd, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix is not symmetric positive-definite.");
159a9603a14SPatrick Farrell 
160a9603a14SPatrick Farrell   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
1611baa6e33SBarry Smith   if (lmP->H0) PetscCall(MatLMVMSetJ0(lmP->M, lmP->H0));
1623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
163a7e14dcfSSatish Balay }
164a7e14dcfSSatish Balay 
165a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
166d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_LMVM(Tao tao)
167d71ae5a4SJacob Faibussowitsch {
168a7e14dcfSSatish Balay   TAO_LMVM *lmP = (TAO_LMVM *)tao->data;
169a7e14dcfSSatish Balay 
170a7e14dcfSSatish Balay   PetscFunctionBegin;
171a7e14dcfSSatish Balay   if (tao->setupcalled) {
1729566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->Xold));
1739566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->Gold));
1749566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->D));
175a7e14dcfSSatish Balay   }
1769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lmP->M));
1771baa6e33SBarry Smith   if (lmP->H0) PetscCall(PetscObjectDereference((PetscObject)lmP->H0));
1789566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
1793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
180a7e14dcfSSatish Balay }
181a7e14dcfSSatish Balay 
182a7e14dcfSSatish Balay /*------------------------------------------------------------*/
183d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetFromOptions_LMVM(Tao tao, PetscOptionItems *PetscOptionsObject)
184d71ae5a4SJacob Faibussowitsch {
185cd929ea3SAlp Dener   TAO_LMVM *lm = (TAO_LMVM *)tao->data;
186a7e14dcfSSatish Balay 
187a7e14dcfSSatish Balay   PetscFunctionBegin;
188d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Limited-memory variable-metric method for unconstrained optimization");
1899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_lmvm_recycle", "enable recycling of the BFGS matrix between subsequent TaoSolve() calls", "", lm->recycle, &lm->recycle, NULL));
1909566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
1919566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(lm->M));
192d0609cedSBarry Smith   PetscOptionsHeadEnd();
1933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
194a7e14dcfSSatish Balay }
195a7e14dcfSSatish Balay 
196a7e14dcfSSatish Balay /*------------------------------------------------------------*/
197d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer)
198d71ae5a4SJacob Faibussowitsch {
199a7e14dcfSSatish Balay   TAO_LMVM *lm = (TAO_LMVM *)tao->data;
200cd929ea3SAlp Dener   PetscBool isascii;
2014d6623b4SAlp Dener   PetscInt  recycled_its;
202a7e14dcfSSatish Balay 
203a7e14dcfSSatish Balay   PetscFunctionBegin;
2049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
205a7e14dcfSSatish Balay   if (isascii) {
20663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Gradient steps: %" PetscInt_FMT "\n", lm->grad));
207cd929ea3SAlp Dener     if (lm->recycle) {
2089566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Recycle: on\n"));
209cd929ea3SAlp Dener       recycled_its = lm->bfgs + lm->grad;
21063a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Total recycled iterations: %" PetscInt_FMT "\n", recycled_its));
211a0bfee83SAlp Dener     }
212a7e14dcfSSatish Balay   }
2133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
214a7e14dcfSSatish Balay }
215a7e14dcfSSatish Balay 
216a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
217a7e14dcfSSatish Balay 
2184aa34175SJason Sarich /*MC
2194aa34175SJason Sarich   TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton
2204aa34175SJason Sarich   optimization solver for unconstrained minimization. It solves
2214aa34175SJason Sarich   the Newton step
2224aa34175SJason Sarich           Hkdk = - gk
2234aa34175SJason Sarich 
2244aa34175SJason Sarich   using an approximation Bk in place of Hk, where Bk is composed using
2254aa34175SJason Sarich   the BFGS update formula. A More-Thuente line search is then used
2264aa34175SJason Sarich   to computed the steplength in the dk direction
2270ad3a497SAlp Dener 
2284aa34175SJason Sarich   Options Database Keys:
229a2b725a8SWilliam Gropp +   -tao_lmvm_recycle - enable recycling LMVM updates between TaoSolve() calls
230a2b725a8SWilliam Gropp -   -tao_lmvm_no_scale - (developer) disables diagonal Broyden scaling on the LMVM approximation
2314aa34175SJason Sarich 
2321eb8069cSJason Sarich   Level: beginner
2334aa34175SJason Sarich M*/
2344aa34175SJason Sarich 
235d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao)
236d71ae5a4SJacob Faibussowitsch {
237a7e14dcfSSatish Balay   TAO_LMVM   *lmP;
2388caf6e8cSBarry Smith   const char *morethuente_type = TAOLINESEARCHMT;
239a7e14dcfSSatish Balay 
240a7e14dcfSSatish Balay   PetscFunctionBegin;
241a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetUp_LMVM;
242a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_LMVM;
243a7e14dcfSSatish Balay   tao->ops->view           = TaoView_LMVM;
244a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_LMVM;
245a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_LMVM;
246a7e14dcfSSatish Balay 
2474dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lmP));
24883c8fe1dSLisandro Dalcin   lmP->D       = NULL;
24983c8fe1dSLisandro Dalcin   lmP->M       = NULL;
25083c8fe1dSLisandro Dalcin   lmP->Xold    = NULL;
25183c8fe1dSLisandro Dalcin   lmP->Gold    = NULL;
252a9603a14SPatrick Farrell   lmP->H0      = NULL;
253cd929ea3SAlp Dener   lmP->recycle = PETSC_FALSE;
254a7e14dcfSSatish Balay 
255a7e14dcfSSatish Balay   tao->data = (void *)lmP;
2566552cf8aSJason Sarich   /* Override default settings (unless already changed) */
2576552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 2000;
2586552cf8aSJason Sarich   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
259a7e14dcfSSatish Balay 
2609566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
2619566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
2629566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
2639566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
2649566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
265cd929ea3SAlp Dener 
2669566063dSJacob Faibussowitsch   PetscCall(KSPInitializePackage());
2679566063dSJacob Faibussowitsch   PetscCall(MatCreate(((PetscObject)tao)->comm, &lmP->M));
2689566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)lmP->M, (PetscObject)tao, 1));
2699566063dSJacob Faibussowitsch   PetscCall(MatSetType(lmP->M, MATLMVMBFGS));
2709566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(lmP->M, "tao_lmvm_"));
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
272a7e14dcfSSatish Balay }
273