xref: /petsc/src/snes/utils/dmdasnes.c (revision 420bcc1b905230dede3c88f397d1a4e60493adde)
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);
1224f572ea9SToby Isaac   PetscAssertPointer(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`:
216f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
2170038884cSBarry Smith . x    - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x)
2180038884cSBarry Smith . f    - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f)
2196cab3a1bSJed Brown - ctx  - optional context passed above
2206cab3a1bSJed Brown 
2216cab3a1bSJed Brown   Level: beginner
2226cab3a1bSJed Brown 
223*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetJacobianLocal()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2246cab3a1bSJed Brown @*/
2250b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetFunctionLocal(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *info, void *x, void *f, void *ctx), void *ctx)
226d71ae5a4SJacob Faibussowitsch {
227942e3340SBarry Smith   DMSNES     sdm;
228942e3340SBarry Smith   DMSNES_DA *dmdasnes;
2296cab3a1bSJed Brown 
2306cab3a1bSJed Brown   PetscFunctionBegin;
2316cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2329566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2339566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
2341aa26658SKarl Rupp 
235709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
2366cab3a1bSJed Brown   dmdasnes->residuallocal      = func;
2376cab3a1bSJed Brown   dmdasnes->residuallocalctx   = ctx;
2381aa26658SKarl Rupp 
2399566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
24022c6f798SBarry Smith   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
2419566063dSJacob Faibussowitsch     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
2422961e153SJed Brown   }
2433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2442961e153SJed Brown }
2452961e153SJed Brown 
2462961e153SJed Brown /*@C
247f6dfbefdSBarry Smith   DMDASNESSetFunctionLocalVec - set a local residual evaluation function that operates on a local vector for `DMDA`
2485505f7afSJunchao Zhang 
2495505f7afSJunchao Zhang   Logically Collective
2505505f7afSJunchao Zhang 
2515505f7afSJunchao Zhang   Input Parameters:
252f6dfbefdSBarry Smith + dm    - `DM` to associate callback with
253f6dfbefdSBarry Smith . imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
2545505f7afSJunchao Zhang . func  - local residual evaluation
2555505f7afSJunchao Zhang - ctx   - optional context for local residual evaluation
2565505f7afSJunchao Zhang 
25720f4b53cSBarry Smith   Calling sequence of `func`:
258f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
2595505f7afSJunchao Zhang . x    - state vector at which to evaluate residual
2605505f7afSJunchao Zhang . f    - residual vector
2615505f7afSJunchao Zhang - ctx  - optional context passed above
2625505f7afSJunchao Zhang 
2635505f7afSJunchao Zhang   Level: beginner
2645505f7afSJunchao Zhang 
265*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetFunctionLocal()`, `DMDASNESSetJacobianLocalVec()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2665505f7afSJunchao Zhang @*/
2670b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetFunctionLocalVec(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *info, Vec x, Vec f, void *ctx), void *ctx)
268d71ae5a4SJacob Faibussowitsch {
2695505f7afSJunchao Zhang   DMSNES     sdm;
2705505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
2715505f7afSJunchao Zhang 
2725505f7afSJunchao Zhang   PetscFunctionBegin;
2735505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2745505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2755505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
2765505f7afSJunchao Zhang 
2775505f7afSJunchao Zhang   dmdasnes->residuallocalimode = imode;
2785505f7afSJunchao Zhang   dmdasnes->residuallocalvec   = func;
2795505f7afSJunchao Zhang   dmdasnes->residuallocalctx   = ctx;
2805505f7afSJunchao Zhang 
2815505f7afSJunchao Zhang   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
2825505f7afSJunchao Zhang   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
2835505f7afSJunchao Zhang     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
2845505f7afSJunchao Zhang   }
2853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2865505f7afSJunchao Zhang }
2875505f7afSJunchao Zhang 
2885505f7afSJunchao Zhang /*@C
289f6dfbefdSBarry Smith   DMDASNESSetJacobianLocal - set a local Jacobian evaluation function for use with `DMDA`
2902961e153SJed Brown 
2912961e153SJed Brown   Logically Collective
2922961e153SJed Brown 
2934165533cSJose E. Roman   Input Parameters:
294f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
2956b7eb5ceSMatthew G. Knepley . func - local Jacobian evaluation
2966b7eb5ceSMatthew G. Knepley - ctx  - optional context for local Jacobian evaluation
2972961e153SJed Brown 
29820f4b53cSBarry Smith   Calling sequence of `func`:
299f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
3000038884cSBarry Smith . x    - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x)
3016b7eb5ceSMatthew G. Knepley . J    - Mat object for the Jacobian
30220f4b53cSBarry Smith . M    - Mat object for the Jacobian preconditioner matrix, often `J`
3032961e153SJed Brown - ctx  - optional context passed above
3042961e153SJed Brown 
3052961e153SJed Brown   Level: beginner
3062961e153SJed Brown 
307*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetFunctionLocal()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3082961e153SJed Brown @*/
3090b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetJacobianLocal(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, void *x, Mat J, Mat M, void *ctx), void *ctx)
310d71ae5a4SJacob Faibussowitsch {
311942e3340SBarry Smith   DMSNES     sdm;
312942e3340SBarry Smith   DMSNES_DA *dmdasnes;
3132961e153SJed Brown 
3142961e153SJed Brown   PetscFunctionBegin;
3152961e153SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3169566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3179566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
3181aa26658SKarl Rupp 
3192961e153SJed Brown   dmdasnes->jacobianlocal    = func;
3202961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
3211aa26658SKarl Rupp 
3229566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
3233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3246cab3a1bSJed Brown }
3252a4ee8f2SPeter Brune 
3262a4ee8f2SPeter Brune /*@C
327f6dfbefdSBarry Smith   DMDASNESSetJacobianLocalVec - set a local Jacobian evaluation function that operates on a local vector with `DMDA`
3285505f7afSJunchao Zhang 
3295505f7afSJunchao Zhang   Logically Collective
3305505f7afSJunchao Zhang 
3315505f7afSJunchao Zhang   Input Parameters:
332f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
3335505f7afSJunchao Zhang . func - local Jacobian evaluation
3345505f7afSJunchao Zhang - ctx  - optional context for local Jacobian evaluation
3355505f7afSJunchao Zhang 
33620f4b53cSBarry Smith   Calling sequence of `func`:
337f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
3385505f7afSJunchao Zhang . x    - state vector at which to evaluate Jacobian
3390b4db180SJacob Faibussowitsch . J    - the Jacobian
3400b4db180SJacob Faibussowitsch . M    - approximate Jacobian from which the preconditioner will be computed, often `J`
3415505f7afSJunchao Zhang - ctx  - optional context passed above
3425505f7afSJunchao Zhang 
3435505f7afSJunchao Zhang   Level: beginner
3445505f7afSJunchao Zhang 
345*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetJacobianLocal()`, `DMDASNESSetFunctionLocalVec()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3465505f7afSJunchao Zhang @*/
3470b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetJacobianLocalVec(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, Vec x, Mat J, Mat M, void *), void *ctx)
348d71ae5a4SJacob Faibussowitsch {
3495505f7afSJunchao Zhang   DMSNES     sdm;
3505505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
3515505f7afSJunchao Zhang 
3525505f7afSJunchao Zhang   PetscFunctionBegin;
3535505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3545505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3555505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
3565505f7afSJunchao Zhang 
3575505f7afSJunchao Zhang   dmdasnes->jacobianlocalvec = func;
3585505f7afSJunchao Zhang   dmdasnes->jacobianlocalctx = ctx;
3595505f7afSJunchao Zhang 
3605505f7afSJunchao Zhang   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
3613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3625505f7afSJunchao Zhang }
3635505f7afSJunchao Zhang 
3640b4db180SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
3655505f7afSJunchao Zhang /*@C
366f6dfbefdSBarry Smith   DMDASNESSetObjectiveLocal - set a local residual evaluation function to used with a `DMDA`
3672a4ee8f2SPeter Brune 
3682a4ee8f2SPeter Brune   Logically Collective
3692a4ee8f2SPeter Brune 
3704165533cSJose E. Roman   Input Parameters:
371f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
372*420bcc1bSBarry Smith . func - local objective evaluation, see `DMDASNESSetObjectiveLocal` for the calling sequence
3732a4ee8f2SPeter Brune - ctx  - optional context for local residual evaluation
3742a4ee8f2SPeter Brune 
3752a4ee8f2SPeter Brune   Level: beginner
3762a4ee8f2SPeter Brune 
377*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDASNESObjective`
3782a4ee8f2SPeter Brune @*/
379d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASNESSetObjectiveLocal(DM dm, DMDASNESObjective func, void *ctx)
380d71ae5a4SJacob Faibussowitsch {
381942e3340SBarry Smith   DMSNES     sdm;
382942e3340SBarry Smith   DMSNES_DA *dmdasnes;
3832a4ee8f2SPeter Brune 
3842a4ee8f2SPeter Brune   PetscFunctionBegin;
3852a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3869566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3879566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
3881aa26658SKarl Rupp 
3892a4ee8f2SPeter Brune   dmdasnes->objectivelocal    = func;
3902a4ee8f2SPeter Brune   dmdasnes->objectivelocalctx = ctx;
3911aa26658SKarl Rupp 
3929566063dSJacob Faibussowitsch   PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes));
3933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3942a4ee8f2SPeter Brune }
395dcad5f1cSBarry Smith 
3960b4db180SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
3975505f7afSJunchao Zhang /*@C
398e4094ef1SJacob Faibussowitsch   DMDASNESSetObjectiveLocalVec - set a local residual evaluation function that operates on a local vector with `DMDA`
3995505f7afSJunchao Zhang 
4005505f7afSJunchao Zhang   Logically Collective
4015505f7afSJunchao Zhang 
4025505f7afSJunchao Zhang   Input Parameters:
403f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
404*420bcc1bSBarry Smith . func - local objective evaluation, see `DMDASNESSetObjectiveLocalVec` for the calling sequence
4055505f7afSJunchao Zhang - ctx  - optional context for local residual evaluation
4065505f7afSJunchao Zhang 
4075505f7afSJunchao Zhang   Level: beginner
4085505f7afSJunchao Zhang 
409*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetObjectiveLocal()`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocalVec()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDASNESObjectiveVec`
4105505f7afSJunchao Zhang @*/
411d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASNESSetObjectiveLocalVec(DM dm, DMDASNESObjectiveVec func, void *ctx)
412d71ae5a4SJacob Faibussowitsch {
4135505f7afSJunchao Zhang   DMSNES     sdm;
4145505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
4155505f7afSJunchao Zhang 
4165505f7afSJunchao Zhang   PetscFunctionBegin;
4175505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4185505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
4195505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
4205505f7afSJunchao Zhang 
4215505f7afSJunchao Zhang   dmdasnes->objectivelocalvec = func;
4225505f7afSJunchao Zhang   dmdasnes->objectivelocalctx = ctx;
4235505f7afSJunchao Zhang 
4245505f7afSJunchao Zhang   PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes));
4253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4265505f7afSJunchao Zhang }
4275505f7afSJunchao Zhang 
428d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputePicard_DMDA(SNES snes, Vec X, Vec F, void *ctx)
429d71ae5a4SJacob Faibussowitsch {
430dcad5f1cSBarry Smith   DM            dm;
431942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
432dcad5f1cSBarry Smith   DMDALocalInfo info;
433dcad5f1cSBarry Smith   Vec           Xloc;
434dcad5f1cSBarry Smith   void         *x, *f;
435dcad5f1cSBarry Smith 
436dcad5f1cSBarry Smith   PetscFunctionBegin;
437dcad5f1cSBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
438dcad5f1cSBarry Smith   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
439dcad5f1cSBarry Smith   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
44028b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->rhsplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
4419566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
4429566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
4439566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
4449566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
4459566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
4469566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
447dcad5f1cSBarry Smith   switch (dmdasnes->residuallocalimode) {
448dcad5f1cSBarry Smith   case INSERT_VALUES: {
4499566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, F, &f));
450792fecdfSBarry Smith     PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx));
4519566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, F, &f));
452dcad5f1cSBarry Smith   } break;
453dcad5f1cSBarry Smith   case ADD_VALUES: {
454dcad5f1cSBarry Smith     Vec Floc;
4559566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Floc));
4569566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
4579566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Floc, &f));
458792fecdfSBarry Smith     PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx));
4599566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
4609566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
4619566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
4629566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
4639566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Floc));
464dcad5f1cSBarry Smith   } break;
465d71ae5a4SJacob Faibussowitsch   default:
466d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdasnes->residuallocalimode);
467dcad5f1cSBarry Smith   }
4689566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
4699566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
471dcad5f1cSBarry Smith }
472dcad5f1cSBarry Smith 
473d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes, Vec X, Mat A, Mat B, void *ctx)
474d71ae5a4SJacob Faibussowitsch {
475dcad5f1cSBarry Smith   DM            dm;
476942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
477dcad5f1cSBarry Smith   DMDALocalInfo info;
478dcad5f1cSBarry Smith   Vec           Xloc;
479dcad5f1cSBarry Smith   void         *x;
480dcad5f1cSBarry Smith 
481dcad5f1cSBarry Smith   PetscFunctionBegin;
48228b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->jacobianplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
4839566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
484dcad5f1cSBarry Smith 
4859566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
4869566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
4879566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
4889566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
4899566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
490792fecdfSBarry Smith   PetscCallBack("SNES Picard DMDA local callback Jacobian", (*dmdasnes->jacobianplocal)(&info, x, A, B, dmdasnes->picardlocalctx));
4919566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
4929566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
49394ab13aaSBarry Smith   if (A != B) {
4949566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4959566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
496dcad5f1cSBarry Smith   }
4973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
498dcad5f1cSBarry Smith }
499dcad5f1cSBarry Smith 
500dcad5f1cSBarry Smith /*@C
501f6dfbefdSBarry Smith   DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration with `DMDA`
502dcad5f1cSBarry Smith 
503dcad5f1cSBarry Smith   Logically Collective
504dcad5f1cSBarry Smith 
5054165533cSJose E. Roman   Input Parameters:
506f6dfbefdSBarry Smith + dm    - `DM` to associate callback with
507f6dfbefdSBarry Smith . imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
508dcad5f1cSBarry Smith . func  - local residual evaluation
509ceaaa498SBarry Smith . jac   - function to compute Jacobian
510dcad5f1cSBarry Smith - ctx   - optional context for local residual evaluation
511dcad5f1cSBarry Smith 
51220f4b53cSBarry Smith   Calling sequence of `func`:
51320f4b53cSBarry Smith + info - defines the subdomain to evaluate the residual on
514dcad5f1cSBarry Smith . x    - dimensional pointer to state at which to evaluate residual
5150b4db180SJacob Faibussowitsch . f    - dimensional pointer to residual, write the residual here
5160b4db180SJacob Faibussowitsch - ctx  - optional context passed above
5170b4db180SJacob Faibussowitsch 
5180b4db180SJacob Faibussowitsch   Calling sequence of `jac`:
5190b4db180SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on
5200b4db180SJacob Faibussowitsch . x    - dimensional pointer to state at which to evaluate residual
5210b4db180SJacob Faibussowitsch . jac  - the Jacobian
5220b4db180SJacob Faibussowitsch . Jp   - approximation to the Jacobian used to compute the preconditioner, often `J`
523dcad5f1cSBarry Smith - ctx  - optional context passed above
524dcad5f1cSBarry Smith 
525dcad5f1cSBarry Smith   Level: beginner
526dcad5f1cSBarry Smith 
52720f4b53cSBarry Smith   Note:
528*420bcc1bSBarry Smith   The user must use `SNESSetFunction`(`snes`,`NULL`,`SNESPicardComputeFunction`,&user));
52920f4b53cSBarry Smith   in their code before calling this routine.
53020f4b53cSBarry Smith 
531*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
532dcad5f1cSBarry Smith @*/
5330b4db180SJacob 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)
534d71ae5a4SJacob Faibussowitsch {
535942e3340SBarry Smith   DMSNES     sdm;
536942e3340SBarry Smith   DMSNES_DA *dmdasnes;
537dcad5f1cSBarry Smith 
538dcad5f1cSBarry Smith   PetscFunctionBegin;
539dcad5f1cSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5409566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
5419566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
5421aa26658SKarl Rupp 
543dcad5f1cSBarry Smith   dmdasnes->residuallocalimode = imode;
544dcad5f1cSBarry Smith   dmdasnes->rhsplocal          = func;
545dcad5f1cSBarry Smith   dmdasnes->jacobianplocal     = jac;
546dcad5f1cSBarry Smith   dmdasnes->picardlocalctx     = ctx;
5471aa26658SKarl Rupp 
5489566063dSJacob Faibussowitsch   PetscCall(DMSNESSetPicard(dm, SNESComputePicard_DMDA, SNESComputePicardJacobian_DMDA, dmdasnes));
5499566063dSJacob Faibussowitsch   PetscCall(DMSNESSetMFFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
5503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
551dcad5f1cSBarry Smith }
552