xref: /petsc/src/tao/unconstrained/impls/lmvm/lmvm.c (revision 03047865b8d8757cf1cf9cda45785c1537b01dc1)
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 
TaoSolve_LMVM(Tao tao)7d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_LMVM(Tao tao)
8d71ae5a4SJacob Faibussowitsch {
9a7e14dcfSSatish Balay   TAO_LMVM                    *lmP = (TAO_LMVM *)tao->data;
10a7e14dcfSSatish Balay   PetscReal                    f, fold, gdx, gnorm;
11a7e14dcfSSatish Balay   PetscReal                    step      = 1.0;
12cd929ea3SAlp Dener   PetscInt                     stepType  = LMVM_STEP_GRAD, nupdates;
13e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
14a7e14dcfSSatish Balay 
15a7e14dcfSSatish Balay   PetscFunctionBegin;
1648a46eb9SPierre Jolivet   if (tao->XL || tao->XU || tao->ops->computebounds) PetscCall(PetscInfo(tao, "WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n"));
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay   /*  Check convergence criteria */
199566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
209566063dSJacob Faibussowitsch   PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
21a9603a14SPatrick Farrell 
22*76c63389SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
23a7e14dcfSSatish Balay 
243ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
259566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
269566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
27dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
283ba16761SJacob Faibussowitsch   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
29a7e14dcfSSatish Balay 
30a7e14dcfSSatish Balay   /*  Set counter for gradient/reset steps */
31cd929ea3SAlp Dener   if (!lmP->recycle) {
32a7e14dcfSSatish Balay     lmP->bfgs = 0;
33a7e14dcfSSatish Balay     lmP->grad = 0;
349566063dSJacob Faibussowitsch     PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
35de6ffafeSAlp Dener   }
36a7e14dcfSSatish Balay 
37a7e14dcfSSatish Balay   /*  Have not converged; continue with Newton method */
383ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
39e1e80dc8SAlp Dener     /* Call general purpose update function */
40270bebe6SStefano Zampini     if (tao->ops->update) {
41270bebe6SStefano Zampini       PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
42270bebe6SStefano Zampini       PetscCall(TaoComputeObjective(tao, tao->solution, &f));
43270bebe6SStefano Zampini     }
44e1e80dc8SAlp Dener 
45a7e14dcfSSatish Balay     /*  Compute direction */
46cd929ea3SAlp Dener     if (lmP->H0) {
479566063dSJacob Faibussowitsch       PetscCall(MatLMVMSetJ0(lmP->M, lmP->H0));
48cd929ea3SAlp Dener       stepType = LMVM_STEP_BFGS;
49cd929ea3SAlp Dener     }
509566063dSJacob Faibussowitsch     PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
519566063dSJacob Faibussowitsch     PetscCall(MatSolve(lmP->M, tao->gradient, lmP->D));
529566063dSJacob Faibussowitsch     PetscCall(MatLMVMGetUpdateCount(lmP->M, &nupdates));
53cd929ea3SAlp Dener     if (nupdates > 0) stepType = LMVM_STEP_BFGS;
54a7e14dcfSSatish Balay 
55a7e14dcfSSatish Balay     /*  Check for success (descent direction) */
56bbb72809SHansol Suh     PetscCall(VecDotRealPart(lmP->D, tao->gradient, &gdx));
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 
669566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
679566063dSJacob Faibussowitsch       PetscCall(MatLMVMClearJ0(lmP->M));
689566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
699566063dSJacob Faibussowitsch       PetscCall(MatSolve(lmP->M, tao->gradient, lmP->D));
70a7e14dcfSSatish Balay 
71a7e14dcfSSatish Balay       /* On a reset, the direction cannot be not a number; it is a
72a7e14dcfSSatish Balay          scaled gradient step.  No need to check for this condition. */
73cd929ea3SAlp Dener       stepType = LMVM_STEP_GRAD;
74a7e14dcfSSatish Balay     }
759566063dSJacob Faibussowitsch     PetscCall(VecScale(lmP->D, -1.0));
76a7e14dcfSSatish Balay 
77a7e14dcfSSatish Balay     /*  Perform the linesearch */
78a7e14dcfSSatish Balay     fold = f;
799566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, lmP->Xold));
809566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, lmP->Gold));
81a7e14dcfSSatish Balay 
829566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status));
839566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
84a7e14dcfSSatish Balay 
85cd929ea3SAlp Dener     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_STEP_GRAD)) {
86a7e14dcfSSatish Balay       /*  Reset factors and use scaled gradient step */
87a7e14dcfSSatish Balay       f = fold;
889566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Xold, tao->solution));
899566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Gold, tao->gradient));
90a7e14dcfSSatish Balay 
91a7e14dcfSSatish Balay       /*  Failed to obtain acceptable iterate with BFGS step */
92a7e14dcfSSatish Balay       /*  Attempt to use the scaled gradient direction */
93a7e14dcfSSatish Balay 
949566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
959566063dSJacob Faibussowitsch       PetscCall(MatLMVMClearJ0(lmP->M));
969566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
979566063dSJacob Faibussowitsch       PetscCall(MatSolve(lmP->M, tao->solution, tao->gradient));
98a7e14dcfSSatish Balay 
99a7e14dcfSSatish Balay       /* On a reset, the direction cannot be not a number; it is a
100a7e14dcfSSatish Balay           scaled gradient step.  No need to check for this condition. */
101cd929ea3SAlp Dener       stepType = LMVM_STEP_GRAD;
1029566063dSJacob Faibussowitsch       PetscCall(VecScale(lmP->D, -1.0));
103a7e14dcfSSatish Balay 
104a7e14dcfSSatish Balay       /*  Perform the linesearch */
1059566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status));
1069566063dSJacob Faibussowitsch       PetscCall(TaoAddLineSearchCounts(tao));
107a7e14dcfSSatish Balay     }
108a7e14dcfSSatish Balay 
10987f595a5SBarry Smith     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
110a7e14dcfSSatish Balay       /*  Failed to find an improving point */
111a7e14dcfSSatish Balay       f = fold;
1129566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Xold, tao->solution));
1139566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Gold, tao->gradient));
114a7e14dcfSSatish Balay       step        = 0.0;
115a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
116cd929ea3SAlp Dener     } else {
117cd929ea3SAlp Dener       /* LS found valid step, so tally up step type */
118cd929ea3SAlp Dener       switch (stepType) {
119d71ae5a4SJacob Faibussowitsch       case LMVM_STEP_BFGS:
120d71ae5a4SJacob Faibussowitsch         ++lmP->bfgs;
121d71ae5a4SJacob Faibussowitsch         break;
122d71ae5a4SJacob Faibussowitsch       case LMVM_STEP_GRAD:
123d71ae5a4SJacob Faibussowitsch         ++lmP->grad;
124d71ae5a4SJacob Faibussowitsch         break;
125d71ae5a4SJacob Faibussowitsch       default:
126d71ae5a4SJacob Faibussowitsch         break;
127cd929ea3SAlp Dener       }
128cd929ea3SAlp Dener       /*  Compute new gradient norm */
1299566063dSJacob Faibussowitsch       PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
130a7e14dcfSSatish Balay     }
131a9603a14SPatrick Farrell 
132cd929ea3SAlp Dener     /* Check convergence */
1338931d482SJason Sarich     tao->niter++;
1349566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
1359566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
136dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
137a7e14dcfSSatish Balay   }
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
139a7e14dcfSSatish Balay }
14087f595a5SBarry Smith 
TaoSetUp_LMVM(Tao tao)141d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_LMVM(Tao tao)
142d71ae5a4SJacob Faibussowitsch {
143a7e14dcfSSatish Balay   TAO_LMVM *lmP = (TAO_LMVM *)tao->data;
144a7e14dcfSSatish Balay   PetscInt  n, N;
145b94d7dedSBarry Smith   PetscBool is_set, is_spd;
146a7e14dcfSSatish Balay 
147a7e14dcfSSatish Balay   PetscFunctionBegin;
148a7e14dcfSSatish Balay   /* Existence of tao->solution checked in TaoSetUp() */
1499566063dSJacob Faibussowitsch   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
1509566063dSJacob Faibussowitsch   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
1519566063dSJacob Faibussowitsch   if (!lmP->D) PetscCall(VecDuplicate(tao->solution, &lmP->D));
1529566063dSJacob Faibussowitsch   if (!lmP->Xold) PetscCall(VecDuplicate(tao->solution, &lmP->Xold));
1539566063dSJacob Faibussowitsch   if (!lmP->Gold) PetscCall(VecDuplicate(tao->solution, &lmP->Gold));
154a7e14dcfSSatish Balay 
155a7e14dcfSSatish Balay   /*  Create matrix for the limited memory approximation */
1569566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(tao->solution, &n));
1579566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tao->solution, &N));
1589566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(lmP->M, n, n, N, N));
1599566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(lmP->M, tao->solution, tao->gradient));
160b94d7dedSBarry Smith   PetscCall(MatIsSPDKnown(lmP->M, &is_set, &is_spd));
161b94d7dedSBarry Smith   PetscCheck(is_set && is_spd, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix is not symmetric positive-definite.");
162a9603a14SPatrick Farrell 
163a9603a14SPatrick Farrell   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
1641baa6e33SBarry Smith   if (lmP->H0) PetscCall(MatLMVMSetJ0(lmP->M, lmP->H0));
1653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
166a7e14dcfSSatish Balay }
167a7e14dcfSSatish Balay 
TaoDestroy_LMVM(Tao tao)168d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_LMVM(Tao tao)
169d71ae5a4SJacob Faibussowitsch {
170a7e14dcfSSatish Balay   TAO_LMVM *lmP = (TAO_LMVM *)tao->data;
171a7e14dcfSSatish Balay 
172a7e14dcfSSatish Balay   PetscFunctionBegin;
173a7e14dcfSSatish Balay   if (tao->setupcalled) {
1749566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->Xold));
1759566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->Gold));
1769566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->D));
177a7e14dcfSSatish Balay   }
1789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lmP->M));
1791baa6e33SBarry Smith   if (lmP->H0) PetscCall(PetscObjectDereference((PetscObject)lmP->H0));
1809566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
182a7e14dcfSSatish Balay }
183a7e14dcfSSatish Balay 
TaoSetFromOptions_LMVM(Tao tao,PetscOptionItems PetscOptionsObject)184ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_LMVM(Tao tao, PetscOptionItems PetscOptionsObject)
185d71ae5a4SJacob Faibussowitsch {
186cd929ea3SAlp Dener   TAO_LMVM *lm = (TAO_LMVM *)tao->data;
187a7e14dcfSSatish Balay 
188a7e14dcfSSatish Balay   PetscFunctionBegin;
189d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Limited-memory variable-metric method for unconstrained optimization");
1909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_lmvm_recycle", "enable recycling of the BFGS matrix between subsequent TaoSolve() calls", "", lm->recycle, &lm->recycle, NULL));
1919566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
1929566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(lm->M));
193d0609cedSBarry Smith   PetscOptionsHeadEnd();
1943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
195a7e14dcfSSatish Balay }
196a7e14dcfSSatish Balay 
TaoView_LMVM(Tao tao,PetscViewer viewer)197d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer)
198d71ae5a4SJacob Faibussowitsch {
199a7e14dcfSSatish Balay   TAO_LMVM *lm = (TAO_LMVM *)tao->data;
200cd929ea3SAlp Dener   PetscBool isascii;
2014d6623b4SAlp Dener   PetscInt  recycled_its;
202a7e14dcfSSatish Balay 
203a7e14dcfSSatish Balay   PetscFunctionBegin;
2049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
205a7e14dcfSSatish Balay   if (isascii) {
206ec9796c4SHansol Suh     PetscCall(PetscViewerASCIIPushTab(viewer));
20763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", lm->grad));
208cd929ea3SAlp Dener     if (lm->recycle) {
2099566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Recycle: on\n"));
210cd929ea3SAlp Dener       recycled_its = lm->bfgs + lm->grad;
21163a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Total recycled iterations: %" PetscInt_FMT "\n", recycled_its));
212a0bfee83SAlp Dener     }
213ec9796c4SHansol Suh     PetscCall(PetscViewerASCIIPrintf(viewer, "LMVM Matrix:\n"));
214ec9796c4SHansol Suh     PetscCall(PetscViewerASCIIPushTab(viewer));
215ec9796c4SHansol Suh     PetscCall(MatView(lm->M, viewer));
216ec9796c4SHansol Suh     PetscCall(PetscViewerASCIIPopTab(viewer));
217ec9796c4SHansol Suh     PetscCall(PetscViewerASCIIPopTab(viewer));
218a7e14dcfSSatish Balay   }
2193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
220a7e14dcfSSatish Balay }
221a7e14dcfSSatish Balay 
2224aa34175SJason Sarich /*MC
2234aa34175SJason Sarich   TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton
2244aa34175SJason Sarich   optimization solver for unconstrained minimization. It solves
2254aa34175SJason Sarich   the Newton step
2264aa34175SJason Sarich           Hkdk = - gk
2274aa34175SJason Sarich 
2284aa34175SJason Sarich   using an approximation Bk in place of Hk, where Bk is composed using
2294aa34175SJason Sarich   the BFGS update formula. A More-Thuente line search is then used
2304aa34175SJason Sarich   to computed the steplength in the dk direction
2310ad3a497SAlp Dener 
2324aa34175SJason Sarich   Options Database Keys:
233a2b725a8SWilliam Gropp +   -tao_lmvm_recycle - enable recycling LMVM updates between TaoSolve() calls
234a2b725a8SWilliam Gropp -   -tao_lmvm_no_scale - (developer) disables diagonal Broyden scaling on the LMVM approximation
2354aa34175SJason Sarich 
2361eb8069cSJason Sarich   Level: beginner
2374aa34175SJason Sarich M*/
2384aa34175SJason Sarich 
TaoCreate_LMVM(Tao tao)239d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao)
240d71ae5a4SJacob Faibussowitsch {
241a7e14dcfSSatish Balay   TAO_LMVM   *lmP;
2428caf6e8cSBarry Smith   const char *morethuente_type = TAOLINESEARCHMT;
243a7e14dcfSSatish Balay 
244a7e14dcfSSatish Balay   PetscFunctionBegin;
245a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetUp_LMVM;
246a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_LMVM;
247a7e14dcfSSatish Balay   tao->ops->view           = TaoView_LMVM;
248a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_LMVM;
249a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_LMVM;
250a7e14dcfSSatish Balay 
2514dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lmP));
25283c8fe1dSLisandro Dalcin   lmP->D       = NULL;
25383c8fe1dSLisandro Dalcin   lmP->M       = NULL;
25483c8fe1dSLisandro Dalcin   lmP->Xold    = NULL;
25583c8fe1dSLisandro Dalcin   lmP->Gold    = NULL;
256a9603a14SPatrick Farrell   lmP->H0      = NULL;
257cd929ea3SAlp Dener   lmP->recycle = PETSC_FALSE;
258a7e14dcfSSatish Balay 
259a7e14dcfSSatish Balay   tao->data = (void *)lmP;
2606552cf8aSJason Sarich   /* Override default settings (unless already changed) */
261606f75f6SBarry Smith   PetscCall(TaoParametersInitialize(tao));
262606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_it, 2000);
263606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_funcs, 4000);
264a7e14dcfSSatish Balay 
2659566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
2669566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
2679566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
2689566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
2699566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
270cd929ea3SAlp Dener 
2719566063dSJacob Faibussowitsch   PetscCall(KSPInitializePackage());
2729566063dSJacob Faibussowitsch   PetscCall(MatCreate(((PetscObject)tao)->comm, &lmP->M));
2739566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)lmP->M, (PetscObject)tao, 1));
2749566063dSJacob Faibussowitsch   PetscCall(MatSetType(lmP->M, MATLMVMBFGS));
2759566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(lmP->M, "tao_lmvm_"));
2763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
277a7e14dcfSSatish Balay }
278