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) { 199dcef436SStefano Zampini ierr = PetscInfo(tao,"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) { 43e1e80dc8SAlp Dener /* Call general purpose update function */ 44e1e80dc8SAlp Dener if (tao->ops->update) { 458fcddce6SStefano Zampini ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr); 46e1e80dc8SAlp Dener } 47e1e80dc8SAlp Dener 48a7e14dcfSSatish Balay /* Compute direction */ 49cd929ea3SAlp Dener if (lmP->H0) { 50cd929ea3SAlp Dener ierr = MatLMVMSetJ0(lmP->M, lmP->H0);CHKERRQ(ierr); 51cd929ea3SAlp Dener stepType = LMVM_STEP_BFGS; 52cd929ea3SAlp Dener } 53a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 549515a401SAlp Dener ierr = MatSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr); 55cd929ea3SAlp Dener ierr = MatLMVMGetUpdateCount(lmP->M, &nupdates); CHKERRQ(ierr); 56cd929ea3SAlp Dener if (nupdates > 0) stepType = LMVM_STEP_BFGS; 57a7e14dcfSSatish Balay 58a7e14dcfSSatish Balay /* Check for success (descent direction) */ 59a7e14dcfSSatish Balay ierr = VecDot(lmP->D, tao->gradient, &gdx);CHKERRQ(ierr); 60a7e14dcfSSatish Balay if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) { 61a7e14dcfSSatish Balay /* Step is not descent or direction produced not a number 62a7e14dcfSSatish Balay We can assert bfgsUpdates > 1 in this case because 63a7e14dcfSSatish Balay the first solve produces the scaled gradient direction, 64a7e14dcfSSatish Balay which is guaranteed to be descent 65a7e14dcfSSatish Balay 66a7e14dcfSSatish Balay Use steepest descent direction (scaled) 67a7e14dcfSSatish Balay */ 68a7e14dcfSSatish Balay 69cd929ea3SAlp Dener ierr = MatLMVMReset(lmP->M, PETSC_FALSE);CHKERRQ(ierr); 700ad3a497SAlp Dener ierr = MatLMVMClearJ0(lmP->M);CHKERRQ(ierr); 71a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 729515a401SAlp Dener ierr = MatSolve(lmP->M,tao->gradient, lmP->D);CHKERRQ(ierr); 73a7e14dcfSSatish Balay 74a7e14dcfSSatish Balay /* On a reset, the direction cannot be not a number; it is a 75a7e14dcfSSatish Balay scaled gradient step. No need to check for this condition. */ 76cd929ea3SAlp Dener stepType = LMVM_STEP_GRAD; 77a7e14dcfSSatish Balay } 78a7e14dcfSSatish Balay ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 79a7e14dcfSSatish Balay 80a7e14dcfSSatish Balay /* Perform the linesearch */ 81a7e14dcfSSatish Balay fold = f; 82a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr); 83a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr); 84a7e14dcfSSatish Balay 85a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status);CHKERRQ(ierr); 86a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 87a7e14dcfSSatish Balay 88cd929ea3SAlp Dener if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_STEP_GRAD)) { 89a7e14dcfSSatish Balay /* Reset factors and use scaled gradient step */ 90a7e14dcfSSatish Balay f = fold; 91a7e14dcfSSatish Balay ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 92a7e14dcfSSatish Balay ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 93a7e14dcfSSatish Balay 94a7e14dcfSSatish Balay /* Failed to obtain acceptable iterate with BFGS step */ 95a7e14dcfSSatish Balay /* Attempt to use the scaled gradient direction */ 96a7e14dcfSSatish Balay 97cd929ea3SAlp Dener ierr = MatLMVMReset(lmP->M, PETSC_FALSE);CHKERRQ(ierr); 980ad3a497SAlp Dener ierr = MatLMVMClearJ0(lmP->M);CHKERRQ(ierr); 99a7e14dcfSSatish Balay ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 1009515a401SAlp Dener ierr = MatSolve(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 101a7e14dcfSSatish Balay 102a7e14dcfSSatish Balay /* On a reset, the direction cannot be not a number; it is a 103a7e14dcfSSatish Balay scaled gradient step. No need to check for this condition. */ 104cd929ea3SAlp Dener stepType = LMVM_STEP_GRAD; 105a7e14dcfSSatish Balay ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 106a7e14dcfSSatish Balay 107a7e14dcfSSatish Balay /* Perform the linesearch */ 108a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status);CHKERRQ(ierr); 109a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 110a7e14dcfSSatish Balay } 111a7e14dcfSSatish Balay 11287f595a5SBarry Smith if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 113a7e14dcfSSatish Balay /* Failed to find an improving point */ 114a7e14dcfSSatish Balay f = fold; 115a7e14dcfSSatish Balay ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 116a7e14dcfSSatish Balay ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 117a7e14dcfSSatish Balay step = 0.0; 118a7e14dcfSSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 119cd929ea3SAlp Dener } else { 120cd929ea3SAlp Dener /* LS found valid step, so tally up step type */ 121cd929ea3SAlp Dener switch (stepType) { 122cd929ea3SAlp Dener case LMVM_STEP_BFGS: 123cd929ea3SAlp Dener ++lmP->bfgs; 124cd929ea3SAlp Dener break; 125cd929ea3SAlp Dener case LMVM_STEP_GRAD: 126cd929ea3SAlp Dener ++lmP->grad; 127cd929ea3SAlp Dener break; 128cd929ea3SAlp Dener default: 129cd929ea3SAlp Dener break; 130cd929ea3SAlp Dener } 131cd929ea3SAlp Dener /* Compute new gradient norm */ 132cd929ea3SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 133a7e14dcfSSatish Balay } 134a9603a14SPatrick Farrell 135cd929ea3SAlp Dener /* Check convergence */ 1368931d482SJason Sarich tao->niter++; 1373ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 1383ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 1393ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 140a7e14dcfSSatish Balay } 141a7e14dcfSSatish Balay PetscFunctionReturn(0); 142a7e14dcfSSatish Balay } 14387f595a5SBarry Smith 144441846f8SBarry Smith static PetscErrorCode TaoSetUp_LMVM(Tao tao) 145a7e14dcfSSatish Balay { 146a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 147a7e14dcfSSatish Balay PetscInt n,N; 148a7e14dcfSSatish Balay PetscErrorCode ierr; 149d5ae2380SAlp Dener PetscBool is_spd; 150a7e14dcfSSatish Balay 151a7e14dcfSSatish Balay PetscFunctionBegin; 152a7e14dcfSSatish Balay /* Existence of tao->solution checked in TaoSetUp() */ 153a7e14dcfSSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); } 154a7e14dcfSSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); } 155a7e14dcfSSatish Balay if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr); } 156a7e14dcfSSatish Balay if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr); } 157a7e14dcfSSatish Balay if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr); } 158a7e14dcfSSatish Balay 159a7e14dcfSSatish Balay /* Create matrix for the limited memory approximation */ 160a7e14dcfSSatish Balay ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 161a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 162cd929ea3SAlp Dener ierr = MatSetSizes(lmP->M, n, n, N, N);CHKERRQ(ierr); 163cd929ea3SAlp Dener ierr = MatLMVMAllocate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 1640ad3a497SAlp Dener ierr = MatGetOption(lmP->M, MAT_SPD, &is_spd);CHKERRQ(ierr); 1650ad3a497SAlp Dener if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix is not symmetric positive-definite."); 166a9603a14SPatrick Farrell 167a9603a14SPatrick Farrell /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 168a9603a14SPatrick Farrell if (lmP->H0) { 169cd929ea3SAlp Dener ierr = MatLMVMSetJ0(lmP->M, lmP->H0);CHKERRQ(ierr); 170a9603a14SPatrick Farrell } 171a9603a14SPatrick Farrell 172a7e14dcfSSatish Balay PetscFunctionReturn(0); 173a7e14dcfSSatish Balay } 174a7e14dcfSSatish Balay 175a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 176441846f8SBarry Smith static PetscErrorCode TaoDestroy_LMVM(Tao tao) 177a7e14dcfSSatish Balay { 178a7e14dcfSSatish Balay TAO_LMVM *lmP = (TAO_LMVM *)tao->data; 179a7e14dcfSSatish Balay PetscErrorCode ierr; 180a7e14dcfSSatish Balay 181a7e14dcfSSatish Balay PetscFunctionBegin; 182a7e14dcfSSatish Balay if (tao->setupcalled) { 183a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr); 184a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr); 185a7e14dcfSSatish Balay ierr = VecDestroy(&lmP->D);CHKERRQ(ierr); 186a7e14dcfSSatish Balay } 187cd929ea3SAlp Dener ierr = MatDestroy(&lmP->M);CHKERRQ(ierr); 188a9603a14SPatrick Farrell if (lmP->H0) { 189a9603a14SPatrick Farrell ierr = PetscObjectDereference((PetscObject)lmP->H0);CHKERRQ(ierr); 190a9603a14SPatrick Farrell } 191a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 192a9603a14SPatrick Farrell 193a7e14dcfSSatish Balay PetscFunctionReturn(0); 194a7e14dcfSSatish Balay } 195a7e14dcfSSatish Balay 196a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 1974416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptionItems *PetscOptionsObject,Tao tao) 198a7e14dcfSSatish Balay { 199cd929ea3SAlp Dener TAO_LMVM *lm = (TAO_LMVM *)tao->data; 200a7e14dcfSSatish Balay PetscErrorCode ierr; 201a7e14dcfSSatish Balay 202a7e14dcfSSatish Balay PetscFunctionBegin; 2031a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");CHKERRQ(ierr); 204cd929ea3SAlp Dener ierr = PetscOptionsBool("-tao_lmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",lm->recycle,&lm->recycle,NULL);CHKERRQ(ierr); 205114d2d62SAlp Dener ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 206cd929ea3SAlp Dener ierr = MatSetFromOptions(lm->M);CHKERRQ(ierr); 207288b7216SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 208a7e14dcfSSatish Balay PetscFunctionReturn(0); 209a7e14dcfSSatish Balay } 210a7e14dcfSSatish Balay 211a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 212441846f8SBarry Smith static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer) 213a7e14dcfSSatish Balay { 214a7e14dcfSSatish Balay TAO_LMVM *lm = (TAO_LMVM *)tao->data; 215cd929ea3SAlp Dener PetscBool isascii; 2164d6623b4SAlp Dener PetscInt recycled_its; 217a7e14dcfSSatish Balay PetscErrorCode ierr; 218a7e14dcfSSatish Balay 219a7e14dcfSSatish Balay PetscFunctionBegin; 220a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 221a7e14dcfSSatish Balay if (isascii) { 222a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer, " Gradient steps: %D\n", lm->grad);CHKERRQ(ierr); 223cd929ea3SAlp Dener if (lm->recycle) { 224288b7216SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " Recycle: on\n");CHKERRQ(ierr); 225cd929ea3SAlp Dener recycled_its = lm->bfgs + lm->grad; 2264d6623b4SAlp Dener ierr = PetscViewerASCIIPrintf(viewer, " Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr); 227a0bfee83SAlp Dener } 228a7e14dcfSSatish Balay } 229a7e14dcfSSatish Balay PetscFunctionReturn(0); 230a7e14dcfSSatish Balay } 231a7e14dcfSSatish Balay 232a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 233a7e14dcfSSatish Balay 2344aa34175SJason Sarich /*MC 2354aa34175SJason Sarich TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton 2364aa34175SJason Sarich optimization solver for unconstrained minimization. It solves 2374aa34175SJason Sarich the Newton step 2384aa34175SJason Sarich Hkdk = - gk 2394aa34175SJason Sarich 2404aa34175SJason Sarich using an approximation Bk in place of Hk, where Bk is composed using 2414aa34175SJason Sarich the BFGS update formula. A More-Thuente line search is then used 2424aa34175SJason Sarich to computed the steplength in the dk direction 2430ad3a497SAlp Dener 2444aa34175SJason Sarich Options Database Keys: 245*a2b725a8SWilliam Gropp + -tao_lmvm_recycle - enable recycling LMVM updates between TaoSolve() calls 246*a2b725a8SWilliam Gropp - -tao_lmvm_no_scale - (developer) disables diagonal Broyden scaling on the LMVM approximation 2474aa34175SJason Sarich 2481eb8069cSJason Sarich Level: beginner 2494aa34175SJason Sarich M*/ 2504aa34175SJason Sarich 251728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao) 252a7e14dcfSSatish Balay { 253a7e14dcfSSatish Balay TAO_LMVM *lmP; 2548caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 255a7e14dcfSSatish Balay PetscErrorCode ierr; 256a7e14dcfSSatish Balay 257a7e14dcfSSatish Balay PetscFunctionBegin; 258a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_LMVM; 259a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_LMVM; 260a7e14dcfSSatish Balay tao->ops->view = TaoView_LMVM; 261a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_LMVM; 262a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_LMVM; 263a7e14dcfSSatish Balay 2643c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr); 265a7e14dcfSSatish Balay lmP->D = 0; 266a7e14dcfSSatish Balay lmP->M = 0; 267a7e14dcfSSatish Balay lmP->Xold = 0; 268a7e14dcfSSatish Balay lmP->Gold = 0; 269a9603a14SPatrick Farrell lmP->H0 = NULL; 270cd929ea3SAlp Dener lmP->recycle = PETSC_FALSE; 271a7e14dcfSSatish Balay 272a7e14dcfSSatish Balay tao->data = (void*)lmP; 2736552cf8aSJason Sarich /* Override default settings (unless already changed) */ 2746552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 2756552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 276a7e14dcfSSatish Balay 277a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 27863b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 279a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 280441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 2815d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 282cd929ea3SAlp Dener 283cd929ea3SAlp Dener ierr = KSPInitializePackage();CHKERRQ(ierr); 284cd929ea3SAlp Dener ierr = MatCreate(((PetscObject)tao)->comm, &lmP->M);CHKERRQ(ierr); 285cd929ea3SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)lmP->M, (PetscObject)tao, 1);CHKERRQ(ierr); 28678e4361aSAlp Dener ierr = MatSetType(lmP->M, MATLMVMBFGS);CHKERRQ(ierr); 287cd929ea3SAlp Dener ierr = MatSetOptionsPrefix(lmP->M, "tao_lmvm_");CHKERRQ(ierr); 288a7e14dcfSSatish Balay PetscFunctionReturn(0); 289a7e14dcfSSatish Balay } 290