xref: /petsc/src/snes/utils/dmdasnes.c (revision e4094ef18e7e53fda86cf35f3a47fda48a8e77d8)
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 
25d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm)
26d71ae5a4SJacob Faibussowitsch {
276cab3a1bSJed Brown   PetscFunctionBegin;
289566063dSJacob Faibussowitsch   PetscCall(PetscFree(sdm->data));
293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
306cab3a1bSJed Brown }
316cab3a1bSJed Brown 
32d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm, DMSNES sdm)
33d71ae5a4SJacob Faibussowitsch {
34e3c0b8a2SPeter Brune   PetscFunctionBegin;
354dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew((DMSNES_DA **)&sdm->data));
361baa6e33SBarry Smith   if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data, oldsdm->data, sizeof(DMSNES_DA)));
373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38e3c0b8a2SPeter Brune }
39e3c0b8a2SPeter Brune 
40d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDASNESGetContext(DM dm, DMSNES sdm, DMSNES_DA **dmdasnes)
41d71ae5a4SJacob Faibussowitsch {
426cab3a1bSJed Brown   PetscFunctionBegin;
430298fd71SBarry Smith   *dmdasnes = NULL;
446cab3a1bSJed Brown   if (!sdm->data) {
454dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((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;
503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
516cab3a1bSJed Brown }
526cab3a1bSJed Brown 
53d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeFunction_DMDA(SNES snes, Vec X, Vec F, void *ctx)
54d71ae5a4SJacob Faibussowitsch {
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));
74792fecdfSBarry Smith     if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocalvec)(&info, Xloc, F, dmdasnes->residuallocalctx));
75ef1023bdSBarry Smith     else {
765505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
775505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, F, &f));
78792fecdfSBarry Smith       PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocal)(&info, x, f, dmdasnes->residuallocalctx));
795505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
805505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, F, &f));
815505f7afSJunchao Zhang     }
829566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_FunctionEval, snes, X, F, 0));
83709ac0d2SJed Brown   } break;
84709ac0d2SJed Brown   case ADD_VALUES: {
85709ac0d2SJed Brown     Vec Floc;
869566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Floc));
879566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
889566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_FunctionEval, snes, X, F, 0));
89792fecdfSBarry Smith     if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocalvec)(&info, Xloc, Floc, dmdasnes->residuallocalctx));
90ef1023bdSBarry Smith     else {
915505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
925505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, Floc, &f));
93792fecdfSBarry Smith       PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocal)(&info, x, f, dmdasnes->residuallocalctx));
945505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
955505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
965505f7afSJunchao Zhang     }
979566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_FunctionEval, snes, X, F, 0));
989566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
999566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
1009566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
1019566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Floc));
102709ac0d2SJed Brown   } break;
103d71ae5a4SJacob Faibussowitsch   default:
104d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdasnes->residuallocalimode);
105709ac0d2SJed Brown   }
1069566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
1071baa6e33SBarry Smith   if (snes->domainerror) PetscCall(VecSetInf(F));
1083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1096cab3a1bSJed Brown }
1106cab3a1bSJed Brown 
111d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeObjective_DMDA(SNES snes, Vec X, PetscReal *ob, void *ctx)
112d71ae5a4SJacob Faibussowitsch {
1132a4ee8f2SPeter Brune   DM            dm;
114942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
1152a4ee8f2SPeter Brune   DMDALocalInfo info;
1162a4ee8f2SPeter Brune   Vec           Xloc;
1172a4ee8f2SPeter Brune   void         *x;
1182a4ee8f2SPeter Brune 
1192a4ee8f2SPeter Brune   PetscFunctionBegin;
1202a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1212a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
122dadcf809SJacob Faibussowitsch   PetscValidRealPointer(ob, 3);
1235505f7afSJunchao Zhang   PetscCheck(dmdasnes->objectivelocal || dmdasnes->objectivelocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
1249566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
1259566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
1269566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1279566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1289566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
129792fecdfSBarry Smith   if (dmdasnes->objectivelocalvec) PetscCallBack("SNES DMDA local callback objective", (*dmdasnes->objectivelocalvec)(&info, Xloc, ob, dmdasnes->objectivelocalctx));
130ef1023bdSBarry Smith   else {
1319566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Xloc, &x));
132792fecdfSBarry Smith     PetscCallBack("SNES DMDA local callback objective", (*dmdasnes->objectivelocal)(&info, x, ob, dmdasnes->objectivelocalctx));
1339566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1345505f7afSJunchao Zhang   }
1359566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
1363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1372a4ee8f2SPeter Brune }
1382a4ee8f2SPeter Brune 
139cc2e6a90SBarry Smith /* Routine is called by example, hence must be labeled PETSC_EXTERN */
140d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes, Vec X, Mat A, Mat B, void *ctx)
141d71ae5a4SJacob Faibussowitsch {
1422961e153SJed Brown   DM            dm;
143942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
1442961e153SJed Brown   DMDALocalInfo info;
1452961e153SJed Brown   Vec           Xloc;
1462961e153SJed Brown   void         *x;
1472961e153SJed Brown 
1482961e153SJed Brown   PetscFunctionBegin;
1495505f7afSJunchao Zhang   PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
1509566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
1512961e153SJed Brown 
1525505f7afSJunchao Zhang   if (dmdasnes->jacobianlocal || dmdasnes->jacobianlocalctx) {
1539566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Xloc));
1549566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1559566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1569566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(dm, &info));
157792fecdfSBarry Smith     if (dmdasnes->jacobianlocalvec) PetscCallBack("SNES DMDA local callback Jacobian", (*dmdasnes->jacobianlocalvec)(&info, Xloc, A, B, dmdasnes->jacobianlocalctx));
158ef1023bdSBarry Smith     else {
1599566063dSJacob Faibussowitsch       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
160792fecdfSBarry Smith       PetscCallBack("SNES DMDA local callback Jacobian", (*dmdasnes->jacobianlocal)(&info, x, A, B, dmdasnes->jacobianlocalctx));
1619566063dSJacob Faibussowitsch       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1625505f7afSJunchao Zhang     }
1639566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Xloc));
1642961e153SJed Brown   } else {
1652961e153SJed Brown     MatFDColoring fdcoloring;
1669566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring));
1672961e153SJed Brown     if (!fdcoloring) {
1682961e153SJed Brown       ISColoring coloring;
1692961e153SJed Brown 
1709566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring));
1719566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring));
1722961e153SJed Brown       switch (dm->coloringtype) {
173d71ae5a4SJacob Faibussowitsch       case IS_COLORING_GLOBAL:
174d71ae5a4SJacob Faibussowitsch         PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))SNESComputeFunction_DMDA, dmdasnes));
175d71ae5a4SJacob Faibussowitsch         break;
176d71ae5a4SJacob Faibussowitsch       default:
177d71ae5a4SJacob Faibussowitsch         SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]);
1782961e153SJed Brown       }
1799566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix));
1809566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1819566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring));
1829566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&coloring));
1839566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring));
1849566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)fdcoloring));
1852961e153SJed Brown 
1862961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
1872961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
1882961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
1892961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
190140e18c1SBarry Smith        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
1912961e153SJed Brown        */
1929566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dm));
1932961e153SJed Brown     }
1949566063dSJacob Faibussowitsch     PetscCall(MatFDColoringApply(B, fdcoloring, X, snes));
1952961e153SJed Brown   }
1962961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
19794ab13aaSBarry Smith   if (A != B) {
1989566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1999566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2002961e153SJed Brown   }
2013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2022961e153SJed Brown }
2032961e153SJed Brown 
2046cab3a1bSJed Brown /*@C
205f6dfbefdSBarry Smith   DMDASNESSetFunctionLocal - set a local residual evaluation function for use with `DMDA`
2066cab3a1bSJed Brown 
2076cab3a1bSJed Brown   Logically Collective
2086cab3a1bSJed Brown 
2094165533cSJose E. Roman   Input Parameters:
210f6dfbefdSBarry Smith + dm    - `DM` to associate callback with
211f6dfbefdSBarry Smith . imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
2126cab3a1bSJed Brown . func  - local residual evaluation
2136cab3a1bSJed Brown - ctx   - optional context for local residual evaluation
2146cab3a1bSJed Brown 
21520f4b53cSBarry Smith   Calling sequence of `func`:
21620f4b53cSBarry Smith $   PetscErrorCode func(DMDALocalInfo *info, void *x, void *f, void *ctx)
217f6dfbefdSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
2180038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x)
2190038884cSBarry Smith .  f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f)
2206cab3a1bSJed Brown - ctx - optional context passed above
2216cab3a1bSJed Brown 
2226cab3a1bSJed Brown   Level: beginner
2236cab3a1bSJed Brown 
224f6dfbefdSBarry Smith .seealso: `DMDA`, `DMDASNESSetJacobianLocal()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2256cab3a1bSJed Brown @*/
226d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASNESSetFunctionLocal(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *, void *, void *, void *), void *ctx)
227d71ae5a4SJacob Faibussowitsch {
228942e3340SBarry Smith   DMSNES     sdm;
229942e3340SBarry Smith   DMSNES_DA *dmdasnes;
2306cab3a1bSJed Brown 
2316cab3a1bSJed Brown   PetscFunctionBegin;
2326cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2339566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2349566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
2351aa26658SKarl Rupp 
236709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
2376cab3a1bSJed Brown   dmdasnes->residuallocal      = func;
2386cab3a1bSJed Brown   dmdasnes->residuallocalctx   = ctx;
2391aa26658SKarl Rupp 
2409566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
24122c6f798SBarry Smith   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
2429566063dSJacob Faibussowitsch     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
2432961e153SJed Brown   }
2443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2452961e153SJed Brown }
2462961e153SJed Brown 
2472961e153SJed Brown /*@C
248f6dfbefdSBarry Smith   DMDASNESSetFunctionLocalVec - set a local residual evaluation function that operates on a local vector for `DMDA`
2495505f7afSJunchao Zhang 
2505505f7afSJunchao Zhang   Logically Collective
2515505f7afSJunchao Zhang 
2525505f7afSJunchao Zhang   Input Parameters:
253f6dfbefdSBarry Smith + dm    - `DM` to associate callback with
254f6dfbefdSBarry Smith . imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
2555505f7afSJunchao Zhang . func  - local residual evaluation
2565505f7afSJunchao Zhang - ctx   - optional context for local residual evaluation
2575505f7afSJunchao Zhang 
25820f4b53cSBarry Smith   Calling sequence of `func`:
25920f4b53cSBarry Smith $   PetscErrorCode func(DMDALocalInfo *info, Vec x, Vec f, void *ctx),
260f6dfbefdSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
2615505f7afSJunchao Zhang .  x - state vector at which to evaluate residual
2625505f7afSJunchao Zhang .  f - residual vector
2635505f7afSJunchao Zhang - ctx - optional context passed above
2645505f7afSJunchao Zhang 
2655505f7afSJunchao Zhang   Level: beginner
2665505f7afSJunchao Zhang 
267f6dfbefdSBarry Smith .seealso: `DMDA`, `DMDASNESSetFunctionLocal()`, `DMDASNESSetJacobianLocalVec()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2685505f7afSJunchao Zhang @*/
269d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASNESSetFunctionLocalVec(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *, Vec, Vec, void *), void *ctx)
270d71ae5a4SJacob Faibussowitsch {
2715505f7afSJunchao Zhang   DMSNES     sdm;
2725505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
2735505f7afSJunchao Zhang 
2745505f7afSJunchao Zhang   PetscFunctionBegin;
2755505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2765505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2775505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
2785505f7afSJunchao Zhang 
2795505f7afSJunchao Zhang   dmdasnes->residuallocalimode = imode;
2805505f7afSJunchao Zhang   dmdasnes->residuallocalvec   = func;
2815505f7afSJunchao Zhang   dmdasnes->residuallocalctx   = ctx;
2825505f7afSJunchao Zhang 
2835505f7afSJunchao Zhang   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
2845505f7afSJunchao Zhang   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
2855505f7afSJunchao Zhang     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
2865505f7afSJunchao Zhang   }
2873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2885505f7afSJunchao Zhang }
2895505f7afSJunchao Zhang 
2905505f7afSJunchao Zhang /*@C
291f6dfbefdSBarry Smith   DMDASNESSetJacobianLocal - set a local Jacobian evaluation function for use with `DMDA`
2922961e153SJed Brown 
2932961e153SJed Brown   Logically Collective
2942961e153SJed Brown 
2954165533cSJose E. Roman   Input Parameters:
296f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
2976b7eb5ceSMatthew G. Knepley . func - local Jacobian evaluation
2986b7eb5ceSMatthew G. Knepley - ctx  - optional context for local Jacobian evaluation
2992961e153SJed Brown 
30020f4b53cSBarry Smith   Calling sequence of `func`:
30120f4b53cSBarry Smith $  PetscErrorCode func(DMDALocalInfo *info, void *x, Mat J, Mat M, void *ctx),
302f6dfbefdSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
3030038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x)
3046b7eb5ceSMatthew G. Knepley .  J - Mat object for the Jacobian
30520f4b53cSBarry Smith .  M - Mat object for the Jacobian preconditioner matrix, often `J`
3062961e153SJed Brown - ctx - optional context passed above
3072961e153SJed Brown 
3082961e153SJed Brown   Level: beginner
3092961e153SJed Brown 
310f6dfbefdSBarry Smith .seealso: `DMDA`, `DMDASNESSetFunctionLocal()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3112961e153SJed Brown @*/
312d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASNESSetJacobianLocal(DM dm, PetscErrorCode (*func)(DMDALocalInfo *, void *, Mat, Mat, void *), void *ctx)
313d71ae5a4SJacob Faibussowitsch {
314942e3340SBarry Smith   DMSNES     sdm;
315942e3340SBarry Smith   DMSNES_DA *dmdasnes;
3162961e153SJed Brown 
3172961e153SJed Brown   PetscFunctionBegin;
3182961e153SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3199566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3209566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
3211aa26658SKarl Rupp 
3222961e153SJed Brown   dmdasnes->jacobianlocal    = func;
3232961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
3241aa26658SKarl Rupp 
3259566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
3263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3276cab3a1bSJed Brown }
3282a4ee8f2SPeter Brune 
3292a4ee8f2SPeter Brune /*@C
330f6dfbefdSBarry Smith   DMDASNESSetJacobianLocalVec - set a local Jacobian evaluation function that operates on a local vector with `DMDA`
3315505f7afSJunchao Zhang 
3325505f7afSJunchao Zhang   Logically Collective
3335505f7afSJunchao Zhang 
3345505f7afSJunchao Zhang   Input Parameters:
335f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
3365505f7afSJunchao Zhang . func - local Jacobian evaluation
3375505f7afSJunchao Zhang - ctx  - optional context for local Jacobian evaluation
3385505f7afSJunchao Zhang 
33920f4b53cSBarry Smith   Calling sequence of `func`:
34020f4b53cSBarry Smith $  PetscErrorCode (*func)(DMDALocalInfo *info, Vec x, Mat J, Mat M, void *ctx),
341f6dfbefdSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
3425505f7afSJunchao Zhang .  x - state vector at which to evaluate Jacobian
3435505f7afSJunchao Zhang .  J - Mat object for the Jacobian
34420f4b53cSBarry Smith .  M - Mat object for the Jacobian preconditioner matrix, often `J`
3455505f7afSJunchao Zhang - ctx - optional context passed above
3465505f7afSJunchao Zhang 
3475505f7afSJunchao Zhang   Level: beginner
3485505f7afSJunchao Zhang 
349f6dfbefdSBarry Smith .seealso: `DMDA`, `DMDASNESSetJacobianLocal()`, `DMDASNESSetFunctionLocalVec()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3505505f7afSJunchao Zhang @*/
351d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASNESSetJacobianLocalVec(DM dm, PetscErrorCode (*func)(DMDALocalInfo *, Vec, Mat, Mat, void *), void *ctx)
352d71ae5a4SJacob Faibussowitsch {
3535505f7afSJunchao Zhang   DMSNES     sdm;
3545505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
3555505f7afSJunchao Zhang 
3565505f7afSJunchao Zhang   PetscFunctionBegin;
3575505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3585505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3595505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
3605505f7afSJunchao Zhang 
3615505f7afSJunchao Zhang   dmdasnes->jacobianlocalvec = func;
3625505f7afSJunchao Zhang   dmdasnes->jacobianlocalctx = ctx;
3635505f7afSJunchao Zhang 
3645505f7afSJunchao Zhang   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
3653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3665505f7afSJunchao Zhang }
3675505f7afSJunchao Zhang 
3685505f7afSJunchao Zhang /*@C
369f6dfbefdSBarry Smith   DMDASNESSetObjectiveLocal - set a local residual evaluation function to used with a `DMDA`
3702a4ee8f2SPeter Brune 
3712a4ee8f2SPeter Brune   Logically Collective
3722a4ee8f2SPeter Brune 
3734165533cSJose E. Roman   Input Parameters:
374f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
3752a4ee8f2SPeter Brune . func - local objective evaluation
3762a4ee8f2SPeter Brune - ctx  - optional context for local residual evaluation
3772a4ee8f2SPeter Brune 
37820f4b53cSBarry Smith   Calling sequence of `func`:
37920f4b53cSBarry Smith $  PetscErrorCode func(DMDALocalInfo *info, void *x, PetscReal obj, void *ctx);
38020f4b53cSBarry Smith +  info - defines the subdomain to evaluate the residual on
3812a4ee8f2SPeter Brune .  x - dimensional pointer to state at which to evaluate residual
3822a4ee8f2SPeter Brune .  ob - eventual objective value
3832a4ee8f2SPeter Brune - ctx - optional context passed above
3842a4ee8f2SPeter Brune 
3852a4ee8f2SPeter Brune   Level: beginner
3862a4ee8f2SPeter Brune 
387f6dfbefdSBarry Smith .seealso: `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3882a4ee8f2SPeter Brune @*/
389d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASNESSetObjectiveLocal(DM dm, DMDASNESObjective func, void *ctx)
390d71ae5a4SJacob Faibussowitsch {
391942e3340SBarry Smith   DMSNES     sdm;
392942e3340SBarry Smith   DMSNES_DA *dmdasnes;
3932a4ee8f2SPeter Brune 
3942a4ee8f2SPeter Brune   PetscFunctionBegin;
3952a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3969566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3979566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
3981aa26658SKarl Rupp 
3992a4ee8f2SPeter Brune   dmdasnes->objectivelocal    = func;
4002a4ee8f2SPeter Brune   dmdasnes->objectivelocalctx = ctx;
4011aa26658SKarl Rupp 
4029566063dSJacob Faibussowitsch   PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes));
4033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4042a4ee8f2SPeter Brune }
405dcad5f1cSBarry Smith 
4065505f7afSJunchao Zhang /*@C
407*e4094ef1SJacob Faibussowitsch   DMDASNESSetObjectiveLocalVec - set a local residual evaluation function that operates on a local vector with `DMDA`
4085505f7afSJunchao Zhang 
4095505f7afSJunchao Zhang   Logically Collective
4105505f7afSJunchao Zhang 
4115505f7afSJunchao Zhang   Input Parameters:
412f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
4135505f7afSJunchao Zhang . func - local objective evaluation
4145505f7afSJunchao Zhang - ctx  - optional context for local residual evaluation
4155505f7afSJunchao Zhang 
41620f4b53cSBarry Smith   Calling sequence of `func`:
41720f4b53cSBarry Smith $  PetscErrorCode func(DMDALocalInfo *info, Vec x, PetscReal *ob, void *ctx);
41820f4b53cSBarry Smith +  info - defines the subdomain to evaluate the residual on
4195505f7afSJunchao Zhang .  x - state vector at which to evaluate residual
4205505f7afSJunchao Zhang .  ob - eventual objective value
4215505f7afSJunchao Zhang - ctx - optional context passed above
4225505f7afSJunchao Zhang 
4235505f7afSJunchao Zhang   Level: beginner
4245505f7afSJunchao Zhang 
425f6dfbefdSBarry Smith .seealso: `DMDA`, `DMDASNESSetObjectiveLocal()`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocalVec()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
4265505f7afSJunchao Zhang @*/
427d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASNESSetObjectiveLocalVec(DM dm, DMDASNESObjectiveVec func, void *ctx)
428d71ae5a4SJacob Faibussowitsch {
4295505f7afSJunchao Zhang   DMSNES     sdm;
4305505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
4315505f7afSJunchao Zhang 
4325505f7afSJunchao Zhang   PetscFunctionBegin;
4335505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4345505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
4355505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
4365505f7afSJunchao Zhang 
4375505f7afSJunchao Zhang   dmdasnes->objectivelocalvec = func;
4385505f7afSJunchao Zhang   dmdasnes->objectivelocalctx = ctx;
4395505f7afSJunchao Zhang 
4405505f7afSJunchao Zhang   PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes));
4413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4425505f7afSJunchao Zhang }
4435505f7afSJunchao Zhang 
444d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputePicard_DMDA(SNES snes, Vec X, Vec F, void *ctx)
445d71ae5a4SJacob Faibussowitsch {
446dcad5f1cSBarry Smith   DM            dm;
447942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
448dcad5f1cSBarry Smith   DMDALocalInfo info;
449dcad5f1cSBarry Smith   Vec           Xloc;
450dcad5f1cSBarry Smith   void         *x, *f;
451dcad5f1cSBarry Smith 
452dcad5f1cSBarry Smith   PetscFunctionBegin;
453dcad5f1cSBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
454dcad5f1cSBarry Smith   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
455dcad5f1cSBarry Smith   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
45628b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->rhsplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
4579566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
4589566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
4599566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
4609566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
4619566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
4629566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
463dcad5f1cSBarry Smith   switch (dmdasnes->residuallocalimode) {
464dcad5f1cSBarry Smith   case INSERT_VALUES: {
4659566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, F, &f));
466792fecdfSBarry Smith     PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx));
4679566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, F, &f));
468dcad5f1cSBarry Smith   } break;
469dcad5f1cSBarry Smith   case ADD_VALUES: {
470dcad5f1cSBarry Smith     Vec Floc;
4719566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Floc));
4729566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
4739566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Floc, &f));
474792fecdfSBarry Smith     PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx));
4759566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
4769566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
4779566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
4789566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
4799566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Floc));
480dcad5f1cSBarry Smith   } break;
481d71ae5a4SJacob Faibussowitsch   default:
482d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdasnes->residuallocalimode);
483dcad5f1cSBarry Smith   }
4849566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
4859566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
4863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
487dcad5f1cSBarry Smith }
488dcad5f1cSBarry Smith 
489d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes, Vec X, Mat A, Mat B, void *ctx)
490d71ae5a4SJacob Faibussowitsch {
491dcad5f1cSBarry Smith   DM            dm;
492942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
493dcad5f1cSBarry Smith   DMDALocalInfo info;
494dcad5f1cSBarry Smith   Vec           Xloc;
495dcad5f1cSBarry Smith   void         *x;
496dcad5f1cSBarry Smith 
497dcad5f1cSBarry Smith   PetscFunctionBegin;
49828b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->jacobianplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
4999566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
500dcad5f1cSBarry Smith 
5019566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
5029566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
5039566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
5049566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
5059566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
506792fecdfSBarry Smith   PetscCallBack("SNES Picard DMDA local callback Jacobian", (*dmdasnes->jacobianplocal)(&info, x, A, B, dmdasnes->picardlocalctx));
5079566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
5089566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
50994ab13aaSBarry Smith   if (A != B) {
5109566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5119566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
512dcad5f1cSBarry Smith   }
5133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
514dcad5f1cSBarry Smith }
515dcad5f1cSBarry Smith 
516dcad5f1cSBarry Smith /*@C
517f6dfbefdSBarry Smith   DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration with `DMDA`
518dcad5f1cSBarry Smith 
519dcad5f1cSBarry Smith   Logically Collective
520dcad5f1cSBarry Smith 
5214165533cSJose E. Roman   Input Parameters:
522f6dfbefdSBarry Smith + dm    - `DM` to associate callback with
523f6dfbefdSBarry Smith . imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
524dcad5f1cSBarry Smith . func  - local residual evaluation
525dcad5f1cSBarry Smith - ctx   - optional context for local residual evaluation
526dcad5f1cSBarry Smith 
52720f4b53cSBarry Smith   Calling sequence of `func`:
52820f4b53cSBarry Smith $  PetscErrorCode func(DMDALocalInfo *info, void *x, void *f, void *ctx);
52920f4b53cSBarry Smith +  info - defines the subdomain to evaluate the residual on
530dcad5f1cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual
531*e4094ef1SJacob Faibussowitsch . jac - dimensional pointer to residual, write the residual here
532dcad5f1cSBarry Smith - ctx - optional context passed above
533dcad5f1cSBarry Smith 
534dcad5f1cSBarry Smith   Level: beginner
535dcad5f1cSBarry Smith 
53620f4b53cSBarry Smith   Note:
53720f4b53cSBarry Smith   The user must use `SNESSetFunction`(snes,`NULL`,`SNESPicardComputeFunction`,&user));
53820f4b53cSBarry Smith   in their code before calling this routine.
53920f4b53cSBarry Smith 
540f6dfbefdSBarry Smith .seealso: `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
541dcad5f1cSBarry Smith @*/
542d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASNESSetPicardLocal(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *, void *, void *, void *), PetscErrorCode (*jac)(DMDALocalInfo *, void *, Mat, Mat, void *), void *ctx)
543d71ae5a4SJacob Faibussowitsch {
544942e3340SBarry Smith   DMSNES     sdm;
545942e3340SBarry Smith   DMSNES_DA *dmdasnes;
546dcad5f1cSBarry Smith 
547dcad5f1cSBarry Smith   PetscFunctionBegin;
548dcad5f1cSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5499566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
5509566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
5511aa26658SKarl Rupp 
552dcad5f1cSBarry Smith   dmdasnes->residuallocalimode = imode;
553dcad5f1cSBarry Smith   dmdasnes->rhsplocal          = func;
554dcad5f1cSBarry Smith   dmdasnes->jacobianplocal     = jac;
555dcad5f1cSBarry Smith   dmdasnes->picardlocalctx     = ctx;
5561aa26658SKarl Rupp 
5579566063dSJacob Faibussowitsch   PetscCall(DMSNESSetPicard(dm, SNESComputePicard_DMDA, SNESComputePicardJacobian_DMDA, dmdasnes));
5589566063dSJacob Faibussowitsch   PetscCall(DMSNESSetMFFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
5593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
560dcad5f1cSBarry Smith }
561