xref: /petsc/src/snes/utils/dmdasnes.c (revision 5505f7af238fa9b78e052818df902f194199a144)
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*);
10*5505f7afSJunchao Zhang 
11*5505f7afSJunchao Zhang   PetscErrorCode (*residuallocalvec)(DMDALocalInfo*,Vec,Vec,void*);
12*5505f7afSJunchao Zhang   PetscErrorCode (*jacobianlocalvec)(DMDALocalInfo*,Vec,Mat,Mat,void*);
13*5505f7afSJunchao 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));
36dcad5f1cSBarry Smith   if (oldsdm->data) {
379566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_DA)));
38dcad5f1cSBarry Smith   }
39e3c0b8a2SPeter Brune   PetscFunctionReturn(0);
40e3c0b8a2SPeter Brune }
41e3c0b8a2SPeter Brune 
42942e3340SBarry Smith static PetscErrorCode DMDASNESGetContext(DM dm,DMSNES sdm,DMSNES_DA  **dmdasnes)
436cab3a1bSJed Brown {
446cab3a1bSJed Brown   PetscFunctionBegin;
450298fd71SBarry Smith   *dmdasnes = NULL;
466cab3a1bSJed Brown   if (!sdm->data) {
479566063dSJacob Faibussowitsch     PetscCall(PetscNewLog(dm,(DMSNES_DA**)&sdm->data));
4822c6f798SBarry Smith     sdm->ops->destroy   = DMSNESDestroy_DMDA;
4922c6f798SBarry Smith     sdm->ops->duplicate = DMSNESDuplicate_DMDA;
506cab3a1bSJed Brown   }
51942e3340SBarry Smith   *dmdasnes = (DMSNES_DA*)sdm->data;
526cab3a1bSJed Brown   PetscFunctionReturn(0);
536cab3a1bSJed Brown }
546cab3a1bSJed Brown 
556cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx)
566cab3a1bSJed Brown {
576cab3a1bSJed Brown   DM             dm;
58942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
596cab3a1bSJed Brown   DMDALocalInfo  info;
606cab3a1bSJed Brown   Vec            Xloc;
616cab3a1bSJed Brown   void           *x,*f;
626cab3a1bSJed Brown 
636cab3a1bSJed Brown   PetscFunctionBegin;
642961e153SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
652961e153SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
662961e153SJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
67*5505f7afSJunchao Zhang   PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
689566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
699566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&Xloc));
709566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
719566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
729566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
73709ac0d2SJed Brown   switch (dmdasnes->residuallocalimode) {
74709ac0d2SJed Brown   case INSERT_VALUES: {
759566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0));
766cab3a1bSJed Brown     CHKMEMQ;
77*5505f7afSJunchao Zhang     if (dmdasnes->residuallocalvec) {
78*5505f7afSJunchao Zhang       PetscCall((*dmdasnes->residuallocalvec)(&info,Xloc,F,dmdasnes->residuallocalctx));
79*5505f7afSJunchao Zhang     } else {
80*5505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm,Xloc,&x));
81*5505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm,F,&f));
829566063dSJacob Faibussowitsch       PetscCall((*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx));
83*5505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
84*5505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm,F,&f));
85*5505f7afSJunchao Zhang     }
866cab3a1bSJed Brown     CHKMEMQ;
879566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0));
88709ac0d2SJed Brown   } break;
89709ac0d2SJed Brown   case ADD_VALUES: {
90709ac0d2SJed Brown     Vec Floc;
919566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm,&Floc));
929566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
939566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0));
94709ac0d2SJed Brown     CHKMEMQ;
95*5505f7afSJunchao Zhang     if (dmdasnes->residuallocalvec) {
96*5505f7afSJunchao Zhang       PetscCall((*dmdasnes->residuallocalvec)(&info,Xloc,Floc,dmdasnes->residuallocalctx));
97*5505f7afSJunchao Zhang     } else {
98*5505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm,Xloc,&x));
99*5505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm,Floc,&f));
1009566063dSJacob Faibussowitsch       PetscCall((*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx));
101*5505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
102*5505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm,Floc,&f));
103*5505f7afSJunchao Zhang     }
104709ac0d2SJed Brown     CHKMEMQ;
1059566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0));
1069566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
1079566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F));
1089566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F));
1099566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm,&Floc));
110709ac0d2SJed Brown   } break;
11198921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
112709ac0d2SJed Brown   }
1139566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&Xloc));
114422a814eSBarry Smith   if (snes->domainerror) {
1159566063dSJacob Faibussowitsch     PetscCall(VecSetInf(F));
116422a814eSBarry Smith   }
1176cab3a1bSJed Brown   PetscFunctionReturn(0);
1186cab3a1bSJed Brown }
1196cab3a1bSJed Brown 
1202a4ee8f2SPeter Brune static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx)
1212a4ee8f2SPeter Brune {
1222a4ee8f2SPeter Brune   DM             dm;
123942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
1242a4ee8f2SPeter Brune   DMDALocalInfo  info;
1252a4ee8f2SPeter Brune   Vec            Xloc;
1262a4ee8f2SPeter Brune   void           *x;
1272a4ee8f2SPeter Brune 
1282a4ee8f2SPeter Brune   PetscFunctionBegin;
1292a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1302a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
131dadcf809SJacob Faibussowitsch   PetscValidRealPointer(ob,3);
132*5505f7afSJunchao Zhang   PetscCheck(dmdasnes->objectivelocal || dmdasnes->objectivelocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
1339566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
1349566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&Xloc));
1359566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
1369566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
1379566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
138*5505f7afSJunchao Zhang   CHKMEMQ;
139*5505f7afSJunchao Zhang   if (dmdasnes->objectivelocalvec) {
140*5505f7afSJunchao Zhang     PetscCall((*dmdasnes->objectivelocalvec)(&info,Xloc,ob,dmdasnes->objectivelocalctx));
141*5505f7afSJunchao Zhang   } else {
1429566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm,Xloc,&x));
1439566063dSJacob Faibussowitsch     PetscCall((*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx));
1449566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
145*5505f7afSJunchao Zhang   }
146*5505f7afSJunchao Zhang   CHKMEMQ;
1479566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&Xloc));
1482a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1492a4ee8f2SPeter Brune }
1502a4ee8f2SPeter Brune 
151cc2e6a90SBarry Smith /* Routine is called by example, hence must be labeled PETSC_EXTERN */
152cc2e6a90SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
1532961e153SJed Brown {
1542961e153SJed Brown   DM             dm;
155942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
1562961e153SJed Brown   DMDALocalInfo  info;
1572961e153SJed Brown   Vec            Xloc;
1582961e153SJed Brown   void           *x;
1592961e153SJed Brown 
1602961e153SJed Brown   PetscFunctionBegin;
161*5505f7afSJunchao Zhang   PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
1629566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
1632961e153SJed Brown 
164*5505f7afSJunchao Zhang   if (dmdasnes->jacobianlocal || dmdasnes->jacobianlocalctx) {
1659566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm,&Xloc));
1669566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
1679566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
1689566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(dm,&info));
169*5505f7afSJunchao Zhang     CHKMEMQ;
170*5505f7afSJunchao Zhang     if (dmdasnes->jacobianlocalvec) {
171*5505f7afSJunchao Zhang       PetscCall((*dmdasnes->jacobianlocalvec)(&info,Xloc,A,B,dmdasnes->jacobianlocalctx));
172*5505f7afSJunchao Zhang     } else {
1739566063dSJacob Faibussowitsch       PetscCall(DMDAVecGetArray(dm,Xloc,&x));
1749566063dSJacob Faibussowitsch       PetscCall((*dmdasnes->jacobianlocal)(&info,x,A,B,dmdasnes->jacobianlocalctx));
1759566063dSJacob Faibussowitsch       PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
176*5505f7afSJunchao Zhang     }
177*5505f7afSJunchao Zhang     CHKMEMQ;
1789566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm,&Xloc));
1792961e153SJed Brown   } else {
1802961e153SJed Brown     MatFDColoring fdcoloring;
1819566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring));
1822961e153SJed Brown     if (!fdcoloring) {
1832961e153SJed Brown       ISColoring coloring;
1842961e153SJed Brown 
1859566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm,dm->coloringtype,&coloring));
1869566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(B,coloring,&fdcoloring));
1872961e153SJed Brown       switch (dm->coloringtype) {
1882961e153SJed Brown       case IS_COLORING_GLOBAL:
1899566063dSJacob Faibussowitsch         PetscCall(MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMDA,dmdasnes));
1902961e153SJed Brown         break;
19198921bdaSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
1922961e153SJed Brown       }
1939566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix));
1949566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1959566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B,coloring,fdcoloring));
1969566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&coloring));
1979566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring));
1989566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)fdcoloring));
1992961e153SJed Brown 
2002961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
2012961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
2022961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
2032961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
204140e18c1SBarry Smith        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
2052961e153SJed Brown        */
2069566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dm));
2072961e153SJed Brown     }
2089566063dSJacob Faibussowitsch     PetscCall(MatFDColoringApply(B,fdcoloring,X,snes));
2092961e153SJed Brown   }
2102961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
21194ab13aaSBarry Smith   if (A != B) {
2129566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
2139566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
2142961e153SJed Brown   }
2152961e153SJed Brown   PetscFunctionReturn(0);
2162961e153SJed Brown }
2172961e153SJed Brown 
2186cab3a1bSJed Brown /*@C
2196cab3a1bSJed Brown    DMDASNESSetFunctionLocal - set a local residual evaluation function
2206cab3a1bSJed Brown 
2216cab3a1bSJed Brown    Logically Collective
2226cab3a1bSJed Brown 
2234165533cSJose E. Roman    Input Parameters:
2246cab3a1bSJed Brown +  dm - DM to associate callback with
225e3c51ba2SJed Brown .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
2266cab3a1bSJed Brown .  func - local residual evaluation
2276cab3a1bSJed Brown -  ctx - optional context for local residual evaluation
2286cab3a1bSJed Brown 
2290038884cSBarry Smith    Calling sequence:
2300038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x, void *f, void *ctx),
2316cab3a1bSJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
2320038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x)
2330038884cSBarry Smith .  f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f)
2346cab3a1bSJed Brown -  ctx - optional context passed above
2356cab3a1bSJed Brown 
2366cab3a1bSJed Brown    Level: beginner
2376cab3a1bSJed Brown 
238db781477SPatrick Sanan .seealso: `DMDASNESSetJacobianLocal()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2396cab3a1bSJed Brown @*/
240709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx)
2416cab3a1bSJed Brown {
242942e3340SBarry Smith   DMSNES         sdm;
243942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
2446cab3a1bSJed Brown 
2456cab3a1bSJed Brown   PetscFunctionBegin;
2466cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2479566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm,&sdm));
2489566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
2491aa26658SKarl Rupp 
250709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
2516cab3a1bSJed Brown   dmdasnes->residuallocal      = func;
2526cab3a1bSJed Brown   dmdasnes->residuallocalctx   = ctx;
2531aa26658SKarl Rupp 
2549566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes));
25522c6f798SBarry Smith   if (!sdm->ops->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
2569566063dSJacob Faibussowitsch     PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes));
2572961e153SJed Brown   }
2582961e153SJed Brown   PetscFunctionReturn(0);
2592961e153SJed Brown }
2602961e153SJed Brown 
2612961e153SJed Brown /*@C
262*5505f7afSJunchao Zhang    DMDASNESSetFunctionLocalVec - set a local residual evaluation function that operates on a local vector
263*5505f7afSJunchao Zhang 
264*5505f7afSJunchao Zhang    Logically Collective
265*5505f7afSJunchao Zhang 
266*5505f7afSJunchao Zhang    Input Parameters:
267*5505f7afSJunchao Zhang +  dm - DM to associate callback with
268*5505f7afSJunchao Zhang .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
269*5505f7afSJunchao Zhang .  func - local residual evaluation
270*5505f7afSJunchao Zhang -  ctx - optional context for local residual evaluation
271*5505f7afSJunchao Zhang 
272*5505f7afSJunchao Zhang    Calling sequence:
273*5505f7afSJunchao Zhang    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x, Vec f, void *ctx),
274*5505f7afSJunchao Zhang +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
275*5505f7afSJunchao Zhang .  x - state vector at which to evaluate residual
276*5505f7afSJunchao Zhang .  f - residual vector
277*5505f7afSJunchao Zhang -  ctx - optional context passed above
278*5505f7afSJunchao Zhang 
279*5505f7afSJunchao Zhang    Level: beginner
280*5505f7afSJunchao Zhang 
281*5505f7afSJunchao Zhang .seealso: `DMDASNESSetFunctionLocal()`, `DMDASNESSetJacobianLocalVec()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
282*5505f7afSJunchao Zhang @*/
283*5505f7afSJunchao Zhang PetscErrorCode DMDASNESSetFunctionLocalVec(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,Vec,Vec,void*),void *ctx)
284*5505f7afSJunchao Zhang {
285*5505f7afSJunchao Zhang   DMSNES         sdm;
286*5505f7afSJunchao Zhang   DMSNES_DA      *dmdasnes;
287*5505f7afSJunchao Zhang 
288*5505f7afSJunchao Zhang   PetscFunctionBegin;
289*5505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
290*5505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm,&sdm));
291*5505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
292*5505f7afSJunchao Zhang 
293*5505f7afSJunchao Zhang   dmdasnes->residuallocalimode = imode;
294*5505f7afSJunchao Zhang   dmdasnes->residuallocalvec   = func;
295*5505f7afSJunchao Zhang   dmdasnes->residuallocalctx   = ctx;
296*5505f7afSJunchao Zhang 
297*5505f7afSJunchao Zhang   PetscCall(DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes));
298*5505f7afSJunchao Zhang   if (!sdm->ops->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
299*5505f7afSJunchao Zhang     PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes));
300*5505f7afSJunchao Zhang   }
301*5505f7afSJunchao Zhang   PetscFunctionReturn(0);
302*5505f7afSJunchao Zhang }
303*5505f7afSJunchao Zhang 
304*5505f7afSJunchao Zhang /*@C
305e79ed307SJed Brown    DMDASNESSetJacobianLocal - set a local Jacobian evaluation function
3062961e153SJed Brown 
3072961e153SJed Brown    Logically Collective
3082961e153SJed Brown 
3094165533cSJose E. Roman    Input Parameters:
3102961e153SJed Brown +  dm - DM to associate callback with
3116b7eb5ceSMatthew G. Knepley .  func - local Jacobian evaluation
3126b7eb5ceSMatthew G. Knepley -  ctx - optional context for local Jacobian evaluation
3132961e153SJed Brown 
3140038884cSBarry Smith    Calling sequence:
3150038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x,Mat J,Mat M,void *ctx),
3166b7eb5ceSMatthew G. Knepley +  info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at
3170038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x)
3186b7eb5ceSMatthew G. Knepley .  J - Mat object for the Jacobian
3196b7eb5ceSMatthew G. Knepley .  M - Mat object for the Jacobian preconditioner matrix
3202961e153SJed Brown -  ctx - optional context passed above
3212961e153SJed Brown 
3222961e153SJed Brown    Level: beginner
3232961e153SJed Brown 
324db781477SPatrick Sanan .seealso: `DMDASNESSetFunctionLocal()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3252961e153SJed Brown @*/
326d1e9a80fSBarry Smith PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
3272961e153SJed Brown {
328942e3340SBarry Smith   DMSNES         sdm;
329942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
3302961e153SJed Brown 
3312961e153SJed Brown   PetscFunctionBegin;
3322961e153SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3339566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm,&sdm));
3349566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
3351aa26658SKarl Rupp 
3362961e153SJed Brown   dmdasnes->jacobianlocal    = func;
3372961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
3381aa26658SKarl Rupp 
3399566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes));
3406cab3a1bSJed Brown   PetscFunctionReturn(0);
3416cab3a1bSJed Brown }
3422a4ee8f2SPeter Brune 
3432a4ee8f2SPeter Brune /*@C
344*5505f7afSJunchao Zhang    DMDASNESSetJacobianLocalVec - set a local Jacobian evaluation function that operates on a local vector
345*5505f7afSJunchao Zhang 
346*5505f7afSJunchao Zhang    Logically Collective
347*5505f7afSJunchao Zhang 
348*5505f7afSJunchao Zhang    Input Parameters:
349*5505f7afSJunchao Zhang +  dm - DM to associate callback with
350*5505f7afSJunchao Zhang .  func - local Jacobian evaluation
351*5505f7afSJunchao Zhang -  ctx - optional context for local Jacobian evaluation
352*5505f7afSJunchao Zhang 
353*5505f7afSJunchao Zhang    Calling sequence:
354*5505f7afSJunchao Zhang    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,Mat J,Mat M,void *ctx),
355*5505f7afSJunchao Zhang +  info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at
356*5505f7afSJunchao Zhang .  x - state vector at which to evaluate Jacobian
357*5505f7afSJunchao Zhang .  J - Mat object for the Jacobian
358*5505f7afSJunchao Zhang .  M - Mat object for the Jacobian preconditioner matrix
359*5505f7afSJunchao Zhang -  ctx - optional context passed above
360*5505f7afSJunchao Zhang 
361*5505f7afSJunchao Zhang    Level: beginner
362*5505f7afSJunchao Zhang 
363*5505f7afSJunchao Zhang .seealso: `DMDASNESSetJacobianLocal()`, `DMDASNESSetFunctionLocalVec()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
364*5505f7afSJunchao Zhang @*/
365*5505f7afSJunchao Zhang PetscErrorCode DMDASNESSetJacobianLocalVec(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,Vec,Mat,Mat,void*),void *ctx)
366*5505f7afSJunchao Zhang {
367*5505f7afSJunchao Zhang   DMSNES         sdm;
368*5505f7afSJunchao Zhang   DMSNES_DA      *dmdasnes;
369*5505f7afSJunchao Zhang 
370*5505f7afSJunchao Zhang   PetscFunctionBegin;
371*5505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
372*5505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm,&sdm));
373*5505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
374*5505f7afSJunchao Zhang 
375*5505f7afSJunchao Zhang   dmdasnes->jacobianlocalvec = func;
376*5505f7afSJunchao Zhang   dmdasnes->jacobianlocalctx = ctx;
377*5505f7afSJunchao Zhang 
378*5505f7afSJunchao Zhang   PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes));
379*5505f7afSJunchao Zhang   PetscFunctionReturn(0);
380*5505f7afSJunchao Zhang }
381*5505f7afSJunchao Zhang 
382*5505f7afSJunchao Zhang /*@C
3832a4ee8f2SPeter Brune    DMDASNESSetObjectiveLocal - set a local residual evaluation function
3842a4ee8f2SPeter Brune 
3852a4ee8f2SPeter Brune    Logically Collective
3862a4ee8f2SPeter Brune 
3874165533cSJose E. Roman    Input Parameters:
3882a4ee8f2SPeter Brune +  dm - DM to associate callback with
3892a4ee8f2SPeter Brune .  func - local objective evaluation
3902a4ee8f2SPeter Brune -  ctx - optional context for local residual evaluation
3912a4ee8f2SPeter Brune 
3922a4ee8f2SPeter Brune    Calling sequence for func:
3932a4ee8f2SPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
3942a4ee8f2SPeter Brune .  x - dimensional pointer to state at which to evaluate residual
3952a4ee8f2SPeter Brune .  ob - eventual objective value
3962a4ee8f2SPeter Brune -  ctx - optional context passed above
3972a4ee8f2SPeter Brune 
3982a4ee8f2SPeter Brune    Level: beginner
3992a4ee8f2SPeter Brune 
400db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
4012a4ee8f2SPeter Brune @*/
4022a4ee8f2SPeter Brune PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx)
4032a4ee8f2SPeter Brune {
404942e3340SBarry Smith   DMSNES         sdm;
405942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
4062a4ee8f2SPeter Brune 
4072a4ee8f2SPeter Brune   PetscFunctionBegin;
4082a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4099566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm,&sdm));
4109566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
4111aa26658SKarl Rupp 
4122a4ee8f2SPeter Brune   dmdasnes->objectivelocal    = func;
4132a4ee8f2SPeter Brune   dmdasnes->objectivelocalctx = ctx;
4141aa26658SKarl Rupp 
4159566063dSJacob Faibussowitsch   PetscCall(DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes));
4162a4ee8f2SPeter Brune   PetscFunctionReturn(0);
4172a4ee8f2SPeter Brune }
418dcad5f1cSBarry Smith 
419*5505f7afSJunchao Zhang /*@C
420*5505f7afSJunchao Zhang    DMDASNESSetObjectiveLocal - set a local residual evaluation function that operates on a local vector
421*5505f7afSJunchao Zhang 
422*5505f7afSJunchao Zhang    Logically Collective
423*5505f7afSJunchao Zhang 
424*5505f7afSJunchao Zhang    Input Parameters:
425*5505f7afSJunchao Zhang +  dm - DM to associate callback with
426*5505f7afSJunchao Zhang .  func - local objective evaluation
427*5505f7afSJunchao Zhang -  ctx - optional context for local residual evaluation
428*5505f7afSJunchao Zhang 
429*5505f7afSJunchao Zhang    Calling sequence
430*5505f7afSJunchao Zhang    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,PetscReal *ob,void *ctx);
431*5505f7afSJunchao Zhang +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
432*5505f7afSJunchao Zhang .  x - state vector at which to evaluate residual
433*5505f7afSJunchao Zhang .  ob - eventual objective value
434*5505f7afSJunchao Zhang -  ctx - optional context passed above
435*5505f7afSJunchao Zhang 
436*5505f7afSJunchao Zhang    Level: beginner
437*5505f7afSJunchao Zhang 
438*5505f7afSJunchao Zhang .seealso: `DMDASNESSetObjectiveLocal()`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocalVec()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
439*5505f7afSJunchao Zhang @*/
440*5505f7afSJunchao Zhang PetscErrorCode DMDASNESSetObjectiveLocalVec(DM dm,DMDASNESObjectiveVec func,void *ctx)
441*5505f7afSJunchao Zhang {
442*5505f7afSJunchao Zhang   DMSNES         sdm;
443*5505f7afSJunchao Zhang   DMSNES_DA      *dmdasnes;
444*5505f7afSJunchao Zhang 
445*5505f7afSJunchao Zhang   PetscFunctionBegin;
446*5505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
447*5505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm,&sdm));
448*5505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
449*5505f7afSJunchao Zhang 
450*5505f7afSJunchao Zhang   dmdasnes->objectivelocalvec = func;
451*5505f7afSJunchao Zhang   dmdasnes->objectivelocalctx = ctx;
452*5505f7afSJunchao Zhang 
453*5505f7afSJunchao Zhang   PetscCall(DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes));
454*5505f7afSJunchao Zhang   PetscFunctionReturn(0);
455*5505f7afSJunchao Zhang }
456*5505f7afSJunchao Zhang 
457dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx)
458dcad5f1cSBarry Smith {
459dcad5f1cSBarry Smith   DM             dm;
460942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
461dcad5f1cSBarry Smith   DMDALocalInfo  info;
462dcad5f1cSBarry Smith   Vec            Xloc;
463dcad5f1cSBarry Smith   void           *x,*f;
464dcad5f1cSBarry Smith 
465dcad5f1cSBarry Smith   PetscFunctionBegin;
466dcad5f1cSBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
467dcad5f1cSBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
468dcad5f1cSBarry Smith   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
46928b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->rhsplocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
4709566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
4719566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&Xloc));
4729566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
4739566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
4749566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
4759566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm,Xloc,&x));
476dcad5f1cSBarry Smith   switch (dmdasnes->residuallocalimode) {
477dcad5f1cSBarry Smith   case INSERT_VALUES: {
4789566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm,F,&f));
479dcad5f1cSBarry Smith     CHKMEMQ;
4809566063dSJacob Faibussowitsch     PetscCall((*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx));
481dcad5f1cSBarry Smith     CHKMEMQ;
4829566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm,F,&f));
483dcad5f1cSBarry Smith   } break;
484dcad5f1cSBarry Smith   case ADD_VALUES: {
485dcad5f1cSBarry Smith     Vec Floc;
4869566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm,&Floc));
4879566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
4889566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm,Floc,&f));
489dcad5f1cSBarry Smith     CHKMEMQ;
4909566063dSJacob Faibussowitsch     PetscCall((*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx));
491dcad5f1cSBarry Smith     CHKMEMQ;
4929566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm,Floc,&f));
4939566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
4949566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F));
4959566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F));
4969566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm,&Floc));
497dcad5f1cSBarry Smith   } break;
49898921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
499dcad5f1cSBarry Smith   }
5009566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
5019566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&Xloc));
502dcad5f1cSBarry Smith   PetscFunctionReturn(0);
503dcad5f1cSBarry Smith }
504dcad5f1cSBarry Smith 
505d1e9a80fSBarry Smith static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
506dcad5f1cSBarry Smith {
507dcad5f1cSBarry Smith   DM             dm;
508942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
509dcad5f1cSBarry Smith   DMDALocalInfo  info;
510dcad5f1cSBarry Smith   Vec            Xloc;
511dcad5f1cSBarry Smith   void           *x;
512dcad5f1cSBarry Smith 
513dcad5f1cSBarry Smith   PetscFunctionBegin;
51428b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->jacobianplocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
5159566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
516dcad5f1cSBarry Smith 
5179566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&Xloc));
5189566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
5199566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
5209566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
5219566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm,Xloc,&x));
522dcad5f1cSBarry Smith   CHKMEMQ;
5239566063dSJacob Faibussowitsch   PetscCall((*dmdasnes->jacobianplocal)(&info,x,A,B,dmdasnes->picardlocalctx));
524dcad5f1cSBarry Smith   CHKMEMQ;
5259566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm,Xloc,&x));
5269566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&Xloc));
52794ab13aaSBarry Smith   if (A != B) {
5289566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
5299566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
530dcad5f1cSBarry Smith   }
531dcad5f1cSBarry Smith   PetscFunctionReturn(0);
532dcad5f1cSBarry Smith }
533dcad5f1cSBarry Smith 
534dcad5f1cSBarry Smith /*@C
535dcad5f1cSBarry Smith    DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration
536dcad5f1cSBarry Smith 
537dcad5f1cSBarry Smith    Logically Collective
538dcad5f1cSBarry Smith 
5394165533cSJose E. Roman    Input Parameters:
540dcad5f1cSBarry Smith +  dm - DM to associate callback with
541dcad5f1cSBarry Smith .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
542dcad5f1cSBarry Smith .  func - local residual evaluation
543dcad5f1cSBarry Smith -  ctx - optional context for local residual evaluation
544dcad5f1cSBarry Smith 
545dcad5f1cSBarry Smith    Calling sequence for func:
546dcad5f1cSBarry Smith +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
547dcad5f1cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual
548dcad5f1cSBarry Smith .  f - dimensional pointer to residual, write the residual here
549dcad5f1cSBarry Smith -  ctx - optional context passed above
550dcad5f1cSBarry Smith 
55195452b02SPatrick Sanan    Notes:
55295452b02SPatrick Sanan     The user must use
5539566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user));
554dcad5f1cSBarry Smith     in their code before calling this routine.
555dcad5f1cSBarry Smith 
556dcad5f1cSBarry Smith    Level: beginner
557dcad5f1cSBarry Smith 
558db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
559dcad5f1cSBarry Smith @*/
560dcad5f1cSBarry Smith PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),
561d1e9a80fSBarry Smith                                       PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
562dcad5f1cSBarry Smith {
563942e3340SBarry Smith   DMSNES         sdm;
564942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
565dcad5f1cSBarry Smith 
566dcad5f1cSBarry Smith   PetscFunctionBegin;
567dcad5f1cSBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5689566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm,&sdm));
5699566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes));
5701aa26658SKarl Rupp 
571dcad5f1cSBarry Smith   dmdasnes->residuallocalimode = imode;
572dcad5f1cSBarry Smith   dmdasnes->rhsplocal          = func;
573dcad5f1cSBarry Smith   dmdasnes->jacobianplocal     = jac;
574dcad5f1cSBarry Smith   dmdasnes->picardlocalctx     = ctx;
5751aa26658SKarl Rupp 
5769566063dSJacob Faibussowitsch   PetscCall(DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes));
5779566063dSJacob Faibussowitsch   PetscCall(DMSNESSetMFFunction(dm,SNESComputeFunction_DMDA,dmdasnes));
578dcad5f1cSBarry Smith   PetscFunctionReturn(0);
579dcad5f1cSBarry Smith }
580