1ba92ff59SBarry Smith #include <petsctaolinesearch.h> 2aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/lmvm/lmvm.h> 3a7e14dcfSSatish Balay 4*cd929ea3SAlp Dener #define LMVM_STEP_BFGS 0 5*cd929ea3SAlp 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 PetscReal delta; 13a7e14dcfSSatish Balay PetscErrorCode ierr; 14*cd929ea3SAlp Dener PetscInt stepType = LMVM_STEP_GRAD, nupdates; 15e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 16a7e14dcfSSatish Balay 17a7e14dcfSSatish Balay PetscFunctionBegin; 18a7e14dcfSSatish Balay 19a7e14dcfSSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 20a7e14dcfSSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n");CHKERRQ(ierr); 21a7e14dcfSSatish Balay } 22a7e14dcfSSatish Balay 23a7e14dcfSSatish Balay /* Check convergence criteria */ 24a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 25a9603a14SPatrick Farrell ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 26a9603a14SPatrick Farrell 2787f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 28a7e14dcfSSatish Balay 293ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 303ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 313ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 323ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 333ecd9318SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 34a7e14dcfSSatish Balay 35a7e14dcfSSatish Balay /* Set initial scaling for the function */ 36*cd929ea3SAlp Dener if (!lmP->H0) { 37*cd929ea3SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm); 38*cd929ea3SAlp Dener ierr = MatLMVMSetJ0Scale(lmP->M, delta);CHKERRQ(ierr); 39a7e14dcfSSatish Balay } 40a7e14dcfSSatish Balay 41a7e14dcfSSatish Balay /* Set counter for gradient/reset steps */ 42*cd929ea3SAlp Dener if (!lmP->recycle) { 43a7e14dcfSSatish Balay lmP->bfgs = 0; 44a7e14dcfSSatish Balay lmP->grad = 0; 45*cd929ea3SAlp Dener ierr = MatLMVMReset(lmP->M, PETSC_FALSE); CHKERRQ(ierr); 46de6ffafeSAlp Dener } 47a7e14dcfSSatish Balay 48a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */ 493ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 50a7e14dcfSSatish Balay /* Compute direction */ 51*cd929ea3SAlp Dener if (lmP->H0) { 52*cd929ea3SAlp Dener ierr = MatLMVMSetJ0(lmP->M, lmP->H0);CHKERRQ(ierr); 53*cd929ea3SAlp Dener stepType = LMVM_STEP_BFGS; 54*cd929ea3SAlp Dener } 55a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 56*cd929ea3SAlp Dener ierr = MatSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 57*cd929ea3SAlp Dener ierr = MatLMVMGetUpdateCount(lmP->M, &nupdates); CHKERRQ(ierr); 58*cd929ea3SAlp Dener if (nupdates > 0) stepType = LMVM_STEP_BFGS; 59a7e14dcfSSatish Balay 60a7e14dcfSSatish Balay /* Check for success (descent direction) */ 61a7e14dcfSSatish Balay ierr = VecDot(lmP->D, tao->gradient, &gdx);CHKERRQ(ierr); 62a7e14dcfSSatish Balay if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) { 63a7e14dcfSSatish Balay /* Step is not descent or direction produced not a number 64a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 65a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 66a7e14dcfSSatish Balay which is guaranteed to be descent 67a7e14dcfSSatish Balay 68a7e14dcfSSatish Balay Use steepest descent direction (scaled) 69a7e14dcfSSatish Balay */ 70a7e14dcfSSatish Balay 71*cd929ea3SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm); 72*cd929ea3SAlp Dener ierr = MatLMVMSetJ0Scale(lmP->M, delta);CHKERRQ(ierr); 73*cd929ea3SAlp Dener ierr = MatLMVMReset(lmP->M, PETSC_FALSE);CHKERRQ(ierr); 74a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 75*cd929ea3SAlp Dener ierr = MatSolve(lmP->M,tao->gradient, lmP->D);CHKERRQ(ierr); 76a7e14dcfSSatish Balay 77a7e14dcfSSatish Balay /* On a reset, the direction cannot be not a number; it is a 78a7e14dcfSSatish Balay scaled gradient step. No need to check for this condition. */ 79*cd929ea3SAlp Dener stepType = LMVM_STEP_GRAD; 80a7e14dcfSSatish Balay } 81a7e14dcfSSatish Balay ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 82a7e14dcfSSatish Balay 83a7e14dcfSSatish Balay /* Perform the linesearch */ 84a7e14dcfSSatish Balay fold = f; 85a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr); 86a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr); 87a7e14dcfSSatish Balay 88a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status);CHKERRQ(ierr); 89a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 90a7e14dcfSSatish Balay 91*cd929ea3SAlp Dener if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_STEP_GRAD)) { 92a7e14dcfSSatish Balay /* Reset factors and use scaled gradient step */ 93a7e14dcfSSatish Balay f = fold; 94a7e14dcfSSatish Balay ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 95a7e14dcfSSatish Balay ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 96a7e14dcfSSatish Balay 97a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with BFGS step */ 98a7e14dcfSSatish Balay /* Attempt to use the scaled gradient direction */ 99a7e14dcfSSatish Balay 100*cd929ea3SAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm); 101*cd929ea3SAlp Dener ierr = MatLMVMSetJ0Scale(lmP->M, delta);CHKERRQ(ierr); 102*cd929ea3SAlp Dener ierr = MatLMVMReset(lmP->M, PETSC_FALSE);CHKERRQ(ierr); 103a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 104*cd929ea3SAlp Dener ierr = MatSolve(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 105a7e14dcfSSatish Balay 106a7e14dcfSSatish Balay /* On a reset, the direction cannot be not a number; it is a 107a7e14dcfSSatish Balay scaled gradient step. No need to check for this condition. */ 108*cd929ea3SAlp Dener stepType = LMVM_STEP_GRAD; 109a7e14dcfSSatish Balay ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 110a7e14dcfSSatish Balay 111a7e14dcfSSatish Balay /* Perform the linesearch */ 112a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status);CHKERRQ(ierr); 113a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 114a7e14dcfSSatish Balay } 115a7e14dcfSSatish Balay 11687f595a5SBarry Smith if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 117a7e14dcfSSatish Balay /* Failed to find an improving point */ 118a7e14dcfSSatish Balay f = fold; 119a7e14dcfSSatish Balay ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 120a7e14dcfSSatish Balay ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 121a7e14dcfSSatish Balay step = 0.0; 122a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 123*cd929ea3SAlp Dener } else { 124*cd929ea3SAlp Dener /* LS found valid step, so tally up step type */ 125*cd929ea3SAlp Dener switch (stepType) { 126*cd929ea3SAlp Dener case LMVM_STEP_BFGS: 127*cd929ea3SAlp Dener ++lmP->bfgs; 128*cd929ea3SAlp Dener break; 129*cd929ea3SAlp Dener case LMVM_STEP_GRAD: 130*cd929ea3SAlp Dener ++lmP->grad; 131*cd929ea3SAlp Dener break; 132*cd929ea3SAlp Dener default: 133*cd929ea3SAlp Dener break; 134*cd929ea3SAlp Dener } 135*cd929ea3SAlp Dener /* Compute new gradient norm */ 136*cd929ea3SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 137a7e14dcfSSatish Balay } 138a9603a14SPatrick Farrell 139*cd929ea3SAlp Dener /* Check convergence */ 1408931d482SJason Sarich tao->niter++; 1413ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 1423ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 1433ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 144a7e14dcfSSatish Balay } 145a7e14dcfSSatish Balay PetscFunctionReturn(0); 146a7e14dcfSSatish Balay } 14787f595a5SBarry Smith 148441846f8SBarry Smith static PetscErrorCode TaoSetUp_LMVM(Tao tao) 149a7e14dcfSSatish Balay { 150a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 151a7e14dcfSSatish Balay PetscInt n,N; 152a7e14dcfSSatish Balay PetscErrorCode ierr; 153a7e14dcfSSatish Balay 154a7e14dcfSSatish Balay PetscFunctionBegin; 155a7e14dcfSSatish Balay /* Existence of tao->solution checked in TaoSetUp() */ 156a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); } 157a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); } 158a7e14dcfSSatish Balay if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr); } 159a7e14dcfSSatish Balay if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr); } 160a7e14dcfSSatish Balay if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr); } 161a7e14dcfSSatish Balay 162a7e14dcfSSatish Balay /* Create matrix for the limited memory approximation */ 163a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 164a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 165*cd929ea3SAlp Dener ierr = MatSetSizes(lmP->M, n, n, N, N);CHKERRQ(ierr); 166*cd929ea3SAlp Dener ierr = MatLMVMAllocate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 167a9603a14SPatrick Farrell 168a9603a14SPatrick Farrell /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 169a9603a14SPatrick Farrell if (lmP->H0) { 170*cd929ea3SAlp Dener ierr = MatLMVMSetJ0(lmP->M, lmP->H0);CHKERRQ(ierr); 171a9603a14SPatrick Farrell } 172a9603a14SPatrick Farrell 173a7e14dcfSSatish Balay PetscFunctionReturn(0); 174a7e14dcfSSatish Balay } 175a7e14dcfSSatish Balay 176a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 177441846f8SBarry Smith static PetscErrorCode TaoDestroy_LMVM(Tao tao) 178a7e14dcfSSatish Balay { 179a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 180a7e14dcfSSatish Balay PetscErrorCode ierr; 181a7e14dcfSSatish Balay 182a7e14dcfSSatish Balay PetscFunctionBegin; 183a7e14dcfSSatish Balay if (tao->setupcalled) { 184a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr); 185a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr); 186a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->D);CHKERRQ(ierr); 187a7e14dcfSSatish Balay } 188*cd929ea3SAlp Dener ierr = MatDestroy(&lmP->M);CHKERRQ(ierr); 189a9603a14SPatrick Farrell if (lmP->H0) { 190a9603a14SPatrick Farrell ierr = PetscObjectDereference((PetscObject)lmP->H0);CHKERRQ(ierr); 191a9603a14SPatrick Farrell } 192a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 193a9603a14SPatrick Farrell 194a7e14dcfSSatish Balay PetscFunctionReturn(0); 195a7e14dcfSSatish Balay } 196a7e14dcfSSatish Balay 197a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 1984416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptionItems *PetscOptionsObject,Tao tao) 199a7e14dcfSSatish Balay { 200*cd929ea3SAlp Dener TAO_LMVM *lm = (TAO_LMVM *)tao->data; 201a7e14dcfSSatish Balay PetscErrorCode ierr; 202a7e14dcfSSatish Balay 203a7e14dcfSSatish Balay PetscFunctionBegin; 2041a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");CHKERRQ(ierr); 205*cd929ea3SAlp Dener ierr = PetscOptionsBool("-tao_lmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",lm->recycle,&lm->recycle,NULL);CHKERRQ(ierr); 206114d2d62SAlp Dener ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 207*cd929ea3SAlp Dener ierr = MatSetFromOptions(lm->M);CHKERRQ(ierr); 208288b7216SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 209a7e14dcfSSatish Balay PetscFunctionReturn(0); 210a7e14dcfSSatish Balay } 211a7e14dcfSSatish Balay 212a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 213441846f8SBarry Smith static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer) 214a7e14dcfSSatish Balay { 215a7e14dcfSSatish Balay TAO_LMVM *lm = (TAO_LMVM *)tao->data; 216*cd929ea3SAlp Dener PetscBool isascii; 2174d6623b4SAlp Dener PetscInt recycled_its; 218a7e14dcfSSatish Balay PetscErrorCode ierr; 219a7e14dcfSSatish Balay 220a7e14dcfSSatish Balay PetscFunctionBegin; 221a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 222a7e14dcfSSatish Balay if (isascii) { 223a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, " Gradient steps: %D\n", lm->grad);CHKERRQ(ierr); 224*cd929ea3SAlp Dener if (lm->recycle) { 225288b7216SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " Recycle: on\n");CHKERRQ(ierr); 226*cd929ea3SAlp Dener recycled_its = lm->bfgs + lm->grad; 2274d6623b4SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr); 228a0bfee83SAlp Dener } 229a7e14dcfSSatish Balay } 230a7e14dcfSSatish Balay PetscFunctionReturn(0); 231a7e14dcfSSatish Balay } 232a7e14dcfSSatish Balay 233a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 234a7e14dcfSSatish Balay 2354aa34175SJason Sarich /*MC 2364aa34175SJason Sarich TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton 2374aa34175SJason Sarich optimization solver for unconstrained minimization. It solves 2384aa34175SJason Sarich the Newton step 2394aa34175SJason Sarich Hkdk = - gk 2404aa34175SJason Sarich 2414aa34175SJason Sarich using an approximation Bk in place of Hk, where Bk is composed using 2424aa34175SJason Sarich the BFGS update formula. A More-Thuente line search is then used 2434aa34175SJason Sarich to computed the steplength in the dk direction 2444aa34175SJason Sarich Options Database Keys: 245*cd929ea3SAlp Dener . -tao_lmvm_recycle - enable recycling LMVM updates between TaoSolve() calls 2464aa34175SJason Sarich 2471eb8069cSJason Sarich Level: beginner 2484aa34175SJason Sarich M*/ 2494aa34175SJason Sarich 250728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao) 251a7e14dcfSSatish Balay { 252a7e14dcfSSatish Balay TAO_LMVM *lmP; 2538caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 254a7e14dcfSSatish Balay PetscErrorCode ierr; 255a7e14dcfSSatish Balay 256a7e14dcfSSatish Balay PetscFunctionBegin; 257a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_LMVM; 258a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_LMVM; 259a7e14dcfSSatish Balay tao->ops->view = TaoView_LMVM; 260a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_LMVM; 261a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_LMVM; 262a7e14dcfSSatish Balay 2633c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr); 264a7e14dcfSSatish Balay lmP->D = 0; 265a7e14dcfSSatish Balay lmP->M = 0; 266a7e14dcfSSatish Balay lmP->Xold = 0; 267a7e14dcfSSatish Balay lmP->Gold = 0; 268a9603a14SPatrick Farrell lmP->H0 = NULL; 269*cd929ea3SAlp Dener lmP->recycle = PETSC_FALSE; 270a7e14dcfSSatish Balay 271a7e14dcfSSatish Balay tao->data = (void*)lmP; 2726552cf8aSJason Sarich /* Override default settings (unless already changed) */ 2736552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 2746552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 275a7e14dcfSSatish Balay 276a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 27763b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 278a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 279441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 2805d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 281*cd929ea3SAlp Dener 282*cd929ea3SAlp Dener ierr = KSPInitializePackage();CHKERRQ(ierr); 283*cd929ea3SAlp Dener ierr = MatCreate(((PetscObject)tao)->comm, &lmP->M);CHKERRQ(ierr); 284*cd929ea3SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)lmP->M, (PetscObject)tao, 1);CHKERRQ(ierr); 285*cd929ea3SAlp Dener ierr = MatSetType(lmP->M, MATLBFGS);CHKERRQ(ierr); 286*cd929ea3SAlp Dener ierr = MatSetOptionsPrefix(lmP->M, "tao_lmvm_");CHKERRQ(ierr); 287a7e14dcfSSatish Balay PetscFunctionReturn(0); 288a7e14dcfSSatish Balay } 289