xref: /petsc/src/snes/utils/dmlocalsnes.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
3ff35dfedSBarry Smith 
4ff35dfedSBarry Smith typedef struct {
5ff35dfedSBarry Smith   PetscErrorCode (*residuallocal)(DM, Vec, Vec, void *);
6d1e9a80fSBarry Smith   PetscErrorCode (*jacobianlocal)(DM, Vec, Mat, Mat, void *);
7bdd6f66aSToby Isaac   PetscErrorCode (*boundarylocal)(DM, Vec, void *);
8ff35dfedSBarry Smith   void *residuallocalctx;
9ff35dfedSBarry Smith   void *jacobianlocalctx;
10bdd6f66aSToby Isaac   void *boundarylocalctx;
11942e3340SBarry Smith } DMSNES_Local;
12ff35dfedSBarry Smith 
13*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDestroy_DMLocal(DMSNES sdm)
14*d71ae5a4SJacob Faibussowitsch {
15ff35dfedSBarry Smith   PetscFunctionBegin;
169566063dSJacob Faibussowitsch   PetscCall(PetscFree(sdm->data));
175ed5e208SMatthew G. Knepley   sdm->data = NULL;
18ff35dfedSBarry Smith   PetscFunctionReturn(0);
19ff35dfedSBarry Smith }
20ff35dfedSBarry Smith 
21*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm, DMSNES sdm)
22*d71ae5a4SJacob Faibussowitsch {
23ff35dfedSBarry Smith   PetscFunctionBegin;
245ed5e208SMatthew G. Knepley   if (sdm->data != oldsdm->data) {
259566063dSJacob Faibussowitsch     PetscCall(PetscFree(sdm->data));
264dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMSNES_Local **)&sdm->data));
279566063dSJacob Faibussowitsch     if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data, oldsdm->data, sizeof(DMSNES_Local)));
28ff35dfedSBarry Smith   }
29ff35dfedSBarry Smith   PetscFunctionReturn(0);
30ff35dfedSBarry Smith }
31ff35dfedSBarry Smith 
32*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMLocalSNESGetContext(DM dm, DMSNES sdm, DMSNES_Local **dmlocalsnes)
33*d71ae5a4SJacob Faibussowitsch {
34ff35dfedSBarry Smith   PetscFunctionBegin;
350298fd71SBarry Smith   *dmlocalsnes = NULL;
36ff35dfedSBarry Smith   if (!sdm->data) {
374dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMSNES_Local **)&sdm->data));
381aa26658SKarl Rupp 
3922c6f798SBarry Smith     sdm->ops->destroy   = DMSNESDestroy_DMLocal;
4022c6f798SBarry Smith     sdm->ops->duplicate = DMSNESDuplicate_DMLocal;
41ff35dfedSBarry Smith   }
42942e3340SBarry Smith   *dmlocalsnes = (DMSNES_Local *)sdm->data;
43ff35dfedSBarry Smith   PetscFunctionReturn(0);
44ff35dfedSBarry Smith }
45ff35dfedSBarry Smith 
46*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes, Vec X, Vec F, void *ctx)
47*d71ae5a4SJacob Faibussowitsch {
48942e3340SBarry Smith   DMSNES_Local *dmlocalsnes = (DMSNES_Local *)ctx;
497a1e7fa1SMatthew G. Knepley   DM            dm;
50ff35dfedSBarry Smith   Vec           Xloc, Floc;
517a1e7fa1SMatthew G. Knepley   PetscBool     transform;
52ff35dfedSBarry Smith 
53ff35dfedSBarry Smith   PetscFunctionBegin;
54ff35dfedSBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
55ff35dfedSBarry Smith   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
56ff35dfedSBarry Smith   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
579566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
589566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
599566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Floc));
609566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Xloc));
619566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Floc));
627a1e7fa1SMatthew G. Knepley   /* Non-conforming routines needs boundary values before G2L */
639566063dSJacob Faibussowitsch   if (dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
649566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
659566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
667a1e7fa1SMatthew G. Knepley   /* Need to reset boundary values if we transformed */
679566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dm, &transform));
689566063dSJacob Faibussowitsch   if (transform && dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
69ff35dfedSBarry Smith   CHKMEMQ;
709566063dSJacob Faibussowitsch   PetscCall((*dmlocalsnes->residuallocal)(dm, Xloc, Floc, dmlocalsnes->residuallocalctx));
71ff35dfedSBarry Smith   CHKMEMQ;
729566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
739566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
749566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
759566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Floc));
769566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
77afa9dec6SMatthew G. Knepley   {
78afa9dec6SMatthew G. Knepley     char        name[PETSC_MAX_PATH_LEN];
79f2f0fc17SMatthew G. Knepley     char        oldname[PETSC_MAX_PATH_LEN];
80f2f0fc17SMatthew G. Knepley     const char *tmp;
81afa9dec6SMatthew G. Knepley     PetscInt    it;
82afa9dec6SMatthew G. Knepley 
839566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(snes, &it));
849566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Solution, Iterate %d", (int)it));
859566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)X, &tmp));
869566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(oldname, tmp, PETSC_MAX_PATH_LEN - 1));
879566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)X, name));
889566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(X, (PetscObject)snes, "-dmsnes_solution_vec_view"));
899566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)X, oldname));
909566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Residual, Iterate %d", (int)it));
919566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)F, name));
929566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(F, (PetscObject)snes, "-dmsnes_residual_vec_view"));
93afa9dec6SMatthew G. Knepley   }
94ff35dfedSBarry Smith   PetscFunctionReturn(0);
95ff35dfedSBarry Smith }
96ff35dfedSBarry Smith 
97*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes, Vec X, Mat A, Mat B, void *ctx)
98*d71ae5a4SJacob Faibussowitsch {
99942e3340SBarry Smith   DMSNES_Local *dmlocalsnes = (DMSNES_Local *)ctx;
1007a1e7fa1SMatthew G. Knepley   DM            dm;
101ff35dfedSBarry Smith   Vec           Xloc;
1027a1e7fa1SMatthew G. Knepley   PetscBool     transform;
103ff35dfedSBarry Smith 
104ff35dfedSBarry Smith   PetscFunctionBegin;
1059566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
106ff35dfedSBarry Smith   if (dmlocalsnes->jacobianlocal) {
1079566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Xloc));
1089566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Xloc));
1097a1e7fa1SMatthew G. Knepley     /* Non-conforming routines needs boundary values before G2L */
1109566063dSJacob Faibussowitsch     if (dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
1119566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1129566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1137a1e7fa1SMatthew G. Knepley     /* Need to reset boundary values if we transformed */
1149566063dSJacob Faibussowitsch     PetscCall(DMHasBasisTransform(dm, &transform));
1159566063dSJacob Faibussowitsch     if (transform && dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
116ff35dfedSBarry Smith     CHKMEMQ;
1179566063dSJacob Faibussowitsch     PetscCall((*dmlocalsnes->jacobianlocal)(dm, Xloc, A, B, dmlocalsnes->jacobianlocalctx));
118ff35dfedSBarry Smith     CHKMEMQ;
1199566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Xloc));
120ff35dfedSBarry Smith   } else {
121ff35dfedSBarry Smith     MatFDColoring fdcoloring;
1229566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring));
123ff35dfedSBarry Smith     if (!fdcoloring) {
124ff35dfedSBarry Smith       ISColoring coloring;
125ff35dfedSBarry Smith 
1269566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring));
1279566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring));
1289566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&coloring));
129ff35dfedSBarry Smith       switch (dm->coloringtype) {
130*d71ae5a4SJacob Faibussowitsch       case IS_COLORING_GLOBAL:
131*d71ae5a4SJacob Faibussowitsch         PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))SNESComputeFunction_DMLocal, dmlocalsnes));
132*d71ae5a4SJacob Faibussowitsch         break;
133*d71ae5a4SJacob Faibussowitsch       default:
134*d71ae5a4SJacob Faibussowitsch         SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]);
135ff35dfedSBarry Smith       }
1369566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix));
1379566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1389566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring));
1399566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring));
1409566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)fdcoloring));
141ff35dfedSBarry Smith 
142ff35dfedSBarry Smith       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
143ff35dfedSBarry Smith        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
144ff35dfedSBarry Smith        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
145ff35dfedSBarry Smith        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
146140e18c1SBarry Smith        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
147ff35dfedSBarry Smith        */
1489566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dm));
149ff35dfedSBarry Smith     }
1509566063dSJacob Faibussowitsch     PetscCall(MatFDColoringApply(B, fdcoloring, X, snes));
151ff35dfedSBarry Smith   }
152ff35dfedSBarry Smith   /* This will be redundant if the user called both, but it's too common to forget. */
15394ab13aaSBarry Smith   if (A != B) {
1549566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1559566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
156ff35dfedSBarry Smith   }
157ff35dfedSBarry Smith   PetscFunctionReturn(0);
158ff35dfedSBarry Smith }
159ff35dfedSBarry Smith 
160ff35dfedSBarry Smith /*@C
161ff35dfedSBarry Smith    DMSNESSetFunctionLocal - set a local residual evaluation function. This function is called with local vector
162ff35dfedSBarry Smith       containing the local vector information PLUS ghost point information. It should compute a result for all local
163f6dfbefdSBarry Smith       elements and `DMSNES` will automatically accumulate the overlapping values.
164ff35dfedSBarry Smith 
165ff35dfedSBarry Smith    Logically Collective
166ff35dfedSBarry Smith 
1674165533cSJose E. Roman    Input Parameters:
168f6dfbefdSBarry Smith +  dm - `DM` to associate callback with
169ff35dfedSBarry Smith .  func - local residual evaluation
170ff35dfedSBarry Smith -  ctx - optional context for local residual evaluation
171ff35dfedSBarry Smith 
172f6dfbefdSBarry Smith    Level: advanced
173ff35dfedSBarry Smith 
174db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
175ff35dfedSBarry Smith @*/
176*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESSetFunctionLocal(DM dm, PetscErrorCode (*func)(DM, Vec, Vec, void *), void *ctx)
177*d71ae5a4SJacob Faibussowitsch {
178942e3340SBarry Smith   DMSNES        sdm;
179942e3340SBarry Smith   DMSNES_Local *dmlocalsnes;
180ff35dfedSBarry Smith 
181ff35dfedSBarry Smith   PetscFunctionBegin;
182ff35dfedSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1839566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
1849566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
1851aa26658SKarl Rupp 
186ff35dfedSBarry Smith   dmlocalsnes->residuallocal    = func;
187ff35dfedSBarry Smith   dmlocalsnes->residuallocalctx = ctx;
1881aa26658SKarl Rupp 
1899566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMLocal, dmlocalsnes));
19022c6f798SBarry Smith   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
1919566063dSJacob Faibussowitsch     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMLocal, dmlocalsnes));
192ff35dfedSBarry Smith   }
193ff35dfedSBarry Smith   PetscFunctionReturn(0);
194ff35dfedSBarry Smith }
195ff35dfedSBarry Smith 
196bdd6f66aSToby Isaac /*@C
197bdd6f66aSToby Isaac    DMSNESSetBoundaryLocal - set a local boundary value function. This function is called with local vector
198bdd6f66aSToby Isaac       containing the local vector information PLUS ghost point information. It should insert values into the local
199bdd6f66aSToby Isaac       vector that do not come from the global vector, such as essential boundary condition data.
200bdd6f66aSToby Isaac 
201bdd6f66aSToby Isaac    Logically Collective
202bdd6f66aSToby Isaac 
2034165533cSJose E. Roman    Input Parameters:
204f6dfbefdSBarry Smith +  dm - `DM` to associate callback with
205bdd6f66aSToby Isaac .  func - local boundary value evaluation
206bdd6f66aSToby Isaac -  ctx - optional context for local boundary value evaluation
207bdd6f66aSToby Isaac 
208f6dfbefdSBarry Smith    Level: advanced
209bdd6f66aSToby Isaac 
210db781477SPatrick Sanan .seealso: `DMSNESSetFunctionLocal()`, `DMDASNESSetJacobianLocal()`
211bdd6f66aSToby Isaac @*/
212*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESSetBoundaryLocal(DM dm, PetscErrorCode (*func)(DM, Vec, void *), void *ctx)
213*d71ae5a4SJacob Faibussowitsch {
214bdd6f66aSToby Isaac   DMSNES        sdm;
215bdd6f66aSToby Isaac   DMSNES_Local *dmlocalsnes;
216bdd6f66aSToby Isaac 
217bdd6f66aSToby Isaac   PetscFunctionBegin;
218bdd6f66aSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2199566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2209566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
221bdd6f66aSToby Isaac 
222bdd6f66aSToby Isaac   dmlocalsnes->boundarylocal    = func;
223bdd6f66aSToby Isaac   dmlocalsnes->boundarylocalctx = ctx;
224bdd6f66aSToby Isaac 
225bdd6f66aSToby Isaac   PetscFunctionReturn(0);
226bdd6f66aSToby Isaac }
227bdd6f66aSToby Isaac 
228ff35dfedSBarry Smith /*@C
229ff35dfedSBarry Smith    DMSNESSetJacobianLocal - set a local Jacobian evaluation function
230ff35dfedSBarry Smith 
231ff35dfedSBarry Smith    Logically Collective
232ff35dfedSBarry Smith 
2334165533cSJose E. Roman    Input Parameters:
234ff35dfedSBarry Smith +  dm - DM to associate callback with
235ff35dfedSBarry Smith .  func - local Jacobian evaluation
236ff35dfedSBarry Smith -  ctx - optional context for local Jacobian evaluation
237ff35dfedSBarry Smith 
238f6dfbefdSBarry Smith    Level: advanced
239ff35dfedSBarry Smith 
240db781477SPatrick Sanan .seealso: `DMSNESSetJacobian()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
241ff35dfedSBarry Smith @*/
242*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESSetJacobianLocal(DM dm, PetscErrorCode (*func)(DM, Vec, Mat, Mat, void *), void *ctx)
243*d71ae5a4SJacob Faibussowitsch {
244942e3340SBarry Smith   DMSNES        sdm;
245942e3340SBarry Smith   DMSNES_Local *dmlocalsnes;
246ff35dfedSBarry Smith 
247ff35dfedSBarry Smith   PetscFunctionBegin;
248ff35dfedSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2499566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2509566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
2511aa26658SKarl Rupp 
252ff35dfedSBarry Smith   dmlocalsnes->jacobianlocal    = func;
253ff35dfedSBarry Smith   dmlocalsnes->jacobianlocalctx = ctx;
2541aa26658SKarl Rupp 
2559566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMLocal, dmlocalsnes));
25628d58a37SPierre Jolivet   PetscFunctionReturn(0);
25728d58a37SPierre Jolivet }
25828d58a37SPierre Jolivet 
25928d58a37SPierre Jolivet /*@C
260f6dfbefdSBarry Smith    DMSNESGetFunctionLocal - get the local residual evaluation function information set with `DMSNESSetFunctionLocal()`.
26128d58a37SPierre Jolivet 
26228d58a37SPierre Jolivet    Not Collective
26328d58a37SPierre Jolivet 
2644165533cSJose E. Roman    Input Parameter:
265f6dfbefdSBarry Smith .  dm - `DM` with the associated callback
26628d58a37SPierre Jolivet 
2674165533cSJose E. Roman    Output Parameters:
26828d58a37SPierre Jolivet +  func - local residual evaluation
26928d58a37SPierre Jolivet -  ctx - context for local residual evaluation
27028d58a37SPierre Jolivet 
27128d58a37SPierre Jolivet    Level: beginner
27228d58a37SPierre Jolivet 
273db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMSNESSetFunctionLocal()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
27428d58a37SPierre Jolivet @*/
275*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESGetFunctionLocal(DM dm, PetscErrorCode (**func)(DM, Vec, Vec, void *), void **ctx)
276*d71ae5a4SJacob Faibussowitsch {
27728d58a37SPierre Jolivet   DMSNES        sdm;
27828d58a37SPierre Jolivet   DMSNES_Local *dmlocalsnes;
27928d58a37SPierre Jolivet 
28028d58a37SPierre Jolivet   PetscFunctionBegin;
28128d58a37SPierre Jolivet   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2829566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
2839566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
28428d58a37SPierre Jolivet   if (func) *func = dmlocalsnes->residuallocal;
28528d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->residuallocalctx;
28628d58a37SPierre Jolivet   PetscFunctionReturn(0);
28728d58a37SPierre Jolivet }
28828d58a37SPierre Jolivet 
28928d58a37SPierre Jolivet /*@C
290f6dfbefdSBarry Smith    DMSNESGetBoundaryLocal - get the local boundary value function set with `DMSNESSetBoundaryLocal()`.
29128d58a37SPierre Jolivet 
29228d58a37SPierre Jolivet    Not Collective
29328d58a37SPierre Jolivet 
2944165533cSJose E. Roman    Input Parameter:
295f6dfbefdSBarry Smith .  dm - `DM` with the associated callback
29628d58a37SPierre Jolivet 
2974165533cSJose E. Roman    Output Parameters:
29828d58a37SPierre Jolivet +  func - local boundary value evaluation
29928d58a37SPierre Jolivet -  ctx - context for local boundary value evaluation
30028d58a37SPierre Jolivet 
30128d58a37SPierre Jolivet    Level: intermediate
30228d58a37SPierre Jolivet 
303db781477SPatrick Sanan .seealso: `DMSNESSetFunctionLocal()`, `DMSNESSetBoundaryLocal()`, `DMDASNESSetJacobianLocal()`
30428d58a37SPierre Jolivet @*/
305*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESGetBoundaryLocal(DM dm, PetscErrorCode (**func)(DM, Vec, void *), void **ctx)
306*d71ae5a4SJacob Faibussowitsch {
30728d58a37SPierre Jolivet   DMSNES        sdm;
30828d58a37SPierre Jolivet   DMSNES_Local *dmlocalsnes;
30928d58a37SPierre Jolivet 
31028d58a37SPierre Jolivet   PetscFunctionBegin;
31128d58a37SPierre Jolivet   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3129566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
3139566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
31428d58a37SPierre Jolivet   if (func) *func = dmlocalsnes->boundarylocal;
31528d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->boundarylocalctx;
31628d58a37SPierre Jolivet   PetscFunctionReturn(0);
31728d58a37SPierre Jolivet }
31828d58a37SPierre Jolivet 
31928d58a37SPierre Jolivet /*@C
320f6dfbefdSBarry Smith    DMSNESGetJacobianLocal - the local Jacobian evaluation function set with `DMSNESSetJacobianLocal()`.
32128d58a37SPierre Jolivet 
32228d58a37SPierre Jolivet    Logically Collective
32328d58a37SPierre Jolivet 
3244165533cSJose E. Roman    Input Parameter:
325f6dfbefdSBarry Smith .  dm - `DM` with the associated callback
32628d58a37SPierre Jolivet 
3274165533cSJose E. Roman    Output Parameters:
32828d58a37SPierre Jolivet +  func - local Jacobian evaluation
32928d58a37SPierre Jolivet -  ctx - context for local Jacobian evaluation
33028d58a37SPierre Jolivet 
33128d58a37SPierre Jolivet    Level: beginner
33228d58a37SPierre Jolivet 
333db781477SPatrick Sanan .seealso: `DMSNESSetJacobianLocal()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
33428d58a37SPierre Jolivet @*/
335*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESGetJacobianLocal(DM dm, PetscErrorCode (**func)(DM, Vec, Mat, Mat, void *), void **ctx)
336*d71ae5a4SJacob Faibussowitsch {
33728d58a37SPierre Jolivet   DMSNES        sdm;
33828d58a37SPierre Jolivet   DMSNES_Local *dmlocalsnes;
33928d58a37SPierre Jolivet 
34028d58a37SPierre Jolivet   PetscFunctionBegin;
34128d58a37SPierre Jolivet   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3429566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
3439566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
34428d58a37SPierre Jolivet   if (func) *func = dmlocalsnes->jacobianlocal;
34528d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->jacobianlocalctx;
346ff35dfedSBarry Smith   PetscFunctionReturn(0);
347ff35dfedSBarry Smith }
348