xref: /petsc/src/snes/utils/dmdasnes.c (revision 0bf52853de8b8fc6e575fbdf9d44e7daab01f4c7)
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 {
7*0bf52853SStefano Zampini   /* array versions for vector data */
8*0bf52853SStefano Zampini   DMDASNESFunction  residuallocal;
9*0bf52853SStefano Zampini   DMDASNESJacobian  jacobianlocal;
10*0bf52853SStefano Zampini   DMDASNESObjective objectivelocal;
115505f7afSJunchao Zhang 
12*0bf52853SStefano Zampini   /* Vec version for vector data */
13*0bf52853SStefano Zampini   DMDASNESFunctionVec  residuallocalvec;
14*0bf52853SStefano Zampini   DMDASNESJacobianVec  jacobianlocalvec;
15*0bf52853SStefano Zampini   DMDASNESObjectiveVec objectivelocalvec;
16*0bf52853SStefano Zampini 
17*0bf52853SStefano Zampini   /* user contexts */
186cab3a1bSJed Brown   void      *residuallocalctx;
196cab3a1bSJed Brown   void      *jacobianlocalctx;
202a4ee8f2SPeter Brune   void      *objectivelocalctx;
21709ac0d2SJed Brown   InsertMode residuallocalimode;
22dcad5f1cSBarry Smith 
23dcad5f1cSBarry Smith   /*   For Picard iteration defined locally */
24dcad5f1cSBarry Smith   PetscErrorCode (*rhsplocal)(DMDALocalInfo *, void *, void *, void *);
25d1e9a80fSBarry Smith   PetscErrorCode (*jacobianplocal)(DMDALocalInfo *, void *, Mat, Mat, void *);
26dcad5f1cSBarry Smith   void *picardlocalctx;
27942e3340SBarry Smith } DMSNES_DA;
286cab3a1bSJed Brown 
29d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm)
30d71ae5a4SJacob Faibussowitsch {
316cab3a1bSJed Brown   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(PetscFree(sdm->data));
333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
346cab3a1bSJed Brown }
356cab3a1bSJed Brown 
36d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm, DMSNES sdm)
37d71ae5a4SJacob Faibussowitsch {
38e3c0b8a2SPeter Brune   PetscFunctionBegin;
394dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew((DMSNES_DA **)&sdm->data));
401baa6e33SBarry Smith   if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data, oldsdm->data, sizeof(DMSNES_DA)));
413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42e3c0b8a2SPeter Brune }
43e3c0b8a2SPeter Brune 
44d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDASNESGetContext(DM dm, DMSNES sdm, DMSNES_DA **dmdasnes)
45d71ae5a4SJacob Faibussowitsch {
466cab3a1bSJed Brown   PetscFunctionBegin;
470298fd71SBarry Smith   *dmdasnes = NULL;
486cab3a1bSJed Brown   if (!sdm->data) {
494dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMSNES_DA **)&sdm->data));
5022c6f798SBarry Smith     sdm->ops->destroy   = DMSNESDestroy_DMDA;
5122c6f798SBarry Smith     sdm->ops->duplicate = DMSNESDuplicate_DMDA;
526cab3a1bSJed Brown   }
53942e3340SBarry Smith   *dmdasnes = (DMSNES_DA *)sdm->data;
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
556cab3a1bSJed Brown }
566cab3a1bSJed Brown 
57d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeFunction_DMDA(SNES snes, Vec X, Vec F, void *ctx)
58d71ae5a4SJacob Faibussowitsch {
596cab3a1bSJed Brown   DM            dm;
60942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
616cab3a1bSJed Brown   DMDALocalInfo info;
626cab3a1bSJed Brown   Vec           Xloc;
63*0bf52853SStefano Zampini   void         *x, *f, *rctx;
646cab3a1bSJed Brown 
656cab3a1bSJed Brown   PetscFunctionBegin;
662961e153SJed Brown   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
672961e153SJed Brown   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
682961e153SJed Brown   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
695505f7afSJunchao Zhang   PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
709566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
719566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
729566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
739566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
749566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
75*0bf52853SStefano Zampini   rctx = dmdasnes->residuallocalctx ? dmdasnes->residuallocalctx : snes->user;
76709ac0d2SJed Brown   switch (dmdasnes->residuallocalimode) {
77709ac0d2SJed Brown   case INSERT_VALUES: {
789566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_FunctionEval, snes, X, F, 0));
79*0bf52853SStefano Zampini     if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocalvec)(&info, Xloc, F, rctx));
80ef1023bdSBarry Smith     else {
815505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
825505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, F, &f));
83*0bf52853SStefano Zampini       PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocal)(&info, x, f, rctx));
845505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
855505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, F, &f));
865505f7afSJunchao Zhang     }
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));
94*0bf52853SStefano Zampini     if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocalvec)(&info, Xloc, Floc, rctx));
95ef1023bdSBarry Smith     else {
965505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
975505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, Floc, &f));
98*0bf52853SStefano Zampini       PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocal)(&info, x, f, rctx));
995505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1005505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
1015505f7afSJunchao Zhang     }
1029566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_FunctionEval, snes, X, F, 0));
1039566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
1049566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
1059566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
1069566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Floc));
107709ac0d2SJed Brown   } break;
108d71ae5a4SJacob Faibussowitsch   default:
109d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdasnes->residuallocalimode);
110709ac0d2SJed Brown   }
1119566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
1121baa6e33SBarry Smith   if (snes->domainerror) PetscCall(VecSetInf(F));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1146cab3a1bSJed Brown }
1156cab3a1bSJed Brown 
116d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeObjective_DMDA(SNES snes, Vec X, PetscReal *ob, void *ctx)
117d71ae5a4SJacob Faibussowitsch {
1182a4ee8f2SPeter Brune   DM            dm;
119942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
1202a4ee8f2SPeter Brune   DMDALocalInfo info;
1212a4ee8f2SPeter Brune   Vec           Xloc;
122*0bf52853SStefano Zampini   void         *x, *octx;
1232a4ee8f2SPeter Brune 
1242a4ee8f2SPeter Brune   PetscFunctionBegin;
1252a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1262a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
1274f572ea9SToby Isaac   PetscAssertPointer(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));
134*0bf52853SStefano Zampini   octx = dmdasnes->objectivelocalctx ? dmdasnes->objectivelocalctx : snes->user;
135*0bf52853SStefano Zampini   if (dmdasnes->objectivelocalvec) PetscCallBack("SNES DMDA local callback objective", (*dmdasnes->objectivelocalvec)(&info, Xloc, ob, octx));
136ef1023bdSBarry Smith   else {
1379566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Xloc, &x));
138*0bf52853SStefano Zampini     PetscCallBack("SNES DMDA local callback objective", (*dmdasnes->objectivelocal)(&info, x, ob, octx));
1399566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1405505f7afSJunchao Zhang   }
1419566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
142*0bf52853SStefano Zampini   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, ob, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
1433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1442a4ee8f2SPeter Brune }
1452a4ee8f2SPeter Brune 
146cc2e6a90SBarry Smith /* Routine is called by example, hence must be labeled PETSC_EXTERN */
147d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes, Vec X, Mat A, Mat B, void *ctx)
148d71ae5a4SJacob Faibussowitsch {
1492961e153SJed Brown   DM            dm;
150942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
1512961e153SJed Brown   DMDALocalInfo info;
1522961e153SJed Brown   Vec           Xloc;
153*0bf52853SStefano Zampini   void         *x, *jctx;
1542961e153SJed Brown 
1552961e153SJed Brown   PetscFunctionBegin;
1565505f7afSJunchao Zhang   PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
1579566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
158*0bf52853SStefano Zampini   jctx = dmdasnes->jacobianlocalctx ? dmdasnes->jacobianlocalctx : snes->user;
159*0bf52853SStefano Zampini   if (dmdasnes->jacobianlocal || dmdasnes->jacobianlocalvec) {
1609566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Xloc));
1619566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1629566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1639566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(dm, &info));
164*0bf52853SStefano Zampini     if (dmdasnes->jacobianlocalvec) PetscCallBack("SNES DMDA local callback Jacobian", (*dmdasnes->jacobianlocalvec)(&info, Xloc, A, B, jctx));
165ef1023bdSBarry Smith     else {
1669566063dSJacob Faibussowitsch       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
167*0bf52853SStefano Zampini       PetscCallBack("SNES DMDA local callback Jacobian", (*dmdasnes->jacobianlocal)(&info, x, A, B, jctx));
1689566063dSJacob Faibussowitsch       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1695505f7afSJunchao Zhang     }
1709566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Xloc));
1712961e153SJed Brown   } else {
1722961e153SJed Brown     MatFDColoring fdcoloring;
1739566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring));
1742961e153SJed Brown     if (!fdcoloring) {
1752961e153SJed Brown       ISColoring coloring;
1762961e153SJed Brown 
1779566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring));
1789566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring));
1792961e153SJed Brown       switch (dm->coloringtype) {
180d71ae5a4SJacob Faibussowitsch       case IS_COLORING_GLOBAL:
181d71ae5a4SJacob Faibussowitsch         PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))SNESComputeFunction_DMDA, dmdasnes));
182d71ae5a4SJacob Faibussowitsch         break;
183d71ae5a4SJacob Faibussowitsch       default:
184d71ae5a4SJacob Faibussowitsch         SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]);
1852961e153SJed Brown       }
1869566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix));
1879566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1889566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring));
1899566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&coloring));
1909566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring));
1919566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)fdcoloring));
1922961e153SJed Brown 
1932961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
1942961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
1952961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
1962961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
197140e18c1SBarry Smith        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
1982961e153SJed Brown        */
1999566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dm));
2002961e153SJed Brown     }
2019566063dSJacob Faibussowitsch     PetscCall(MatFDColoringApply(B, fdcoloring, X, snes));
2022961e153SJed Brown   }
2032961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
20494ab13aaSBarry Smith   if (A != B) {
2059566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2069566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2072961e153SJed Brown   }
2083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2092961e153SJed Brown }
2102961e153SJed Brown 
2116cab3a1bSJed Brown /*@C
212f6dfbefdSBarry Smith   DMDASNESSetFunctionLocal - set a local residual evaluation function for use with `DMDA`
2136cab3a1bSJed Brown 
2146cab3a1bSJed Brown   Logically Collective
2156cab3a1bSJed Brown 
2164165533cSJose E. Roman   Input Parameters:
217f6dfbefdSBarry Smith + dm    - `DM` to associate callback with
218f6dfbefdSBarry Smith . imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
2196cab3a1bSJed Brown . func  - local residual evaluation
2206cab3a1bSJed Brown - ctx   - optional context for local residual evaluation
2216cab3a1bSJed Brown 
22220f4b53cSBarry Smith   Calling sequence of `func`:
223f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
2240038884cSBarry Smith . x    - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x)
2250038884cSBarry Smith . f    - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f)
2266cab3a1bSJed Brown - ctx  - optional context passed above
2276cab3a1bSJed Brown 
2286cab3a1bSJed Brown   Level: beginner
2296cab3a1bSJed Brown 
230420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetJacobianLocal()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2316cab3a1bSJed Brown @*/
2320b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetFunctionLocal(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *info, void *x, void *f, void *ctx), void *ctx)
233d71ae5a4SJacob Faibussowitsch {
234942e3340SBarry Smith   DMSNES     sdm;
235942e3340SBarry Smith   DMSNES_DA *dmdasnes;
2366cab3a1bSJed Brown 
2376cab3a1bSJed Brown   PetscFunctionBegin;
2386cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2399566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2409566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
2411aa26658SKarl Rupp 
242709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
2436cab3a1bSJed Brown   dmdasnes->residuallocal      = func;
2446cab3a1bSJed Brown   dmdasnes->residuallocalctx   = ctx;
2451aa26658SKarl Rupp 
2469566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
24722c6f798SBarry Smith   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
2489566063dSJacob Faibussowitsch     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
2492961e153SJed Brown   }
2503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2512961e153SJed Brown }
2522961e153SJed Brown 
2532961e153SJed Brown /*@C
254f6dfbefdSBarry Smith   DMDASNESSetFunctionLocalVec - set a local residual evaluation function that operates on a local vector for `DMDA`
2555505f7afSJunchao Zhang 
2565505f7afSJunchao Zhang   Logically Collective
2575505f7afSJunchao Zhang 
2585505f7afSJunchao Zhang   Input Parameters:
259f6dfbefdSBarry Smith + dm    - `DM` to associate callback with
260f6dfbefdSBarry Smith . imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
2615505f7afSJunchao Zhang . func  - local residual evaluation
2625505f7afSJunchao Zhang - ctx   - optional context for local residual evaluation
2635505f7afSJunchao Zhang 
26420f4b53cSBarry Smith   Calling sequence of `func`:
265f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
2665505f7afSJunchao Zhang . x    - state vector at which to evaluate residual
2675505f7afSJunchao Zhang . f    - residual vector
2685505f7afSJunchao Zhang - ctx  - optional context passed above
2695505f7afSJunchao Zhang 
2705505f7afSJunchao Zhang   Level: beginner
2715505f7afSJunchao Zhang 
272420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetFunctionLocal()`, `DMDASNESSetJacobianLocalVec()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2735505f7afSJunchao Zhang @*/
2740b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetFunctionLocalVec(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *info, Vec x, Vec f, void *ctx), void *ctx)
275d71ae5a4SJacob Faibussowitsch {
2765505f7afSJunchao Zhang   DMSNES     sdm;
2775505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
2785505f7afSJunchao Zhang 
2795505f7afSJunchao Zhang   PetscFunctionBegin;
2805505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2815505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2825505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
2835505f7afSJunchao Zhang 
2845505f7afSJunchao Zhang   dmdasnes->residuallocalimode = imode;
2855505f7afSJunchao Zhang   dmdasnes->residuallocalvec   = func;
2865505f7afSJunchao Zhang   dmdasnes->residuallocalctx   = ctx;
2875505f7afSJunchao Zhang 
2885505f7afSJunchao Zhang   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
2895505f7afSJunchao Zhang   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
2905505f7afSJunchao Zhang     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
2915505f7afSJunchao Zhang   }
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2935505f7afSJunchao Zhang }
2945505f7afSJunchao Zhang 
2955505f7afSJunchao Zhang /*@C
296f6dfbefdSBarry Smith   DMDASNESSetJacobianLocal - set a local Jacobian evaluation function for use with `DMDA`
2972961e153SJed Brown 
2982961e153SJed Brown   Logically Collective
2992961e153SJed Brown 
3004165533cSJose E. Roman   Input Parameters:
301f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
3026b7eb5ceSMatthew G. Knepley . func - local Jacobian evaluation
3036b7eb5ceSMatthew G. Knepley - ctx  - optional context for local Jacobian evaluation
3042961e153SJed Brown 
30520f4b53cSBarry Smith   Calling sequence of `func`:
306f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
3070038884cSBarry Smith . x    - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x)
3086b7eb5ceSMatthew G. Knepley . J    - Mat object for the Jacobian
30920f4b53cSBarry Smith . M    - Mat object for the Jacobian preconditioner matrix, often `J`
3102961e153SJed Brown - ctx  - optional context passed above
3112961e153SJed Brown 
3122961e153SJed Brown   Level: beginner
3132961e153SJed Brown 
314420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetFunctionLocal()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3152961e153SJed Brown @*/
3160b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetJacobianLocal(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, void *x, Mat J, Mat M, void *ctx), void *ctx)
317d71ae5a4SJacob Faibussowitsch {
318942e3340SBarry Smith   DMSNES     sdm;
319942e3340SBarry Smith   DMSNES_DA *dmdasnes;
3202961e153SJed Brown 
3212961e153SJed Brown   PetscFunctionBegin;
3222961e153SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3239566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3249566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
3251aa26658SKarl Rupp 
3262961e153SJed Brown   dmdasnes->jacobianlocal    = func;
3272961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
3281aa26658SKarl Rupp 
3299566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
3303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3316cab3a1bSJed Brown }
3322a4ee8f2SPeter Brune 
3332a4ee8f2SPeter Brune /*@C
334f6dfbefdSBarry Smith   DMDASNESSetJacobianLocalVec - set a local Jacobian evaluation function that operates on a local vector with `DMDA`
3355505f7afSJunchao Zhang 
3365505f7afSJunchao Zhang   Logically Collective
3375505f7afSJunchao Zhang 
3385505f7afSJunchao Zhang   Input Parameters:
339f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
3405505f7afSJunchao Zhang . func - local Jacobian evaluation
3415505f7afSJunchao Zhang - ctx  - optional context for local Jacobian evaluation
3425505f7afSJunchao Zhang 
34320f4b53cSBarry Smith   Calling sequence of `func`:
344f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
3455505f7afSJunchao Zhang . x    - state vector at which to evaluate Jacobian
3460b4db180SJacob Faibussowitsch . J    - the Jacobian
3470b4db180SJacob Faibussowitsch . M    - approximate Jacobian from which the preconditioner will be computed, often `J`
3485505f7afSJunchao Zhang - ctx  - optional context passed above
3495505f7afSJunchao Zhang 
3505505f7afSJunchao Zhang   Level: beginner
3515505f7afSJunchao Zhang 
352420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetJacobianLocal()`, `DMDASNESSetFunctionLocalVec()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3535505f7afSJunchao Zhang @*/
3540b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetJacobianLocalVec(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, Vec x, Mat J, Mat M, void *), void *ctx)
355d71ae5a4SJacob Faibussowitsch {
3565505f7afSJunchao Zhang   DMSNES     sdm;
3575505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
3585505f7afSJunchao Zhang 
3595505f7afSJunchao Zhang   PetscFunctionBegin;
3605505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3615505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3625505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
3635505f7afSJunchao Zhang 
3645505f7afSJunchao Zhang   dmdasnes->jacobianlocalvec = func;
3655505f7afSJunchao Zhang   dmdasnes->jacobianlocalctx = ctx;
3665505f7afSJunchao Zhang 
3675505f7afSJunchao Zhang   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
3683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3695505f7afSJunchao Zhang }
3705505f7afSJunchao Zhang 
3715505f7afSJunchao Zhang /*@C
372f6dfbefdSBarry Smith   DMDASNESSetObjectiveLocal - set a local residual evaluation function to used with a `DMDA`
3732a4ee8f2SPeter Brune 
3742a4ee8f2SPeter Brune   Logically Collective
3752a4ee8f2SPeter Brune 
3764165533cSJose E. Roman   Input Parameters:
377f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
378420bcc1bSBarry Smith . func - local objective evaluation, see `DMDASNESSetObjectiveLocal` for the calling sequence
3792a4ee8f2SPeter Brune - ctx  - optional context for local residual evaluation
3802a4ee8f2SPeter Brune 
381*0bf52853SStefano Zampini   Calling sequence of `func`:
382*0bf52853SStefano Zampini + info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
383*0bf52853SStefano Zampini . x    - dimensional pointer to state at which to evaluate the objective (e.g. PetscScalar *x or **x or ***x)
384*0bf52853SStefano Zampini . obj  - returned objective value for the local subdomain
385*0bf52853SStefano Zampini - ctx  - optional context passed above
386*0bf52853SStefano Zampini 
3872a4ee8f2SPeter Brune   Level: beginner
3882a4ee8f2SPeter Brune 
389420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDASNESObjective`
3902a4ee8f2SPeter Brune @*/
391*0bf52853SStefano Zampini PetscErrorCode DMDASNESSetObjectiveLocal(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, void *x, PetscReal *obj, void *), void *ctx)
392d71ae5a4SJacob Faibussowitsch {
393942e3340SBarry Smith   DMSNES     sdm;
394942e3340SBarry Smith   DMSNES_DA *dmdasnes;
3952a4ee8f2SPeter Brune 
3962a4ee8f2SPeter Brune   PetscFunctionBegin;
3972a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3989566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3999566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
4001aa26658SKarl Rupp 
4012a4ee8f2SPeter Brune   dmdasnes->objectivelocal    = func;
4022a4ee8f2SPeter Brune   dmdasnes->objectivelocalctx = ctx;
4031aa26658SKarl Rupp 
4049566063dSJacob Faibussowitsch   PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes));
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4062a4ee8f2SPeter Brune }
407dcad5f1cSBarry Smith 
4085505f7afSJunchao Zhang /*@C
409e4094ef1SJacob Faibussowitsch   DMDASNESSetObjectiveLocalVec - set a local residual evaluation function that operates on a local vector with `DMDA`
4105505f7afSJunchao Zhang 
4115505f7afSJunchao Zhang   Logically Collective
4125505f7afSJunchao Zhang 
4135505f7afSJunchao Zhang   Input Parameters:
414f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
415420bcc1bSBarry Smith . func - local objective evaluation, see `DMDASNESSetObjectiveLocalVec` for the calling sequence
4165505f7afSJunchao Zhang - ctx  - optional context for local residual evaluation
4175505f7afSJunchao Zhang 
418*0bf52853SStefano Zampini   Calling sequence of `func`:
419*0bf52853SStefano Zampini + info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
420*0bf52853SStefano Zampini . x    - state vector at which to evaluate the objective
421*0bf52853SStefano Zampini . obj  - returned objective value for the local subdomain
422*0bf52853SStefano Zampini - ctx  - optional context passed above
423*0bf52853SStefano Zampini 
4245505f7afSJunchao Zhang   Level: beginner
4255505f7afSJunchao Zhang 
426420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetObjectiveLocal()`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocalVec()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDASNESObjectiveVec`
4275505f7afSJunchao Zhang @*/
428*0bf52853SStefano Zampini PetscErrorCode DMDASNESSetObjectiveLocalVec(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, Vec x, PetscReal *obj, void *), void *ctx)
429d71ae5a4SJacob Faibussowitsch {
4305505f7afSJunchao Zhang   DMSNES     sdm;
4315505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
4325505f7afSJunchao Zhang 
4335505f7afSJunchao Zhang   PetscFunctionBegin;
4345505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4355505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
4365505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
4375505f7afSJunchao Zhang 
4385505f7afSJunchao Zhang   dmdasnes->objectivelocalvec = func;
4395505f7afSJunchao Zhang   dmdasnes->objectivelocalctx = ctx;
4405505f7afSJunchao Zhang 
4415505f7afSJunchao Zhang   PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes));
4423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4435505f7afSJunchao Zhang }
4445505f7afSJunchao Zhang 
445d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputePicard_DMDA(SNES snes, Vec X, Vec F, void *ctx)
446d71ae5a4SJacob Faibussowitsch {
447dcad5f1cSBarry Smith   DM            dm;
448942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
449dcad5f1cSBarry Smith   DMDALocalInfo info;
450dcad5f1cSBarry Smith   Vec           Xloc;
451dcad5f1cSBarry Smith   void         *x, *f;
452dcad5f1cSBarry Smith 
453dcad5f1cSBarry Smith   PetscFunctionBegin;
454dcad5f1cSBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
455dcad5f1cSBarry Smith   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
456dcad5f1cSBarry Smith   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
45728b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->rhsplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
4589566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
4599566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
4609566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
4619566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
4629566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
4639566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
464dcad5f1cSBarry Smith   switch (dmdasnes->residuallocalimode) {
465dcad5f1cSBarry Smith   case INSERT_VALUES: {
4669566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, F, &f));
467792fecdfSBarry Smith     PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx));
4689566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, F, &f));
469dcad5f1cSBarry Smith   } break;
470dcad5f1cSBarry Smith   case ADD_VALUES: {
471dcad5f1cSBarry Smith     Vec Floc;
4729566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Floc));
4739566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
4749566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Floc, &f));
475792fecdfSBarry Smith     PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx));
4769566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
4779566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
4789566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
4799566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
4809566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Floc));
481dcad5f1cSBarry Smith   } break;
482d71ae5a4SJacob Faibussowitsch   default:
483d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdasnes->residuallocalimode);
484dcad5f1cSBarry Smith   }
4859566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
4869566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
4873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
488dcad5f1cSBarry Smith }
489dcad5f1cSBarry Smith 
490d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes, Vec X, Mat A, Mat B, void *ctx)
491d71ae5a4SJacob Faibussowitsch {
492dcad5f1cSBarry Smith   DM            dm;
493942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
494dcad5f1cSBarry Smith   DMDALocalInfo info;
495dcad5f1cSBarry Smith   Vec           Xloc;
496dcad5f1cSBarry Smith   void         *x;
497dcad5f1cSBarry Smith 
498dcad5f1cSBarry Smith   PetscFunctionBegin;
49928b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->jacobianplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
5009566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
501dcad5f1cSBarry Smith 
5029566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
5039566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
5049566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
5059566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
5069566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
507792fecdfSBarry Smith   PetscCallBack("SNES Picard DMDA local callback Jacobian", (*dmdasnes->jacobianplocal)(&info, x, A, B, dmdasnes->picardlocalctx));
5089566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
5099566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
51094ab13aaSBarry Smith   if (A != B) {
5119566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5129566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
513dcad5f1cSBarry Smith   }
5143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
515dcad5f1cSBarry Smith }
516dcad5f1cSBarry Smith 
517dcad5f1cSBarry Smith /*@C
518f6dfbefdSBarry Smith   DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration with `DMDA`
519dcad5f1cSBarry Smith 
520dcad5f1cSBarry Smith   Logically Collective
521dcad5f1cSBarry Smith 
5224165533cSJose E. Roman   Input Parameters:
523f6dfbefdSBarry Smith + dm    - `DM` to associate callback with
524f6dfbefdSBarry Smith . imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
525dcad5f1cSBarry Smith . func  - local residual evaluation
526ceaaa498SBarry Smith . jac   - function to compute Jacobian
527dcad5f1cSBarry Smith - ctx   - optional context for local residual evaluation
528dcad5f1cSBarry Smith 
52920f4b53cSBarry Smith   Calling sequence of `func`:
53020f4b53cSBarry Smith + info - defines the subdomain to evaluate the residual on
531dcad5f1cSBarry Smith . x    - dimensional pointer to state at which to evaluate residual
5320b4db180SJacob Faibussowitsch . f    - dimensional pointer to residual, write the residual here
5330b4db180SJacob Faibussowitsch - ctx  - optional context passed above
5340b4db180SJacob Faibussowitsch 
5350b4db180SJacob Faibussowitsch   Calling sequence of `jac`:
5360b4db180SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
5370b4db180SJacob Faibussowitsch . x    - dimensional pointer to state at which to evaluate residual
5380b4db180SJacob Faibussowitsch . jac  - the Jacobian
5390b4db180SJacob Faibussowitsch . Jp   - approximation to the Jacobian used to compute the preconditioner, often `J`
540dcad5f1cSBarry Smith - ctx  - optional context passed above
541dcad5f1cSBarry Smith 
542dcad5f1cSBarry Smith   Level: beginner
543dcad5f1cSBarry Smith 
54420f4b53cSBarry Smith   Note:
545420bcc1bSBarry Smith   The user must use `SNESSetFunction`(`snes`,`NULL`,`SNESPicardComputeFunction`,&user));
54620f4b53cSBarry Smith   in their code before calling this routine.
54720f4b53cSBarry Smith 
548420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
549dcad5f1cSBarry Smith @*/
5500b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetPicardLocal(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *info, void *x, void *f, void *ctx), PetscErrorCode (*jac)(DMDALocalInfo *info, void *x, Mat jac, Mat Jp, void *ctx), void *ctx)
551d71ae5a4SJacob Faibussowitsch {
552942e3340SBarry Smith   DMSNES     sdm;
553942e3340SBarry Smith   DMSNES_DA *dmdasnes;
554dcad5f1cSBarry Smith 
555dcad5f1cSBarry Smith   PetscFunctionBegin;
556dcad5f1cSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5579566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
5589566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
5591aa26658SKarl Rupp 
560dcad5f1cSBarry Smith   dmdasnes->residuallocalimode = imode;
561dcad5f1cSBarry Smith   dmdasnes->rhsplocal          = func;
562dcad5f1cSBarry Smith   dmdasnes->jacobianplocal     = jac;
563dcad5f1cSBarry Smith   dmdasnes->picardlocalctx     = ctx;
5641aa26658SKarl Rupp 
5659566063dSJacob Faibussowitsch   PetscCall(DMSNESSetPicard(dm, SNESComputePicard_DMDA, SNESComputePicardJacobian_DMDA, dmdasnes));
5669566063dSJacob Faibussowitsch   PetscCall(DMSNESSetMFFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
5673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
568dcad5f1cSBarry Smith }
569