xref: /petsc/src/snes/utils/dmdasnes.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
16cab3a1bSJed Brown #include <petscdmda.h>          /*I "petscdmda.h" I*/
2af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
3af0996ceSBarry Smith #include <petsc/private/snesimpl.h>   /*I "petscsnes.h" I*/
46cab3a1bSJed Brown 
56cab3a1bSJed Brown /* This structure holds the user-provided DMDA callbacks */
66cab3a1bSJed Brown typedef struct {
76cab3a1bSJed Brown   PetscErrorCode (*residuallocal)(DMDALocalInfo*,void*,void*,void*);
8d1e9a80fSBarry Smith   PetscErrorCode (*jacobianlocal)(DMDALocalInfo*,void*,Mat,Mat,void*);
92a4ee8f2SPeter Brune   PetscErrorCode (*objectivelocal)(DMDALocalInfo*,void*,PetscReal*,void*);
105505f7afSJunchao Zhang 
115505f7afSJunchao Zhang   PetscErrorCode (*residuallocalvec)(DMDALocalInfo*,Vec,Vec,void*);
125505f7afSJunchao Zhang   PetscErrorCode (*jacobianlocalvec)(DMDALocalInfo*,Vec,Mat,Mat,void*);
135505f7afSJunchao Zhang   PetscErrorCode (*objectivelocalvec)(DMDALocalInfo*,Vec,PetscReal*,void*);
146cab3a1bSJed Brown   void       *residuallocalctx;
156cab3a1bSJed Brown   void       *jacobianlocalctx;
162a4ee8f2SPeter Brune   void       *objectivelocalctx;
17709ac0d2SJed Brown   InsertMode residuallocalimode;
18dcad5f1cSBarry Smith 
19dcad5f1cSBarry Smith   /*   For Picard iteration defined locally */
20dcad5f1cSBarry Smith   PetscErrorCode (*rhsplocal)(DMDALocalInfo*,void*,void*,void*);
21d1e9a80fSBarry Smith   PetscErrorCode (*jacobianplocal)(DMDALocalInfo*,void*,Mat,Mat,void*);
22dcad5f1cSBarry Smith   void *picardlocalctx;
23942e3340SBarry Smith } DMSNES_DA;
246cab3a1bSJed Brown 
25942e3340SBarry Smith static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm)
266cab3a1bSJed Brown {
276cab3a1bSJed Brown   PetscFunctionBegin;
289566063dSJacob Faibussowitsch   PetscCall(PetscFree(sdm->data));
296cab3a1bSJed Brown   PetscFunctionReturn(0);
306cab3a1bSJed Brown }
316cab3a1bSJed Brown 
3222c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm,DMSNES sdm)
33e3c0b8a2SPeter Brune {
34e3c0b8a2SPeter Brune   PetscFunctionBegin;
359566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(sdm,(DMSNES_DA**)&sdm->data));
36*1baa6e33SBarry Smith   if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_DA)));
37e3c0b8a2SPeter Brune   PetscFunctionReturn(0);
38e3c0b8a2SPeter Brune }
39e3c0b8a2SPeter Brune 
40942e3340SBarry Smith static PetscErrorCode DMDASNESGetContext(DM dm,DMSNES sdm,DMSNES_DA  **dmdasnes)
416cab3a1bSJed Brown {
426cab3a1bSJed Brown   PetscFunctionBegin;
430298fd71SBarry Smith   *dmdasnes = NULL;
446cab3a1bSJed Brown   if (!sdm->data) {
459566063dSJacob Faibussowitsch     PetscCall(PetscNewLog(dm,(DMSNES_DA**)&sdm->data));
4622c6f798SBarry Smith     sdm->ops->destroy   = DMSNESDestroy_DMDA;
4722c6f798SBarry Smith     sdm->ops->duplicate = DMSNESDuplicate_DMDA;
486cab3a1bSJed Brown   }
49942e3340SBarry Smith   *dmdasnes = (DMSNES_DA*)sdm->data;
506cab3a1bSJed Brown   PetscFunctionReturn(0);
516cab3a1bSJed Brown }
526cab3a1bSJed Brown 
536cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx)
546cab3a1bSJed Brown {
556cab3a1bSJed Brown   DM             dm;
56942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
576cab3a1bSJed Brown   DMDALocalInfo  info;
586cab3a1bSJed Brown   Vec            Xloc;
596cab3a1bSJed Brown   void           *x,*f;
606cab3a1bSJed Brown 
616cab3a1bSJed Brown   PetscFunctionBegin;
622961e153SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
632961e153SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
642961e153SJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
655505f7afSJunchao Zhang   PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
669566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
679566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&Xloc));
689566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
699566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
709566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
71709ac0d2SJed Brown   switch (dmdasnes->residuallocalimode) {
72709ac0d2SJed Brown   case INSERT_VALUES: {
739566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0));
746cab3a1bSJed Brown     CHKMEMQ;
755505f7afSJunchao Zhang     if (dmdasnes->residuallocalvec) {
765505f7afSJunchao Zhang       PetscCall((*dmdasnes->residuallocalvec)(&info,Xloc,F,dmdasnes->residuallocalctx));
775505f7afSJunchao Zhang     } else {
785505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm,Xloc,&x));
795505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm,F,&f));
809566063dSJacob Faibussowitsch       PetscCall((*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx));
815505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
825505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm,F,&f));
835505f7afSJunchao Zhang     }
846cab3a1bSJed Brown     CHKMEMQ;
859566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0));
86709ac0d2SJed Brown   } break;
87709ac0d2SJed Brown   case ADD_VALUES: {
88709ac0d2SJed Brown     Vec Floc;
899566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm,&Floc));
909566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
919566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0));
92709ac0d2SJed Brown     CHKMEMQ;
935505f7afSJunchao Zhang     if (dmdasnes->residuallocalvec) {
945505f7afSJunchao Zhang       PetscCall((*dmdasnes->residuallocalvec)(&info,Xloc,Floc,dmdasnes->residuallocalctx));
955505f7afSJunchao Zhang     } else {
965505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm,Xloc,&x));
975505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm,Floc,&f));
989566063dSJacob Faibussowitsch       PetscCall((*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx));
995505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
1005505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm,Floc,&f));
1015505f7afSJunchao Zhang     }
102709ac0d2SJed Brown     CHKMEMQ;
1039566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0));
1049566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
1059566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F));
1069566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F));
1079566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm,&Floc));
108709ac0d2SJed Brown   } break;
10998921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
110709ac0d2SJed Brown   }
1119566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&Xloc));
112*1baa6e33SBarry Smith   if (snes->domainerror) PetscCall(VecSetInf(F));
1136cab3a1bSJed Brown   PetscFunctionReturn(0);
1146cab3a1bSJed Brown }
1156cab3a1bSJed Brown 
1162a4ee8f2SPeter Brune static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx)
1172a4ee8f2SPeter Brune {
1182a4ee8f2SPeter Brune   DM             dm;
119942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
1202a4ee8f2SPeter Brune   DMDALocalInfo  info;
1212a4ee8f2SPeter Brune   Vec            Xloc;
1222a4ee8f2SPeter Brune   void           *x;
1232a4ee8f2SPeter Brune 
1242a4ee8f2SPeter Brune   PetscFunctionBegin;
1252a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1262a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
127dadcf809SJacob Faibussowitsch   PetscValidRealPointer(ob,3);
1285505f7afSJunchao Zhang   PetscCheck(dmdasnes->objectivelocal || dmdasnes->objectivelocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
1299566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
1309566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&Xloc));
1319566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
1329566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
1339566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
1345505f7afSJunchao Zhang   CHKMEMQ;
1355505f7afSJunchao Zhang   if (dmdasnes->objectivelocalvec) {
1365505f7afSJunchao Zhang     PetscCall((*dmdasnes->objectivelocalvec)(&info,Xloc,ob,dmdasnes->objectivelocalctx));
1375505f7afSJunchao Zhang   } else {
1389566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm,Xloc,&x));
1399566063dSJacob Faibussowitsch     PetscCall((*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx));
1409566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
1415505f7afSJunchao Zhang   }
1425505f7afSJunchao Zhang   CHKMEMQ;
1439566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&Xloc));
1442a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1452a4ee8f2SPeter Brune }
1462a4ee8f2SPeter Brune 
147cc2e6a90SBarry Smith /* Routine is called by example, hence must be labeled PETSC_EXTERN */
148cc2e6a90SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
1492961e153SJed Brown {
1502961e153SJed Brown   DM             dm;
151942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
1522961e153SJed Brown   DMDALocalInfo  info;
1532961e153SJed Brown   Vec            Xloc;
1542961e153SJed Brown   void           *x;
1552961e153SJed Brown 
1562961e153SJed Brown   PetscFunctionBegin;
1575505f7afSJunchao Zhang   PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
1589566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
1592961e153SJed Brown 
1605505f7afSJunchao Zhang   if (dmdasnes->jacobianlocal || dmdasnes->jacobianlocalctx) {
1619566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm,&Xloc));
1629566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
1639566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
1649566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(dm,&info));
1655505f7afSJunchao Zhang     CHKMEMQ;
1665505f7afSJunchao Zhang     if (dmdasnes->jacobianlocalvec) {
1675505f7afSJunchao Zhang       PetscCall((*dmdasnes->jacobianlocalvec)(&info,Xloc,A,B,dmdasnes->jacobianlocalctx));
1685505f7afSJunchao Zhang     } else {
1699566063dSJacob Faibussowitsch       PetscCall(DMDAVecGetArray(dm,Xloc,&x));
1709566063dSJacob Faibussowitsch       PetscCall((*dmdasnes->jacobianlocal)(&info,x,A,B,dmdasnes->jacobianlocalctx));
1719566063dSJacob Faibussowitsch       PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
1725505f7afSJunchao Zhang     }
1735505f7afSJunchao Zhang     CHKMEMQ;
1749566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm,&Xloc));
1752961e153SJed Brown   } else {
1762961e153SJed Brown     MatFDColoring fdcoloring;
1779566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring));
1782961e153SJed Brown     if (!fdcoloring) {
1792961e153SJed Brown       ISColoring coloring;
1802961e153SJed Brown 
1819566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm,dm->coloringtype,&coloring));
1829566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(B,coloring,&fdcoloring));
1832961e153SJed Brown       switch (dm->coloringtype) {
1842961e153SJed Brown       case IS_COLORING_GLOBAL:
1859566063dSJacob Faibussowitsch         PetscCall(MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMDA,dmdasnes));
1862961e153SJed Brown         break;
18798921bdaSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
1882961e153SJed Brown       }
1899566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix));
1909566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1919566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B,coloring,fdcoloring));
1929566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&coloring));
1939566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring));
1949566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)fdcoloring));
1952961e153SJed Brown 
1962961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
1972961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
1982961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
1992961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
200140e18c1SBarry Smith        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
2012961e153SJed Brown        */
2029566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dm));
2032961e153SJed Brown     }
2049566063dSJacob Faibussowitsch     PetscCall(MatFDColoringApply(B,fdcoloring,X,snes));
2052961e153SJed Brown   }
2062961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
20794ab13aaSBarry Smith   if (A != B) {
2089566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
2099566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
2102961e153SJed Brown   }
2112961e153SJed Brown   PetscFunctionReturn(0);
2122961e153SJed Brown }
2132961e153SJed Brown 
2146cab3a1bSJed Brown /*@C
2156cab3a1bSJed Brown    DMDASNESSetFunctionLocal - set a local residual evaluation function
2166cab3a1bSJed Brown 
2176cab3a1bSJed Brown    Logically Collective
2186cab3a1bSJed Brown 
2194165533cSJose E. Roman    Input Parameters:
2206cab3a1bSJed Brown +  dm - DM to associate callback with
221e3c51ba2SJed Brown .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
2226cab3a1bSJed Brown .  func - local residual evaluation
2236cab3a1bSJed Brown -  ctx - optional context for local residual evaluation
2246cab3a1bSJed Brown 
2250038884cSBarry Smith    Calling sequence:
2260038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x, void *f, void *ctx),
2276cab3a1bSJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
2280038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x)
2290038884cSBarry Smith .  f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f)
2306cab3a1bSJed Brown -  ctx - optional context passed above
2316cab3a1bSJed Brown 
2326cab3a1bSJed Brown    Level: beginner
2336cab3a1bSJed Brown 
234db781477SPatrick Sanan .seealso: `DMDASNESSetJacobianLocal()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2356cab3a1bSJed Brown @*/
236709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx)
2376cab3a1bSJed Brown {
238942e3340SBarry Smith   DMSNES         sdm;
239942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
2406cab3a1bSJed Brown 
2416cab3a1bSJed Brown   PetscFunctionBegin;
2426cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2439566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm,&sdm));
2449566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
2451aa26658SKarl Rupp 
246709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
2476cab3a1bSJed Brown   dmdasnes->residuallocal      = func;
2486cab3a1bSJed Brown   dmdasnes->residuallocalctx   = ctx;
2491aa26658SKarl Rupp 
2509566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes));
25122c6f798SBarry Smith   if (!sdm->ops->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
2529566063dSJacob Faibussowitsch     PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes));
2532961e153SJed Brown   }
2542961e153SJed Brown   PetscFunctionReturn(0);
2552961e153SJed Brown }
2562961e153SJed Brown 
2572961e153SJed Brown /*@C
2585505f7afSJunchao Zhang    DMDASNESSetFunctionLocalVec - set a local residual evaluation function that operates on a local vector
2595505f7afSJunchao Zhang 
2605505f7afSJunchao Zhang    Logically Collective
2615505f7afSJunchao Zhang 
2625505f7afSJunchao Zhang    Input Parameters:
2635505f7afSJunchao Zhang +  dm - DM to associate callback with
2645505f7afSJunchao Zhang .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
2655505f7afSJunchao Zhang .  func - local residual evaluation
2665505f7afSJunchao Zhang -  ctx - optional context for local residual evaluation
2675505f7afSJunchao Zhang 
2685505f7afSJunchao Zhang    Calling sequence:
2695505f7afSJunchao Zhang    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x, Vec f, void *ctx),
2705505f7afSJunchao Zhang +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
2715505f7afSJunchao Zhang .  x - state vector at which to evaluate residual
2725505f7afSJunchao Zhang .  f - residual vector
2735505f7afSJunchao Zhang -  ctx - optional context passed above
2745505f7afSJunchao Zhang 
2755505f7afSJunchao Zhang    Level: beginner
2765505f7afSJunchao Zhang 
2775505f7afSJunchao Zhang .seealso: `DMDASNESSetFunctionLocal()`, `DMDASNESSetJacobianLocalVec()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2785505f7afSJunchao Zhang @*/
2795505f7afSJunchao Zhang PetscErrorCode DMDASNESSetFunctionLocalVec(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,Vec,Vec,void*),void *ctx)
2805505f7afSJunchao Zhang {
2815505f7afSJunchao Zhang   DMSNES         sdm;
2825505f7afSJunchao Zhang   DMSNES_DA      *dmdasnes;
2835505f7afSJunchao Zhang 
2845505f7afSJunchao Zhang   PetscFunctionBegin;
2855505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2865505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm,&sdm));
2875505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
2885505f7afSJunchao Zhang 
2895505f7afSJunchao Zhang   dmdasnes->residuallocalimode = imode;
2905505f7afSJunchao Zhang   dmdasnes->residuallocalvec   = func;
2915505f7afSJunchao Zhang   dmdasnes->residuallocalctx   = ctx;
2925505f7afSJunchao Zhang 
2935505f7afSJunchao Zhang   PetscCall(DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes));
2945505f7afSJunchao Zhang   if (!sdm->ops->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
2955505f7afSJunchao Zhang     PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes));
2965505f7afSJunchao Zhang   }
2975505f7afSJunchao Zhang   PetscFunctionReturn(0);
2985505f7afSJunchao Zhang }
2995505f7afSJunchao Zhang 
3005505f7afSJunchao Zhang /*@C
301e79ed307SJed Brown    DMDASNESSetJacobianLocal - set a local Jacobian evaluation function
3022961e153SJed Brown 
3032961e153SJed Brown    Logically Collective
3042961e153SJed Brown 
3054165533cSJose E. Roman    Input Parameters:
3062961e153SJed Brown +  dm - DM to associate callback with
3076b7eb5ceSMatthew G. Knepley .  func - local Jacobian evaluation
3086b7eb5ceSMatthew G. Knepley -  ctx - optional context for local Jacobian evaluation
3092961e153SJed Brown 
3100038884cSBarry Smith    Calling sequence:
3110038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x,Mat J,Mat M,void *ctx),
3126b7eb5ceSMatthew G. Knepley +  info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at
3130038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x)
3146b7eb5ceSMatthew G. Knepley .  J - Mat object for the Jacobian
3156b7eb5ceSMatthew G. Knepley .  M - Mat object for the Jacobian preconditioner matrix
3162961e153SJed Brown -  ctx - optional context passed above
3172961e153SJed Brown 
3182961e153SJed Brown    Level: beginner
3192961e153SJed Brown 
320db781477SPatrick Sanan .seealso: `DMDASNESSetFunctionLocal()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3212961e153SJed Brown @*/
322d1e9a80fSBarry Smith PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
3232961e153SJed Brown {
324942e3340SBarry Smith   DMSNES         sdm;
325942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
3262961e153SJed Brown 
3272961e153SJed Brown   PetscFunctionBegin;
3282961e153SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3299566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm,&sdm));
3309566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
3311aa26658SKarl Rupp 
3322961e153SJed Brown   dmdasnes->jacobianlocal    = func;
3332961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
3341aa26658SKarl Rupp 
3359566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes));
3366cab3a1bSJed Brown   PetscFunctionReturn(0);
3376cab3a1bSJed Brown }
3382a4ee8f2SPeter Brune 
3392a4ee8f2SPeter Brune /*@C
3405505f7afSJunchao Zhang    DMDASNESSetJacobianLocalVec - set a local Jacobian evaluation function that operates on a local vector
3415505f7afSJunchao Zhang 
3425505f7afSJunchao Zhang    Logically Collective
3435505f7afSJunchao Zhang 
3445505f7afSJunchao Zhang    Input Parameters:
3455505f7afSJunchao Zhang +  dm - DM to associate callback with
3465505f7afSJunchao Zhang .  func - local Jacobian evaluation
3475505f7afSJunchao Zhang -  ctx - optional context for local Jacobian evaluation
3485505f7afSJunchao Zhang 
3495505f7afSJunchao Zhang    Calling sequence:
3505505f7afSJunchao Zhang    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,Mat J,Mat M,void *ctx),
3515505f7afSJunchao Zhang +  info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at
3525505f7afSJunchao Zhang .  x - state vector at which to evaluate Jacobian
3535505f7afSJunchao Zhang .  J - Mat object for the Jacobian
3545505f7afSJunchao Zhang .  M - Mat object for the Jacobian preconditioner matrix
3555505f7afSJunchao Zhang -  ctx - optional context passed above
3565505f7afSJunchao Zhang 
3575505f7afSJunchao Zhang    Level: beginner
3585505f7afSJunchao Zhang 
3595505f7afSJunchao Zhang .seealso: `DMDASNESSetJacobianLocal()`, `DMDASNESSetFunctionLocalVec()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3605505f7afSJunchao Zhang @*/
3615505f7afSJunchao Zhang PetscErrorCode DMDASNESSetJacobianLocalVec(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,Vec,Mat,Mat,void*),void *ctx)
3625505f7afSJunchao Zhang {
3635505f7afSJunchao Zhang   DMSNES         sdm;
3645505f7afSJunchao Zhang   DMSNES_DA      *dmdasnes;
3655505f7afSJunchao Zhang 
3665505f7afSJunchao Zhang   PetscFunctionBegin;
3675505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3685505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm,&sdm));
3695505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
3705505f7afSJunchao Zhang 
3715505f7afSJunchao Zhang   dmdasnes->jacobianlocalvec = func;
3725505f7afSJunchao Zhang   dmdasnes->jacobianlocalctx = ctx;
3735505f7afSJunchao Zhang 
3745505f7afSJunchao Zhang   PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes));
3755505f7afSJunchao Zhang   PetscFunctionReturn(0);
3765505f7afSJunchao Zhang }
3775505f7afSJunchao Zhang 
3785505f7afSJunchao Zhang /*@C
3792a4ee8f2SPeter Brune    DMDASNESSetObjectiveLocal - set a local residual evaluation function
3802a4ee8f2SPeter Brune 
3812a4ee8f2SPeter Brune    Logically Collective
3822a4ee8f2SPeter Brune 
3834165533cSJose E. Roman    Input Parameters:
3842a4ee8f2SPeter Brune +  dm - DM to associate callback with
3852a4ee8f2SPeter Brune .  func - local objective evaluation
3862a4ee8f2SPeter Brune -  ctx - optional context for local residual evaluation
3872a4ee8f2SPeter Brune 
3882a4ee8f2SPeter Brune    Calling sequence for func:
3892a4ee8f2SPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
3902a4ee8f2SPeter Brune .  x - dimensional pointer to state at which to evaluate residual
3912a4ee8f2SPeter Brune .  ob - eventual objective value
3922a4ee8f2SPeter Brune -  ctx - optional context passed above
3932a4ee8f2SPeter Brune 
3942a4ee8f2SPeter Brune    Level: beginner
3952a4ee8f2SPeter Brune 
396db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3972a4ee8f2SPeter Brune @*/
3982a4ee8f2SPeter Brune PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx)
3992a4ee8f2SPeter Brune {
400942e3340SBarry Smith   DMSNES         sdm;
401942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
4022a4ee8f2SPeter Brune 
4032a4ee8f2SPeter Brune   PetscFunctionBegin;
4042a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4059566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm,&sdm));
4069566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
4071aa26658SKarl Rupp 
4082a4ee8f2SPeter Brune   dmdasnes->objectivelocal    = func;
4092a4ee8f2SPeter Brune   dmdasnes->objectivelocalctx = ctx;
4101aa26658SKarl Rupp 
4119566063dSJacob Faibussowitsch   PetscCall(DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes));
4122a4ee8f2SPeter Brune   PetscFunctionReturn(0);
4132a4ee8f2SPeter Brune }
414dcad5f1cSBarry Smith 
4155505f7afSJunchao Zhang /*@C
4165505f7afSJunchao Zhang    DMDASNESSetObjectiveLocal - set a local residual evaluation function that operates on a local vector
4175505f7afSJunchao Zhang 
4185505f7afSJunchao Zhang    Logically Collective
4195505f7afSJunchao Zhang 
4205505f7afSJunchao Zhang    Input Parameters:
4215505f7afSJunchao Zhang +  dm - DM to associate callback with
4225505f7afSJunchao Zhang .  func - local objective evaluation
4235505f7afSJunchao Zhang -  ctx - optional context for local residual evaluation
4245505f7afSJunchao Zhang 
4255505f7afSJunchao Zhang    Calling sequence
4265505f7afSJunchao Zhang    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,PetscReal *ob,void *ctx);
4275505f7afSJunchao Zhang +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
4285505f7afSJunchao Zhang .  x - state vector at which to evaluate residual
4295505f7afSJunchao Zhang .  ob - eventual objective value
4305505f7afSJunchao Zhang -  ctx - optional context passed above
4315505f7afSJunchao Zhang 
4325505f7afSJunchao Zhang    Level: beginner
4335505f7afSJunchao Zhang 
4345505f7afSJunchao Zhang .seealso: `DMDASNESSetObjectiveLocal()`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocalVec()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
4355505f7afSJunchao Zhang @*/
4365505f7afSJunchao Zhang PetscErrorCode DMDASNESSetObjectiveLocalVec(DM dm,DMDASNESObjectiveVec func,void *ctx)
4375505f7afSJunchao Zhang {
4385505f7afSJunchao Zhang   DMSNES         sdm;
4395505f7afSJunchao Zhang   DMSNES_DA      *dmdasnes;
4405505f7afSJunchao Zhang 
4415505f7afSJunchao Zhang   PetscFunctionBegin;
4425505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4435505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm,&sdm));
4445505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
4455505f7afSJunchao Zhang 
4465505f7afSJunchao Zhang   dmdasnes->objectivelocalvec = func;
4475505f7afSJunchao Zhang   dmdasnes->objectivelocalctx = ctx;
4485505f7afSJunchao Zhang 
4495505f7afSJunchao Zhang   PetscCall(DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes));
4505505f7afSJunchao Zhang   PetscFunctionReturn(0);
4515505f7afSJunchao Zhang }
4525505f7afSJunchao Zhang 
453dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx)
454dcad5f1cSBarry Smith {
455dcad5f1cSBarry Smith   DM             dm;
456942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
457dcad5f1cSBarry Smith   DMDALocalInfo  info;
458dcad5f1cSBarry Smith   Vec            Xloc;
459dcad5f1cSBarry Smith   void           *x,*f;
460dcad5f1cSBarry Smith 
461dcad5f1cSBarry Smith   PetscFunctionBegin;
462dcad5f1cSBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
463dcad5f1cSBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
464dcad5f1cSBarry Smith   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
46528b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->rhsplocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
4669566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
4679566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&Xloc));
4689566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
4699566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
4709566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
4719566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm,Xloc,&x));
472dcad5f1cSBarry Smith   switch (dmdasnes->residuallocalimode) {
473dcad5f1cSBarry Smith   case INSERT_VALUES: {
4749566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm,F,&f));
475dcad5f1cSBarry Smith     CHKMEMQ;
4769566063dSJacob Faibussowitsch     PetscCall((*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx));
477dcad5f1cSBarry Smith     CHKMEMQ;
4789566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm,F,&f));
479dcad5f1cSBarry Smith   } break;
480dcad5f1cSBarry Smith   case ADD_VALUES: {
481dcad5f1cSBarry Smith     Vec Floc;
4829566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm,&Floc));
4839566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
4849566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm,Floc,&f));
485dcad5f1cSBarry Smith     CHKMEMQ;
4869566063dSJacob Faibussowitsch     PetscCall((*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx));
487dcad5f1cSBarry Smith     CHKMEMQ;
4889566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm,Floc,&f));
4899566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
4909566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F));
4919566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F));
4929566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm,&Floc));
493dcad5f1cSBarry Smith   } break;
49498921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
495dcad5f1cSBarry Smith   }
4969566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
4979566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&Xloc));
498dcad5f1cSBarry Smith   PetscFunctionReturn(0);
499dcad5f1cSBarry Smith }
500dcad5f1cSBarry Smith 
501d1e9a80fSBarry Smith static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
502dcad5f1cSBarry Smith {
503dcad5f1cSBarry Smith   DM             dm;
504942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
505dcad5f1cSBarry Smith   DMDALocalInfo  info;
506dcad5f1cSBarry Smith   Vec            Xloc;
507dcad5f1cSBarry Smith   void           *x;
508dcad5f1cSBarry Smith 
509dcad5f1cSBarry Smith   PetscFunctionBegin;
51028b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->jacobianplocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
5119566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
512dcad5f1cSBarry Smith 
5139566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&Xloc));
5149566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
5159566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
5169566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
5179566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm,Xloc,&x));
518dcad5f1cSBarry Smith   CHKMEMQ;
5199566063dSJacob Faibussowitsch   PetscCall((*dmdasnes->jacobianplocal)(&info,x,A,B,dmdasnes->picardlocalctx));
520dcad5f1cSBarry Smith   CHKMEMQ;
5219566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
5229566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&Xloc));
52394ab13aaSBarry Smith   if (A != B) {
5249566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
5259566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
526dcad5f1cSBarry Smith   }
527dcad5f1cSBarry Smith   PetscFunctionReturn(0);
528dcad5f1cSBarry Smith }
529dcad5f1cSBarry Smith 
530dcad5f1cSBarry Smith /*@C
531dcad5f1cSBarry Smith    DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration
532dcad5f1cSBarry Smith 
533dcad5f1cSBarry Smith    Logically Collective
534dcad5f1cSBarry Smith 
5354165533cSJose E. Roman    Input Parameters:
536dcad5f1cSBarry Smith +  dm - DM to associate callback with
537dcad5f1cSBarry Smith .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
538dcad5f1cSBarry Smith .  func - local residual evaluation
539dcad5f1cSBarry Smith -  ctx - optional context for local residual evaluation
540dcad5f1cSBarry Smith 
541dcad5f1cSBarry Smith    Calling sequence for func:
542dcad5f1cSBarry Smith +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
543dcad5f1cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual
544dcad5f1cSBarry Smith .  f - dimensional pointer to residual, write the residual here
545dcad5f1cSBarry Smith -  ctx - optional context passed above
546dcad5f1cSBarry Smith 
54795452b02SPatrick Sanan    Notes:
54895452b02SPatrick Sanan     The user must use
5499566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user));
550dcad5f1cSBarry Smith     in their code before calling this routine.
551dcad5f1cSBarry Smith 
552dcad5f1cSBarry Smith    Level: beginner
553dcad5f1cSBarry Smith 
554db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
555dcad5f1cSBarry Smith @*/
556dcad5f1cSBarry Smith PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),
557d1e9a80fSBarry Smith                                       PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
558dcad5f1cSBarry Smith {
559942e3340SBarry Smith   DMSNES         sdm;
560942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
561dcad5f1cSBarry Smith 
562dcad5f1cSBarry Smith   PetscFunctionBegin;
563dcad5f1cSBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5649566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm,&sdm));
5659566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
5661aa26658SKarl Rupp 
567dcad5f1cSBarry Smith   dmdasnes->residuallocalimode = imode;
568dcad5f1cSBarry Smith   dmdasnes->rhsplocal          = func;
569dcad5f1cSBarry Smith   dmdasnes->jacobianplocal     = jac;
570dcad5f1cSBarry Smith   dmdasnes->picardlocalctx     = ctx;
5711aa26658SKarl Rupp 
5729566063dSJacob Faibussowitsch   PetscCall(DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes));
5739566063dSJacob Faibussowitsch   PetscCall(DMSNESSetMFFunction(dm,SNESComputeFunction_DMDA,dmdasnes));
574dcad5f1cSBarry Smith   PetscFunctionReturn(0);
575dcad5f1cSBarry Smith }
576