1ba92ff59SBarry Smith #include <petsctaolinesearch.h> 2aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h> 3aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/lmvm/lmvm.h> 4a7e14dcfSSatish Balay 5a7e14dcfSSatish Balay #define LMVM_BFGS 0 6a7e14dcfSSatish Balay #define LMVM_SCALED_GRADIENT 1 7a7e14dcfSSatish Balay #define LMVM_GRADIENT 2 8a7e14dcfSSatish Balay 9a7e14dcfSSatish Balay #undef __FUNCT__ 10a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_LMVM" 11441846f8SBarry Smith static PetscErrorCode TaoSolve_LMVM(Tao tao) 12a7e14dcfSSatish Balay { 13a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 14a7e14dcfSSatish Balay PetscReal f, fold, gdx, gnorm; 15a7e14dcfSSatish Balay PetscReal step = 1.0; 16a7e14dcfSSatish Balay PetscReal delta; 17a7e14dcfSSatish Balay PetscErrorCode ierr; 18a7e14dcfSSatish Balay PetscInt stepType; 19a7e14dcfSSatish Balay PetscInt iter = 0; 20*e4cb33bbSBarry Smith TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 21*e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 22a7e14dcfSSatish Balay 23a7e14dcfSSatish Balay PetscFunctionBegin; 24a7e14dcfSSatish Balay 25a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 26a7e14dcfSSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n");CHKERRQ(ierr); 27a7e14dcfSSatish Balay } 28a7e14dcfSSatish Balay 29a7e14dcfSSatish Balay /* Check convergence criteria */ 30a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 31a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 3287f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 33a7e14dcfSSatish Balay 34a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr); 3587f595a5SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 36a7e14dcfSSatish Balay 37a7e14dcfSSatish Balay /* Set initial scaling for the function */ 38a7e14dcfSSatish Balay if (f != 0.0) { 39a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 4087f595a5SBarry Smith } else { 41a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 42a7e14dcfSSatish Balay } 43a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(lmP->M,delta);CHKERRQ(ierr); 44a7e14dcfSSatish Balay 45a7e14dcfSSatish Balay /* Set counter for gradient/reset steps */ 46a7e14dcfSSatish Balay lmP->bfgs = 0; 47a7e14dcfSSatish Balay lmP->sgrad = 0; 48a7e14dcfSSatish Balay lmP->grad = 0; 49a7e14dcfSSatish Balay 50a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 51a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 52a7e14dcfSSatish Balay /* Compute direction */ 53a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 54a7e14dcfSSatish Balay ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 55a7e14dcfSSatish Balay ++lmP->bfgs; 56a7e14dcfSSatish Balay 57a7e14dcfSSatish Balay /* Check for success (descent direction) */ 58a7e14dcfSSatish Balay ierr = VecDot(lmP->D, tao->gradient, &gdx);CHKERRQ(ierr); 59a7e14dcfSSatish Balay if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) { 60a7e14dcfSSatish Balay /* Step is not descent or direction produced not a number 61a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 62a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 63a7e14dcfSSatish Balay which is guaranteed to be descent 64a7e14dcfSSatish Balay 65a7e14dcfSSatish Balay Use steepest descent direction (scaled) 66a7e14dcfSSatish Balay */ 67a7e14dcfSSatish Balay 68a7e14dcfSSatish Balay ++lmP->grad; 69a7e14dcfSSatish Balay 70a7e14dcfSSatish Balay if (f != 0.0) { 71a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 7287f595a5SBarry Smith } else { 73a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 74a7e14dcfSSatish Balay } 75a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr); 76a7e14dcfSSatish Balay ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 77a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 78a7e14dcfSSatish Balay ierr = MatLMVMSolve(lmP->M,tao->gradient, lmP->D);CHKERRQ(ierr); 79a7e14dcfSSatish Balay 80a7e14dcfSSatish Balay /* On a reset, the direction cannot be not a number; it is a 81a7e14dcfSSatish Balay scaled gradient step. No need to check for this condition. */ 82a7e14dcfSSatish Balay 83a7e14dcfSSatish Balay lmP->bfgs = 1; 84a7e14dcfSSatish Balay ++lmP->sgrad; 85a7e14dcfSSatish Balay stepType = LMVM_SCALED_GRADIENT; 8687f595a5SBarry Smith } else { 87a7e14dcfSSatish Balay if (1 == lmP->bfgs) { 88a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 89a7e14dcfSSatish Balay ++lmP->sgrad; 90a7e14dcfSSatish Balay stepType = LMVM_SCALED_GRADIENT; 9187f595a5SBarry Smith } else { 92a7e14dcfSSatish Balay ++lmP->bfgs; 93a7e14dcfSSatish Balay stepType = LMVM_BFGS; 94a7e14dcfSSatish Balay } 95a7e14dcfSSatish Balay } 96a7e14dcfSSatish Balay ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 97a7e14dcfSSatish Balay 98a7e14dcfSSatish Balay /* Perform the linesearch */ 99a7e14dcfSSatish Balay fold = f; 100a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr); 101a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr); 102a7e14dcfSSatish Balay 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 10687f595a5SBarry Smith while (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_GRADIENT)) { 107a7e14dcfSSatish Balay /* Linesearch failed */ 108a7e14dcfSSatish Balay /* Reset factors and use scaled gradient step */ 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 113a7e14dcfSSatish Balay switch(stepType) { 114a7e14dcfSSatish Balay case LMVM_BFGS: 115a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with BFGS step */ 116a7e14dcfSSatish Balay /* Attempt to use the scaled gradient direction */ 117a7e14dcfSSatish Balay 118a7e14dcfSSatish Balay if (f != 0.0) { 119a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 12087f595a5SBarry Smith } else { 121a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 122a7e14dcfSSatish Balay } 123a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr); 124a7e14dcfSSatish Balay ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 125a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 126a7e14dcfSSatish Balay ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 127a7e14dcfSSatish Balay 128a7e14dcfSSatish Balay /* On a reset, the direction cannot be not a number; it is a 129a7e14dcfSSatish Balay scaled gradient step. No need to check for this condition. */ 130a7e14dcfSSatish Balay 131a7e14dcfSSatish Balay lmP->bfgs = 1; 132a7e14dcfSSatish Balay ++lmP->sgrad; 133a7e14dcfSSatish Balay stepType = LMVM_SCALED_GRADIENT; 134a7e14dcfSSatish Balay break; 135a7e14dcfSSatish Balay 136a7e14dcfSSatish Balay case LMVM_SCALED_GRADIENT: 137a7e14dcfSSatish Balay /* The scaled gradient step did not produce a new iterate; 138a7e14dcfSSatish Balay attempt to use the gradient direction. 139a7e14dcfSSatish Balay Need to make sure we are not using a different diagonal scaling */ 140a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(lmP->M, 1.0);CHKERRQ(ierr); 141a7e14dcfSSatish Balay ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 142a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 143a7e14dcfSSatish Balay ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 144a7e14dcfSSatish Balay 145a7e14dcfSSatish Balay lmP->bfgs = 1; 146a7e14dcfSSatish Balay ++lmP->grad; 147a7e14dcfSSatish Balay stepType = LMVM_GRADIENT; 148a7e14dcfSSatish Balay break; 149a7e14dcfSSatish Balay } 150a7e14dcfSSatish Balay ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 151a7e14dcfSSatish Balay 152a7e14dcfSSatish Balay /* Perform the linesearch */ 153a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status);CHKERRQ(ierr); 154a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 155a7e14dcfSSatish Balay } 156a7e14dcfSSatish Balay 15787f595a5SBarry Smith if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 158a7e14dcfSSatish Balay /* Failed to find an improving point */ 159a7e14dcfSSatish Balay f = fold; 160a7e14dcfSSatish Balay ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 161a7e14dcfSSatish Balay ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 162a7e14dcfSSatish Balay step = 0.0; 163a7e14dcfSSatish Balay reason = TAO_DIVERGED_LS_FAILURE; 164a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 165a7e14dcfSSatish Balay } 166a7e14dcfSSatish Balay /* Check for termination */ 167a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 168a7e14dcfSSatish Balay iter++; 169a7e14dcfSSatish Balay ierr = TaoMonitor(tao,iter,f,gnorm,0.0,step,&reason);CHKERRQ(ierr); 170a7e14dcfSSatish Balay } 171a7e14dcfSSatish Balay PetscFunctionReturn(0); 172a7e14dcfSSatish Balay } 17387f595a5SBarry Smith 174a7e14dcfSSatish Balay #undef __FUNCT__ 175a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_LMVM" 176441846f8SBarry Smith static PetscErrorCode TaoSetUp_LMVM(Tao tao) 177a7e14dcfSSatish Balay { 178a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 179a7e14dcfSSatish Balay PetscInt n,N; 180a7e14dcfSSatish Balay PetscErrorCode ierr; 181a7e14dcfSSatish Balay 182a7e14dcfSSatish Balay PetscFunctionBegin; 183a7e14dcfSSatish Balay /* Existence of tao->solution checked in TaoSetUp() */ 184a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); } 185a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); } 186a7e14dcfSSatish Balay if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr); } 187a7e14dcfSSatish Balay if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr); } 188a7e14dcfSSatish Balay if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr); } 189a7e14dcfSSatish Balay 190a7e14dcfSSatish Balay /* Create matrix for the limited memory approximation */ 191a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 192a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 193a7e14dcfSSatish Balay ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&lmP->M);CHKERRQ(ierr); 194a7e14dcfSSatish Balay ierr = MatLMVMAllocateVectors(lmP->M,tao->solution);CHKERRQ(ierr); 195a7e14dcfSSatish Balay PetscFunctionReturn(0); 196a7e14dcfSSatish Balay } 197a7e14dcfSSatish Balay 198a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 199a7e14dcfSSatish Balay #undef __FUNCT__ 200a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_LMVM" 201441846f8SBarry Smith static PetscErrorCode TaoDestroy_LMVM(Tao tao) 202a7e14dcfSSatish Balay { 203a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 204a7e14dcfSSatish Balay PetscErrorCode ierr; 205a7e14dcfSSatish Balay 206a7e14dcfSSatish Balay PetscFunctionBegin; 207a7e14dcfSSatish Balay if (tao->setupcalled) { 208a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr); 209a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr); 210a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->D);CHKERRQ(ierr); 211a7e14dcfSSatish Balay ierr = MatDestroy(&lmP->M);CHKERRQ(ierr); 212a7e14dcfSSatish Balay } 213a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 214a7e14dcfSSatish Balay PetscFunctionReturn(0); 215a7e14dcfSSatish Balay } 216a7e14dcfSSatish Balay 217a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 218a7e14dcfSSatish Balay #undef __FUNCT__ 219a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_LMVM" 220441846f8SBarry Smith static PetscErrorCode TaoSetFromOptions_LMVM(Tao tao) 221a7e14dcfSSatish Balay { 222a7e14dcfSSatish Balay PetscErrorCode ierr; 223a7e14dcfSSatish Balay 224a7e14dcfSSatish Balay PetscFunctionBegin; 225a7e14dcfSSatish Balay ierr = PetscOptionsHead("Limited-memory variable-metric method for unconstrained optimization");CHKERRQ(ierr); 226a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 227a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 228a7e14dcfSSatish Balay PetscFunctionReturn(0); 229a7e14dcfSSatish Balay PetscFunctionReturn(0); 230a7e14dcfSSatish Balay } 231a7e14dcfSSatish Balay 232a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 233a7e14dcfSSatish Balay #undef __FUNCT__ 234a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_LMVM" 235441846f8SBarry Smith static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer) 236a7e14dcfSSatish Balay { 237a7e14dcfSSatish Balay TAO_LMVM *lm = (TAO_LMVM *)tao->data; 238a7e14dcfSSatish Balay PetscBool isascii; 239a7e14dcfSSatish Balay PetscErrorCode ierr; 240a7e14dcfSSatish Balay 241a7e14dcfSSatish Balay PetscFunctionBegin; 242a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 243a7e14dcfSSatish Balay if (isascii) { 244a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 245a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", lm->bfgs);CHKERRQ(ierr); 246a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", lm->sgrad);CHKERRQ(ierr); 247a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lm->grad);CHKERRQ(ierr); 248a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 249a7e14dcfSSatish Balay } 250a7e14dcfSSatish Balay PetscFunctionReturn(0); 251a7e14dcfSSatish Balay } 252a7e14dcfSSatish Balay 253a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 254a7e14dcfSSatish Balay 255a7e14dcfSSatish Balay EXTERN_C_BEGIN 256a7e14dcfSSatish Balay #undef __FUNCT__ 257a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_LMVM" 258441846f8SBarry Smith PetscErrorCode TaoCreate_LMVM(Tao tao) 259a7e14dcfSSatish Balay { 260a7e14dcfSSatish Balay TAO_LMVM *lmP; 261a7e14dcfSSatish Balay const char *morethuente_type = TAOLINESEARCH_MT; 262a7e14dcfSSatish Balay PetscErrorCode ierr; 263a7e14dcfSSatish Balay 264a7e14dcfSSatish Balay PetscFunctionBegin; 265a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_LMVM; 266a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_LMVM; 267a7e14dcfSSatish Balay tao->ops->view = TaoView_LMVM; 268a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_LMVM; 269a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_LMVM; 270a7e14dcfSSatish Balay 2713c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr); 272a7e14dcfSSatish Balay lmP->D = 0; 273a7e14dcfSSatish Balay lmP->M = 0; 274a7e14dcfSSatish Balay lmP->Xold = 0; 275a7e14dcfSSatish Balay lmP->Gold = 0; 276a7e14dcfSSatish Balay 277a7e14dcfSSatish Balay tao->data = (void*)lmP; 278a7e14dcfSSatish Balay tao->max_it = 2000; 279a7e14dcfSSatish Balay tao->max_funcs = 4000; 280a7e14dcfSSatish Balay tao->fatol = 1e-4; 281a7e14dcfSSatish Balay tao->frtol = 1e-4; 282a7e14dcfSSatish Balay 283a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 284a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 285441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 286a7e14dcfSSatish Balay PetscFunctionReturn(0); 287a7e14dcfSSatish Balay } 288a7e14dcfSSatish Balay EXTERN_C_END 289