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 9441846f8SBarry Smith static PetscErrorCode TaoSolve_LMVM(Tao tao) 10a7e14dcfSSatish Balay { 11a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 12a7e14dcfSSatish Balay PetscReal f, fold, gdx, gnorm; 13a7e14dcfSSatish Balay PetscReal step = 1.0; 14a7e14dcfSSatish Balay PetscReal delta; 15a7e14dcfSSatish Balay PetscErrorCode ierr; 16a7e14dcfSSatish Balay PetscInt stepType; 17e4cb33bbSBarry Smith TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 18e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 19a7e14dcfSSatish Balay 20a7e14dcfSSatish Balay PetscFunctionBegin; 21a7e14dcfSSatish Balay 22a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 23a7e14dcfSSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n");CHKERRQ(ierr); 24a7e14dcfSSatish Balay } 25a7e14dcfSSatish Balay 26a7e14dcfSSatish Balay /* Check convergence criteria */ 27a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 28a9603a14SPatrick Farrell ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 29a9603a14SPatrick Farrell 3087f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 31a7e14dcfSSatish Balay 328931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr); 3387f595a5SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 34a7e14dcfSSatish Balay 35a7e14dcfSSatish Balay /* Set initial scaling for the function */ 36a7e14dcfSSatish Balay if (f != 0.0) { 37a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 3887f595a5SBarry Smith } else { 39a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 40a7e14dcfSSatish Balay } 41a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(lmP->M,delta);CHKERRQ(ierr); 42a7e14dcfSSatish Balay 43a7e14dcfSSatish Balay /* Set counter for gradient/reset steps */ 44a7e14dcfSSatish Balay lmP->bfgs = 0; 45a7e14dcfSSatish Balay lmP->sgrad = 0; 46a7e14dcfSSatish Balay lmP->grad = 0; 47a7e14dcfSSatish Balay 48a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 49a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 50a7e14dcfSSatish Balay /* Compute direction */ 51a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 52a7e14dcfSSatish Balay ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 53a7e14dcfSSatish Balay ++lmP->bfgs; 54a7e14dcfSSatish Balay 55a7e14dcfSSatish Balay /* Check for success (descent direction) */ 56a7e14dcfSSatish Balay ierr = VecDot(lmP->D, tao->gradient, &gdx);CHKERRQ(ierr); 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 66a7e14dcfSSatish Balay ++lmP->grad; 67a7e14dcfSSatish Balay 68a7e14dcfSSatish Balay if (f != 0.0) { 69a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 7087f595a5SBarry Smith } else { 71a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 72a7e14dcfSSatish Balay } 73a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr); 74a7e14dcfSSatish Balay ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 75a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 76a7e14dcfSSatish Balay ierr = MatLMVMSolve(lmP->M,tao->gradient, lmP->D);CHKERRQ(ierr); 77a7e14dcfSSatish Balay 78a7e14dcfSSatish Balay /* On a reset, the direction cannot be not a number; it is a 79a7e14dcfSSatish Balay scaled gradient step. No need to check for this condition. */ 80a7e14dcfSSatish Balay 81a7e14dcfSSatish Balay lmP->bfgs = 1; 82a7e14dcfSSatish Balay ++lmP->sgrad; 83a7e14dcfSSatish Balay stepType = LMVM_SCALED_GRADIENT; 8487f595a5SBarry Smith } else { 85*a0bfee83SAlp Dener if (1 == lmP->bfgs && !lmP->recycle) { 86a7e14dcfSSatish Balay /* The first BFGS direction is always the scaled gradient */ 87a7e14dcfSSatish Balay ++lmP->sgrad; 88a7e14dcfSSatish Balay stepType = LMVM_SCALED_GRADIENT; 8987f595a5SBarry Smith } else { 90a7e14dcfSSatish Balay ++lmP->bfgs; 91a7e14dcfSSatish Balay stepType = LMVM_BFGS; 92a7e14dcfSSatish Balay } 93a7e14dcfSSatish Balay } 94a7e14dcfSSatish Balay ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 95a7e14dcfSSatish Balay 96a7e14dcfSSatish Balay /* Perform the linesearch */ 97a7e14dcfSSatish Balay fold = f; 98a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr); 99a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr); 100a7e14dcfSSatish Balay 101a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status);CHKERRQ(ierr); 102a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 103a7e14dcfSSatish Balay 10487f595a5SBarry Smith while (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_GRADIENT)) { 105a7e14dcfSSatish Balay /* Linesearch failed */ 106a7e14dcfSSatish Balay /* Reset factors and use scaled gradient step */ 107a7e14dcfSSatish Balay f = fold; 108a7e14dcfSSatish Balay ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 109a7e14dcfSSatish Balay ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 110a7e14dcfSSatish Balay 111a7e14dcfSSatish Balay switch(stepType) { 112a7e14dcfSSatish Balay case LMVM_BFGS: 113a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with BFGS step */ 114a7e14dcfSSatish Balay /* Attempt to use the scaled gradient direction */ 115a7e14dcfSSatish Balay 116a7e14dcfSSatish Balay if (f != 0.0) { 117a7e14dcfSSatish Balay delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 11887f595a5SBarry Smith } else { 119a7e14dcfSSatish Balay delta = 2.0 / (gnorm*gnorm); 120a7e14dcfSSatish Balay } 121a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr); 122a7e14dcfSSatish Balay ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 123a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 124a7e14dcfSSatish Balay ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 125a7e14dcfSSatish Balay 126a7e14dcfSSatish Balay /* On a reset, the direction cannot be not a number; it is a 127a7e14dcfSSatish Balay scaled gradient step. No need to check for this condition. */ 128a7e14dcfSSatish Balay 129a7e14dcfSSatish Balay lmP->bfgs = 1; 130a7e14dcfSSatish Balay ++lmP->sgrad; 131a7e14dcfSSatish Balay stepType = LMVM_SCALED_GRADIENT; 132a7e14dcfSSatish Balay break; 133a7e14dcfSSatish Balay 134a7e14dcfSSatish Balay case LMVM_SCALED_GRADIENT: 135a7e14dcfSSatish Balay /* The scaled gradient step did not produce a new iterate; 136a7e14dcfSSatish Balay attempt to use the gradient direction. 137a7e14dcfSSatish Balay Need to make sure we are not using a different diagonal scaling */ 138a7e14dcfSSatish Balay ierr = MatLMVMSetDelta(lmP->M, 1.0);CHKERRQ(ierr); 139a7e14dcfSSatish Balay ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 140a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 141a7e14dcfSSatish Balay ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 142a7e14dcfSSatish Balay 143a7e14dcfSSatish Balay lmP->bfgs = 1; 144a7e14dcfSSatish Balay ++lmP->grad; 145a7e14dcfSSatish Balay stepType = LMVM_GRADIENT; 146a7e14dcfSSatish Balay break; 147a7e14dcfSSatish Balay } 148a7e14dcfSSatish Balay ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 149a7e14dcfSSatish Balay 150a7e14dcfSSatish Balay /* Perform the linesearch */ 151a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status);CHKERRQ(ierr); 152a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 153a7e14dcfSSatish Balay } 154a7e14dcfSSatish Balay 15587f595a5SBarry Smith if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 156a7e14dcfSSatish Balay /* Failed to find an improving point */ 157a7e14dcfSSatish Balay f = fold; 158a7e14dcfSSatish Balay ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 159a7e14dcfSSatish Balay ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 160a7e14dcfSSatish Balay step = 0.0; 161a7e14dcfSSatish Balay reason = TAO_DIVERGED_LS_FAILURE; 162a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 163a7e14dcfSSatish Balay } 164a9603a14SPatrick Farrell 165a7e14dcfSSatish Balay /* Check for termination */ 166a9603a14SPatrick Farrell ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 167a9603a14SPatrick Farrell 1688931d482SJason Sarich tao->niter++; 1698931d482SJason Sarich ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step,&reason);CHKERRQ(ierr); 170a7e14dcfSSatish Balay } 171a7e14dcfSSatish Balay PetscFunctionReturn(0); 172a7e14dcfSSatish Balay } 17387f595a5SBarry Smith 174441846f8SBarry Smith static PetscErrorCode TaoSetUp_LMVM(Tao tao) 175a7e14dcfSSatish Balay { 176a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 177a7e14dcfSSatish Balay PetscInt n,N; 178a7e14dcfSSatish Balay PetscErrorCode ierr; 179a9603a14SPatrick Farrell KSP H0ksp; 180a7e14dcfSSatish Balay 181a7e14dcfSSatish Balay PetscFunctionBegin; 182a7e14dcfSSatish Balay /* Existence of tao->solution checked in TaoSetUp() */ 183a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); } 184a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); } 185a7e14dcfSSatish Balay if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr); } 186a7e14dcfSSatish Balay if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr); } 187a7e14dcfSSatish Balay if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr); } 188*a0bfee83SAlp Dener lmP->recycle = PETSC_FALSE; 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); 195a9603a14SPatrick Farrell 196a9603a14SPatrick Farrell /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 197a9603a14SPatrick Farrell if (lmP->H0) { 198a9603a14SPatrick Farrell const char *prefix; 199a9603a14SPatrick Farrell PC H0pc; 200a9603a14SPatrick Farrell 201a9603a14SPatrick Farrell ierr = MatLMVMSetH0(lmP->M, lmP->H0);CHKERRQ(ierr); 202a9603a14SPatrick Farrell ierr = MatLMVMGetH0KSP(lmP->M, &H0ksp);CHKERRQ(ierr); 203a9603a14SPatrick Farrell 204a9603a14SPatrick Farrell ierr = TaoGetOptionsPrefix(tao, &prefix);CHKERRQ(ierr); 205a9603a14SPatrick Farrell ierr = KSPSetOptionsPrefix(H0ksp, prefix);CHKERRQ(ierr); 206a9603a14SPatrick Farrell ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0ksp, "tao_h0_");CHKERRQ(ierr); 207a9603a14SPatrick Farrell ierr = KSPGetPC(H0ksp, &H0pc);CHKERRQ(ierr); 208a9603a14SPatrick Farrell ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0pc, "tao_h0_");CHKERRQ(ierr); 209a9603a14SPatrick Farrell 210a9603a14SPatrick Farrell ierr = KSPSetFromOptions(H0ksp);CHKERRQ(ierr); 211a9603a14SPatrick Farrell ierr = KSPSetUp(H0ksp);CHKERRQ(ierr); 212a9603a14SPatrick Farrell } 213a9603a14SPatrick Farrell 214a7e14dcfSSatish Balay PetscFunctionReturn(0); 215a7e14dcfSSatish Balay } 216a7e14dcfSSatish Balay 217a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 218441846f8SBarry Smith static PetscErrorCode TaoDestroy_LMVM(Tao tao) 219a7e14dcfSSatish Balay { 220a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 221a7e14dcfSSatish Balay PetscErrorCode ierr; 222a7e14dcfSSatish Balay 223a7e14dcfSSatish Balay PetscFunctionBegin; 224a7e14dcfSSatish Balay if (tao->setupcalled) { 225a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr); 226a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr); 227a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->D);CHKERRQ(ierr); 228a7e14dcfSSatish Balay ierr = MatDestroy(&lmP->M);CHKERRQ(ierr); 229a7e14dcfSSatish Balay } 230a9603a14SPatrick Farrell 231a9603a14SPatrick Farrell if (lmP->H0) { 232a9603a14SPatrick Farrell ierr = PetscObjectDereference((PetscObject)lmP->H0);CHKERRQ(ierr); 233a9603a14SPatrick Farrell } 234a9603a14SPatrick Farrell 235a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 236a9603a14SPatrick Farrell 237a7e14dcfSSatish Balay PetscFunctionReturn(0); 238a7e14dcfSSatish Balay } 239a7e14dcfSSatish Balay 240a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 2414416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptionItems *PetscOptionsObject,Tao tao) 242a7e14dcfSSatish Balay { 243a7e14dcfSSatish Balay PetscErrorCode ierr; 244*a0bfee83SAlp Dener TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 245a7e14dcfSSatish Balay 246a7e14dcfSSatish Balay PetscFunctionBegin; 2471a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");CHKERRQ(ierr); 248*a0bfee83SAlp Dener ierr = PetscOptionsBool("-tao_lmvm_recycle","Recycle the BFGS information between subsequent TaoSolve() calls.","TaoSolve_LMVM",lmP->recycle,&lmP->recycle,NULL);CHKERRQ(ierr); 249a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 250a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 251a7e14dcfSSatish Balay PetscFunctionReturn(0); 252a7e14dcfSSatish Balay } 253a7e14dcfSSatish Balay 254a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 255441846f8SBarry Smith static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer) 256a7e14dcfSSatish Balay { 257a7e14dcfSSatish Balay TAO_LMVM *lm = (TAO_LMVM *)tao->data; 258a7e14dcfSSatish Balay PetscBool isascii; 259a7e14dcfSSatish Balay PetscErrorCode ierr; 260a7e14dcfSSatish Balay 261a7e14dcfSSatish Balay PetscFunctionBegin; 262a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 263a7e14dcfSSatish Balay if (isascii) { 264a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 265a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", lm->bfgs);CHKERRQ(ierr); 266a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", lm->sgrad);CHKERRQ(ierr); 267a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lm->grad);CHKERRQ(ierr); 268*a0bfee83SAlp Dener if (lm->recycle) { 269*a0bfee83SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "recycle: on\n");CHKERRQ(ierr); 270*a0bfee83SAlp Dener } else { 271*a0bfee83SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "recycle: off\n");CHKERRQ(ierr); 272*a0bfee83SAlp Dener } 273a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 274a7e14dcfSSatish Balay } 275a7e14dcfSSatish Balay PetscFunctionReturn(0); 276a7e14dcfSSatish Balay } 277a7e14dcfSSatish Balay 278a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 279a7e14dcfSSatish Balay 2804aa34175SJason Sarich /*MC 2814aa34175SJason Sarich TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton 2824aa34175SJason Sarich optimization solver for unconstrained minimization. It solves 2834aa34175SJason Sarich the Newton step 2844aa34175SJason Sarich Hkdk = - gk 2854aa34175SJason Sarich 2864aa34175SJason Sarich using an approximation Bk in place of Hk, where Bk is composed using 2874aa34175SJason Sarich the BFGS update formula. A More-Thuente line search is then used 2884aa34175SJason Sarich to computed the steplength in the dk direction 2894aa34175SJason Sarich Options Database Keys: 2904aa34175SJason Sarich + -tao_lmm_vectors - number of vectors to use for approximation 2914aa34175SJason Sarich . -tao_lmm_scale_type - "none","scalar","broyden" 2924aa34175SJason Sarich . -tao_lmm_limit_type - "none","average","relative","absolute" 2934aa34175SJason Sarich . -tao_lmm_rescale_type - "none","scalar","gl" 2944aa34175SJason Sarich . -tao_lmm_limit_mu - mu limiting factor 2954aa34175SJason Sarich . -tao_lmm_limit_nu - nu limiting factor 2964aa34175SJason Sarich . -tao_lmm_delta_min - minimum delta value 2974aa34175SJason Sarich . -tao_lmm_delta_max - maximum delta value 2984aa34175SJason Sarich . -tao_lmm_broyden_phi - phi factor for Broyden scaling 2994aa34175SJason Sarich . -tao_lmm_scalar_alpha - alpha factor for scalar scaling 3004aa34175SJason Sarich . -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal 3014aa34175SJason Sarich . -tao_lmm_rescale_beta - beta factor for rescaling diagonal 3024aa34175SJason Sarich . -tao_lmm_scalar_history - amount of history for scalar scaling 3034aa34175SJason Sarich . -tao_lmm_rescale_history - amount of history for rescaling diagonal 3044aa34175SJason Sarich - -tao_lmm_eps - rejection tolerance 3054aa34175SJason Sarich 3061eb8069cSJason Sarich Level: beginner 3074aa34175SJason Sarich M*/ 3084aa34175SJason Sarich 309728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao) 310a7e14dcfSSatish Balay { 311a7e14dcfSSatish Balay TAO_LMVM *lmP; 3128caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 313a7e14dcfSSatish Balay PetscErrorCode ierr; 314a7e14dcfSSatish Balay 315a7e14dcfSSatish Balay PetscFunctionBegin; 316a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_LMVM; 317a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_LMVM; 318a7e14dcfSSatish Balay tao->ops->view = TaoView_LMVM; 319a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_LMVM; 320a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_LMVM; 321a7e14dcfSSatish Balay 3223c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr); 323a7e14dcfSSatish Balay lmP->D = 0; 324a7e14dcfSSatish Balay lmP->M = 0; 325a7e14dcfSSatish Balay lmP->Xold = 0; 326a7e14dcfSSatish Balay lmP->Gold = 0; 327a9603a14SPatrick Farrell lmP->H0 = NULL; 328a7e14dcfSSatish Balay 329a7e14dcfSSatish Balay tao->data = (void*)lmP; 3306552cf8aSJason Sarich /* Override default settings (unless already changed) */ 3316552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 3326552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 333a7e14dcfSSatish Balay 334a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 33563b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 336a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 337441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 3385d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 339a7e14dcfSSatish Balay PetscFunctionReturn(0); 340a7e14dcfSSatish Balay } 341