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 79371c9d4SSatish Balay static PetscErrorCode TaoSolve_LMVM(Tao tao) { 8a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 9a7e14dcfSSatish Balay PetscReal f, fold, gdx, gnorm; 10a7e14dcfSSatish Balay PetscReal step = 1.0; 11cd929ea3SAlp Dener PetscInt stepType = LMVM_STEP_GRAD, nupdates; 12e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 13a7e14dcfSSatish Balay 14a7e14dcfSSatish Balay PetscFunctionBegin; 15a7e14dcfSSatish Balay 16*48a46eb9SPierre 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); 283ecd9318SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 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) */ 539566063dSJacob Faibussowitsch PetscCall(VecDot(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) { 1169371c9d4SSatish Balay case LMVM_STEP_BFGS: ++lmP->bfgs; break; 1179371c9d4SSatish Balay case LMVM_STEP_GRAD: ++lmP->grad; break; 1189371c9d4SSatish Balay default: break; 119cd929ea3SAlp Dener } 120cd929ea3SAlp Dener /* Compute new gradient norm */ 1219566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm)); 122a7e14dcfSSatish Balay } 123a9603a14SPatrick Farrell 124cd929ea3SAlp Dener /* Check convergence */ 1258931d482SJason Sarich tao->niter++; 1269566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its)); 1279566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step)); 128dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 129a7e14dcfSSatish Balay } 130a7e14dcfSSatish Balay PetscFunctionReturn(0); 131a7e14dcfSSatish Balay } 13287f595a5SBarry Smith 1339371c9d4SSatish Balay static PetscErrorCode TaoSetUp_LMVM(Tao tao) { 134a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 135a7e14dcfSSatish Balay PetscInt n, N; 136b94d7dedSBarry Smith PetscBool is_set, is_spd; 137a7e14dcfSSatish Balay 138a7e14dcfSSatish Balay PetscFunctionBegin; 139a7e14dcfSSatish Balay /* Existence of tao->solution checked in TaoSetUp() */ 1409566063dSJacob Faibussowitsch if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 1419566063dSJacob Faibussowitsch if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 1429566063dSJacob Faibussowitsch if (!lmP->D) PetscCall(VecDuplicate(tao->solution, &lmP->D)); 1439566063dSJacob Faibussowitsch if (!lmP->Xold) PetscCall(VecDuplicate(tao->solution, &lmP->Xold)); 1449566063dSJacob Faibussowitsch if (!lmP->Gold) PetscCall(VecDuplicate(tao->solution, &lmP->Gold)); 145a7e14dcfSSatish Balay 146a7e14dcfSSatish Balay /* Create matrix for the limited memory approximation */ 1479566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution, &n)); 1489566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &N)); 1499566063dSJacob Faibussowitsch PetscCall(MatSetSizes(lmP->M, n, n, N, N)); 1509566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(lmP->M, tao->solution, tao->gradient)); 151b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(lmP->M, &is_set, &is_spd)); 152b94d7dedSBarry Smith PetscCheck(is_set && is_spd, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix is not symmetric positive-definite."); 153a9603a14SPatrick Farrell 154a9603a14SPatrick Farrell /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 1551baa6e33SBarry Smith if (lmP->H0) PetscCall(MatLMVMSetJ0(lmP->M, lmP->H0)); 156a7e14dcfSSatish Balay PetscFunctionReturn(0); 157a7e14dcfSSatish Balay } 158a7e14dcfSSatish Balay 159a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 1609371c9d4SSatish Balay static PetscErrorCode TaoDestroy_LMVM(Tao tao) { 161a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 162a7e14dcfSSatish Balay 163a7e14dcfSSatish Balay PetscFunctionBegin; 164a7e14dcfSSatish Balay if (tao->setupcalled) { 1659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lmP->Xold)); 1669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lmP->Gold)); 1679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lmP->D)); 168a7e14dcfSSatish Balay } 1699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&lmP->M)); 1701baa6e33SBarry Smith if (lmP->H0) PetscCall(PetscObjectDereference((PetscObject)lmP->H0)); 1719566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 172a7e14dcfSSatish Balay PetscFunctionReturn(0); 173a7e14dcfSSatish Balay } 174a7e14dcfSSatish Balay 175a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 1769371c9d4SSatish Balay static PetscErrorCode TaoSetFromOptions_LMVM(Tao tao, PetscOptionItems *PetscOptionsObject) { 177cd929ea3SAlp Dener TAO_LMVM *lm = (TAO_LMVM *)tao->data; 178a7e14dcfSSatish Balay 179a7e14dcfSSatish Balay PetscFunctionBegin; 180d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Limited-memory variable-metric method for unconstrained optimization"); 1819566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_lmvm_recycle", "enable recycling of the BFGS matrix between subsequent TaoSolve() calls", "", lm->recycle, &lm->recycle, NULL)); 1829566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 1839566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(lm->M)); 184d0609cedSBarry Smith PetscOptionsHeadEnd(); 185a7e14dcfSSatish Balay PetscFunctionReturn(0); 186a7e14dcfSSatish Balay } 187a7e14dcfSSatish Balay 188a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 1899371c9d4SSatish Balay static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer) { 190a7e14dcfSSatish Balay TAO_LMVM *lm = (TAO_LMVM *)tao->data; 191cd929ea3SAlp Dener PetscBool isascii; 1924d6623b4SAlp Dener PetscInt recycled_its; 193a7e14dcfSSatish Balay 194a7e14dcfSSatish Balay PetscFunctionBegin; 1959566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 196a7e14dcfSSatish Balay if (isascii) { 19763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Gradient steps: %" PetscInt_FMT "\n", lm->grad)); 198cd929ea3SAlp Dener if (lm->recycle) { 1999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Recycle: on\n")); 200cd929ea3SAlp Dener recycled_its = lm->bfgs + lm->grad; 20163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Total recycled iterations: %" PetscInt_FMT "\n", recycled_its)); 202a0bfee83SAlp Dener } 203a7e14dcfSSatish Balay } 204a7e14dcfSSatish Balay PetscFunctionReturn(0); 205a7e14dcfSSatish Balay } 206a7e14dcfSSatish Balay 207a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 208a7e14dcfSSatish Balay 2094aa34175SJason Sarich /*MC 2104aa34175SJason Sarich TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton 2114aa34175SJason Sarich optimization solver for unconstrained minimization. It solves 2124aa34175SJason Sarich the Newton step 2134aa34175SJason Sarich Hkdk = - gk 2144aa34175SJason Sarich 2154aa34175SJason Sarich using an approximation Bk in place of Hk, where Bk is composed using 2164aa34175SJason Sarich the BFGS update formula. A More-Thuente line search is then used 2174aa34175SJason Sarich to computed the steplength in the dk direction 2180ad3a497SAlp Dener 2194aa34175SJason Sarich Options Database Keys: 220a2b725a8SWilliam Gropp + -tao_lmvm_recycle - enable recycling LMVM updates between TaoSolve() calls 221a2b725a8SWilliam Gropp - -tao_lmvm_no_scale - (developer) disables diagonal Broyden scaling on the LMVM approximation 2224aa34175SJason Sarich 2231eb8069cSJason Sarich Level: beginner 2244aa34175SJason Sarich M*/ 2254aa34175SJason Sarich 2269371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao) { 227a7e14dcfSSatish Balay TAO_LMVM *lmP; 2288caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 229a7e14dcfSSatish Balay 230a7e14dcfSSatish Balay PetscFunctionBegin; 231a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_LMVM; 232a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_LMVM; 233a7e14dcfSSatish Balay tao->ops->view = TaoView_LMVM; 234a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_LMVM; 235a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_LMVM; 236a7e14dcfSSatish Balay 2379566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao, &lmP)); 23883c8fe1dSLisandro Dalcin lmP->D = NULL; 23983c8fe1dSLisandro Dalcin lmP->M = NULL; 24083c8fe1dSLisandro Dalcin lmP->Xold = NULL; 24183c8fe1dSLisandro Dalcin lmP->Gold = NULL; 242a9603a14SPatrick Farrell lmP->H0 = NULL; 243cd929ea3SAlp Dener lmP->recycle = PETSC_FALSE; 244a7e14dcfSSatish Balay 245a7e14dcfSSatish Balay tao->data = (void *)lmP; 2466552cf8aSJason Sarich /* Override default settings (unless already changed) */ 2476552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 2486552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 249a7e14dcfSSatish Balay 2509566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 2529566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type)); 2539566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao)); 2549566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix)); 255cd929ea3SAlp Dener 2569566063dSJacob Faibussowitsch PetscCall(KSPInitializePackage()); 2579566063dSJacob Faibussowitsch PetscCall(MatCreate(((PetscObject)tao)->comm, &lmP->M)); 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)lmP->M, (PetscObject)tao, 1)); 2599566063dSJacob Faibussowitsch PetscCall(MatSetType(lmP->M, MATLMVMBFGS)); 2609566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(lmP->M, "tao_lmvm_")); 261a7e14dcfSSatish Balay PetscFunctionReturn(0); 262a7e14dcfSSatish Balay } 263