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; 12a7e14dcfSSatish Balay PetscErrorCode ierr; 13cd929ea3SAlp Dener PetscInt stepType = LMVM_STEP_GRAD, nupdates; 14e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 15a7e14dcfSSatish Balay 16a7e14dcfSSatish Balay PetscFunctionBegin; 17a7e14dcfSSatish Balay 18a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 19a7e14dcfSSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n");CHKERRQ(ierr); 20a7e14dcfSSatish Balay } 21a7e14dcfSSatish Balay 22a7e14dcfSSatish Balay /* Check convergence criteria */ 23a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 24a9603a14SPatrick Farrell ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 25a9603a14SPatrick Farrell 2687f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 27a7e14dcfSSatish Balay 283ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 293ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 303ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 313ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 323ecd9318SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 33a7e14dcfSSatish Balay 34a7e14dcfSSatish Balay /* Set counter for gradient/reset steps */ 35cd929ea3SAlp Dener if (!lmP->recycle) { 36a7e14dcfSSatish Balay lmP->bfgs = 0; 37a7e14dcfSSatish Balay lmP->grad = 0; 38cd929ea3SAlp Dener ierr = MatLMVMReset(lmP->M, PETSC_FALSE); CHKERRQ(ierr); 39de6ffafeSAlp Dener } 40a7e14dcfSSatish Balay 41a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 423ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 43a7e14dcfSSatish Balay /* Compute direction */ 44cd929ea3SAlp Dener if (lmP->H0) { 45cd929ea3SAlp Dener ierr = MatLMVMSetJ0(lmP->M, lmP->H0);CHKERRQ(ierr); 46cd929ea3SAlp Dener stepType = LMVM_STEP_BFGS; 47cd929ea3SAlp Dener } 48a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 49cd929ea3SAlp Dener ierr = MatSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 50cd929ea3SAlp Dener ierr = MatLMVMGetUpdateCount(lmP->M, &nupdates); CHKERRQ(ierr); 51cd929ea3SAlp Dener if (nupdates > 0) stepType = LMVM_STEP_BFGS; 52a7e14dcfSSatish Balay 53a7e14dcfSSatish Balay /* Check for success (descent direction) */ 54a7e14dcfSSatish Balay ierr = VecDot(lmP->D, tao->gradient, &gdx);CHKERRQ(ierr); 55a7e14dcfSSatish Balay if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) { 56a7e14dcfSSatish Balay /* Step is not descent or direction produced not a number 57a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 58a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 59a7e14dcfSSatish Balay which is guaranteed to be descent 60a7e14dcfSSatish Balay 61a7e14dcfSSatish Balay Use steepest descent direction (scaled) 62a7e14dcfSSatish Balay */ 63a7e14dcfSSatish Balay 64*b2d8c577SAlp Dener ierr = MatLMVMResetJ0(lmP->M);CHKERRQ(ierr); 65cd929ea3SAlp Dener ierr = MatLMVMReset(lmP->M, PETSC_FALSE);CHKERRQ(ierr); 66a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 67cd929ea3SAlp Dener ierr = MatSolve(lmP->M,tao->gradient, lmP->D);CHKERRQ(ierr); 68a7e14dcfSSatish Balay 69a7e14dcfSSatish Balay /* On a reset, the direction cannot be not a number; it is a 70a7e14dcfSSatish Balay scaled gradient step. No need to check for this condition. */ 71cd929ea3SAlp Dener stepType = LMVM_STEP_GRAD; 72a7e14dcfSSatish Balay } 73a7e14dcfSSatish Balay ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 74a7e14dcfSSatish Balay 75a7e14dcfSSatish Balay /* Perform the linesearch */ 76a7e14dcfSSatish Balay fold = f; 77a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr); 78a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr); 79a7e14dcfSSatish Balay 80a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status);CHKERRQ(ierr); 81a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 82a7e14dcfSSatish Balay 83cd929ea3SAlp Dener if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_STEP_GRAD)) { 84a7e14dcfSSatish Balay /* Reset factors and use scaled gradient step */ 85a7e14dcfSSatish Balay f = fold; 86a7e14dcfSSatish Balay ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 87a7e14dcfSSatish Balay ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 88a7e14dcfSSatish Balay 89a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with BFGS step */ 90a7e14dcfSSatish Balay /* Attempt to use the scaled gradient direction */ 91a7e14dcfSSatish Balay 92*b2d8c577SAlp Dener ierr = MatLMVMResetJ0(lmP->M);CHKERRQ(ierr); 93cd929ea3SAlp Dener ierr = MatLMVMReset(lmP->M, PETSC_FALSE);CHKERRQ(ierr); 94a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 95cd929ea3SAlp Dener ierr = MatSolve(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 96a7e14dcfSSatish Balay 97a7e14dcfSSatish Balay /* On a reset, the direction cannot be not a number; it is a 98a7e14dcfSSatish Balay scaled gradient step. No need to check for this condition. */ 99cd929ea3SAlp Dener stepType = LMVM_STEP_GRAD; 100a7e14dcfSSatish Balay ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 101a7e14dcfSSatish Balay 102a7e14dcfSSatish Balay /* Perform the linesearch */ 103a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status);CHKERRQ(ierr); 104a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 105a7e14dcfSSatish Balay } 106a7e14dcfSSatish Balay 10787f595a5SBarry Smith if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 108a7e14dcfSSatish Balay /* Failed to find an improving point */ 109a7e14dcfSSatish Balay f = fold; 110a7e14dcfSSatish Balay ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 111a7e14dcfSSatish Balay ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 112a7e14dcfSSatish Balay step = 0.0; 113a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 114cd929ea3SAlp Dener } else { 115cd929ea3SAlp Dener /* LS found valid step, so tally up step type */ 116cd929ea3SAlp Dener switch (stepType) { 117cd929ea3SAlp Dener case LMVM_STEP_BFGS: 118cd929ea3SAlp Dener ++lmP->bfgs; 119cd929ea3SAlp Dener break; 120cd929ea3SAlp Dener case LMVM_STEP_GRAD: 121cd929ea3SAlp Dener ++lmP->grad; 122cd929ea3SAlp Dener break; 123cd929ea3SAlp Dener default: 124cd929ea3SAlp Dener break; 125cd929ea3SAlp Dener } 126cd929ea3SAlp Dener /* Compute new gradient norm */ 127cd929ea3SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 128a7e14dcfSSatish Balay } 129a9603a14SPatrick Farrell 130cd929ea3SAlp Dener /* Check convergence */ 1318931d482SJason Sarich tao->niter++; 1323ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 1333ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 1343ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 135a7e14dcfSSatish Balay } 136a7e14dcfSSatish Balay PetscFunctionReturn(0); 137a7e14dcfSSatish Balay } 13887f595a5SBarry Smith 139441846f8SBarry Smith static PetscErrorCode TaoSetUp_LMVM(Tao tao) 140a7e14dcfSSatish Balay { 141a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 142a7e14dcfSSatish Balay PetscInt n,N; 143a7e14dcfSSatish Balay PetscErrorCode ierr; 144a7e14dcfSSatish Balay 145a7e14dcfSSatish Balay PetscFunctionBegin; 146a7e14dcfSSatish Balay /* Existence of tao->solution checked in TaoSetUp() */ 147a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); } 148a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); } 149a7e14dcfSSatish Balay if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr); } 150a7e14dcfSSatish Balay if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr); } 151a7e14dcfSSatish Balay if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr); } 152a7e14dcfSSatish Balay 153a7e14dcfSSatish Balay /* Create matrix for the limited memory approximation */ 154a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 155a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 156cd929ea3SAlp Dener ierr = MatSetSizes(lmP->M, n, n, N, N);CHKERRQ(ierr); 157cd929ea3SAlp Dener ierr = MatLMVMAllocate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 158a9603a14SPatrick Farrell 159a9603a14SPatrick Farrell /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 160a9603a14SPatrick Farrell if (lmP->H0) { 161cd929ea3SAlp Dener ierr = MatLMVMSetJ0(lmP->M, lmP->H0);CHKERRQ(ierr); 162a9603a14SPatrick Farrell } 163a9603a14SPatrick Farrell 164a7e14dcfSSatish Balay PetscFunctionReturn(0); 165a7e14dcfSSatish Balay } 166a7e14dcfSSatish Balay 167a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 168441846f8SBarry Smith static PetscErrorCode TaoDestroy_LMVM(Tao tao) 169a7e14dcfSSatish Balay { 170a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 171a7e14dcfSSatish Balay PetscErrorCode ierr; 172a7e14dcfSSatish Balay 173a7e14dcfSSatish Balay PetscFunctionBegin; 174a7e14dcfSSatish Balay if (tao->setupcalled) { 175a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr); 176a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr); 177a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->D);CHKERRQ(ierr); 178a7e14dcfSSatish Balay } 179cd929ea3SAlp Dener ierr = MatDestroy(&lmP->M);CHKERRQ(ierr); 180a9603a14SPatrick Farrell if (lmP->H0) { 181a9603a14SPatrick Farrell ierr = PetscObjectDereference((PetscObject)lmP->H0);CHKERRQ(ierr); 182a9603a14SPatrick Farrell } 183a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 184a9603a14SPatrick Farrell 185a7e14dcfSSatish Balay PetscFunctionReturn(0); 186a7e14dcfSSatish Balay } 187a7e14dcfSSatish Balay 188a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 1894416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptionItems *PetscOptionsObject,Tao tao) 190a7e14dcfSSatish Balay { 191cd929ea3SAlp Dener TAO_LMVM *lm = (TAO_LMVM *)tao->data; 192a7e14dcfSSatish Balay PetscErrorCode ierr; 193a7e14dcfSSatish Balay 194a7e14dcfSSatish Balay PetscFunctionBegin; 1951a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");CHKERRQ(ierr); 196cd929ea3SAlp Dener ierr = PetscOptionsBool("-tao_lmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",lm->recycle,&lm->recycle,NULL);CHKERRQ(ierr); 197114d2d62SAlp Dener ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 198cd929ea3SAlp Dener ierr = MatSetFromOptions(lm->M);CHKERRQ(ierr); 199288b7216SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 200a7e14dcfSSatish Balay PetscFunctionReturn(0); 201a7e14dcfSSatish Balay } 202a7e14dcfSSatish Balay 203a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 204441846f8SBarry Smith static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer) 205a7e14dcfSSatish Balay { 206a7e14dcfSSatish Balay TAO_LMVM *lm = (TAO_LMVM *)tao->data; 207cd929ea3SAlp Dener PetscBool isascii; 2084d6623b4SAlp Dener PetscInt recycled_its; 209a7e14dcfSSatish Balay PetscErrorCode ierr; 210a7e14dcfSSatish Balay 211a7e14dcfSSatish Balay PetscFunctionBegin; 212a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 213a7e14dcfSSatish Balay if (isascii) { 214a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, " Gradient steps: %D\n", lm->grad);CHKERRQ(ierr); 215cd929ea3SAlp Dener if (lm->recycle) { 216288b7216SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " Recycle: on\n");CHKERRQ(ierr); 217cd929ea3SAlp Dener recycled_its = lm->bfgs + lm->grad; 2184d6623b4SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr); 219a0bfee83SAlp Dener } 220a7e14dcfSSatish Balay } 221a7e14dcfSSatish Balay PetscFunctionReturn(0); 222a7e14dcfSSatish Balay } 223a7e14dcfSSatish Balay 224a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 225a7e14dcfSSatish Balay 2264aa34175SJason Sarich /*MC 2274aa34175SJason Sarich TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton 2284aa34175SJason Sarich optimization solver for unconstrained minimization. It solves 2294aa34175SJason Sarich the Newton step 2304aa34175SJason Sarich Hkdk = - gk 2314aa34175SJason Sarich 2324aa34175SJason Sarich using an approximation Bk in place of Hk, where Bk is composed using 2334aa34175SJason Sarich the BFGS update formula. A More-Thuente line search is then used 2344aa34175SJason Sarich to computed the steplength in the dk direction 2354aa34175SJason Sarich Options Database Keys: 236cd929ea3SAlp Dener . -tao_lmvm_recycle - enable recycling LMVM updates between TaoSolve() calls 2374aa34175SJason Sarich 2381eb8069cSJason Sarich Level: beginner 2394aa34175SJason Sarich M*/ 2404aa34175SJason Sarich 241728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao) 242a7e14dcfSSatish Balay { 243a7e14dcfSSatish Balay TAO_LMVM *lmP; 2448caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 245a7e14dcfSSatish Balay PetscErrorCode ierr; 246a7e14dcfSSatish Balay 247a7e14dcfSSatish Balay PetscFunctionBegin; 248a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_LMVM; 249a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_LMVM; 250a7e14dcfSSatish Balay tao->ops->view = TaoView_LMVM; 251a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_LMVM; 252a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_LMVM; 253a7e14dcfSSatish Balay 2543c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr); 255a7e14dcfSSatish Balay lmP->D = 0; 256a7e14dcfSSatish Balay lmP->M = 0; 257a7e14dcfSSatish Balay lmP->Xold = 0; 258a7e14dcfSSatish Balay lmP->Gold = 0; 259a9603a14SPatrick Farrell lmP->H0 = NULL; 260cd929ea3SAlp Dener lmP->recycle = PETSC_FALSE; 261a7e14dcfSSatish Balay 262a7e14dcfSSatish Balay tao->data = (void*)lmP; 2636552cf8aSJason Sarich /* Override default settings (unless already changed) */ 2646552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 2656552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 266a7e14dcfSSatish Balay 267a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 26863b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 269a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 270441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 2715d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 272cd929ea3SAlp Dener 273cd929ea3SAlp Dener ierr = KSPInitializePackage();CHKERRQ(ierr); 274cd929ea3SAlp Dener ierr = MatCreate(((PetscObject)tao)->comm, &lmP->M);CHKERRQ(ierr); 275cd929ea3SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)lmP->M, (PetscObject)tao, 1);CHKERRQ(ierr); 276cd929ea3SAlp Dener ierr = MatSetType(lmP->M, MATLBFGS);CHKERRQ(ierr); 277cd929ea3SAlp Dener ierr = MatSetOptionsPrefix(lmP->M, "tao_lmvm_");CHKERRQ(ierr); 278a7e14dcfSSatish Balay PetscFunctionReturn(0); 279a7e14dcfSSatish Balay } 280