xref: /petsc/src/snes/utils/dmlocalsnes.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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 
139371c9d4SSatish Balay static PetscErrorCode DMSNESDestroy_DMLocal(DMSNES sdm) {
14ff35dfedSBarry Smith   PetscFunctionBegin;
159566063dSJacob Faibussowitsch   PetscCall(PetscFree(sdm->data));
165ed5e208SMatthew G. Knepley   sdm->data = NULL;
17ff35dfedSBarry Smith   PetscFunctionReturn(0);
18ff35dfedSBarry Smith }
19ff35dfedSBarry Smith 
209371c9d4SSatish Balay static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm, DMSNES sdm) {
21ff35dfedSBarry Smith   PetscFunctionBegin;
225ed5e208SMatthew G. Knepley   if (sdm->data != oldsdm->data) {
239566063dSJacob Faibussowitsch     PetscCall(PetscFree(sdm->data));
24*4dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMSNES_Local **)&sdm->data));
259566063dSJacob Faibussowitsch     if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data, oldsdm->data, sizeof(DMSNES_Local)));
26ff35dfedSBarry Smith   }
27ff35dfedSBarry Smith   PetscFunctionReturn(0);
28ff35dfedSBarry Smith }
29ff35dfedSBarry Smith 
309371c9d4SSatish Balay static PetscErrorCode DMLocalSNESGetContext(DM dm, DMSNES sdm, DMSNES_Local **dmlocalsnes) {
31ff35dfedSBarry Smith   PetscFunctionBegin;
320298fd71SBarry Smith   *dmlocalsnes = NULL;
33ff35dfedSBarry Smith   if (!sdm->data) {
34*4dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMSNES_Local **)&sdm->data));
351aa26658SKarl Rupp 
3622c6f798SBarry Smith     sdm->ops->destroy   = DMSNESDestroy_DMLocal;
3722c6f798SBarry Smith     sdm->ops->duplicate = DMSNESDuplicate_DMLocal;
38ff35dfedSBarry Smith   }
39942e3340SBarry Smith   *dmlocalsnes = (DMSNES_Local *)sdm->data;
40ff35dfedSBarry Smith   PetscFunctionReturn(0);
41ff35dfedSBarry Smith }
42ff35dfedSBarry Smith 
439371c9d4SSatish Balay static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes, Vec X, Vec F, void *ctx) {
44942e3340SBarry Smith   DMSNES_Local *dmlocalsnes = (DMSNES_Local *)ctx;
457a1e7fa1SMatthew G. Knepley   DM            dm;
46ff35dfedSBarry Smith   Vec           Xloc, Floc;
477a1e7fa1SMatthew G. Knepley   PetscBool     transform;
48ff35dfedSBarry Smith 
49ff35dfedSBarry Smith   PetscFunctionBegin;
50ff35dfedSBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
51ff35dfedSBarry Smith   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
52ff35dfedSBarry Smith   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
539566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
549566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
559566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Floc));
569566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Xloc));
579566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Floc));
587a1e7fa1SMatthew G. Knepley   /* Non-conforming routines needs boundary values before G2L */
599566063dSJacob Faibussowitsch   if (dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
609566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
619566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
627a1e7fa1SMatthew G. Knepley   /* Need to reset boundary values if we transformed */
639566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dm, &transform));
649566063dSJacob Faibussowitsch   if (transform && dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
65ff35dfedSBarry Smith   CHKMEMQ;
669566063dSJacob Faibussowitsch   PetscCall((*dmlocalsnes->residuallocal)(dm, Xloc, Floc, dmlocalsnes->residuallocalctx));
67ff35dfedSBarry Smith   CHKMEMQ;
689566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
699566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
709566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
719566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Floc));
729566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
73afa9dec6SMatthew G. Knepley   {
74afa9dec6SMatthew G. Knepley     char        name[PETSC_MAX_PATH_LEN];
75f2f0fc17SMatthew G. Knepley     char        oldname[PETSC_MAX_PATH_LEN];
76f2f0fc17SMatthew G. Knepley     const char *tmp;
77afa9dec6SMatthew G. Knepley     PetscInt    it;
78afa9dec6SMatthew G. Knepley 
799566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(snes, &it));
809566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Solution, Iterate %d", (int)it));
819566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)X, &tmp));
829566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(oldname, tmp, PETSC_MAX_PATH_LEN - 1));
839566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)X, name));
849566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(X, (PetscObject)snes, "-dmsnes_solution_vec_view"));
859566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)X, oldname));
869566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Residual, Iterate %d", (int)it));
879566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)F, name));
889566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(F, (PetscObject)snes, "-dmsnes_residual_vec_view"));
89afa9dec6SMatthew G. Knepley   }
90ff35dfedSBarry Smith   PetscFunctionReturn(0);
91ff35dfedSBarry Smith }
92ff35dfedSBarry Smith 
939371c9d4SSatish Balay static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes, Vec X, Mat A, Mat B, void *ctx) {
94942e3340SBarry Smith   DMSNES_Local *dmlocalsnes = (DMSNES_Local *)ctx;
957a1e7fa1SMatthew G. Knepley   DM            dm;
96ff35dfedSBarry Smith   Vec           Xloc;
977a1e7fa1SMatthew G. Knepley   PetscBool     transform;
98ff35dfedSBarry Smith 
99ff35dfedSBarry Smith   PetscFunctionBegin;
1009566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
101ff35dfedSBarry Smith   if (dmlocalsnes->jacobianlocal) {
1029566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Xloc));
1039566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Xloc));
1047a1e7fa1SMatthew G. Knepley     /* Non-conforming routines needs boundary values before G2L */
1059566063dSJacob Faibussowitsch     if (dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
1069566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1079566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1087a1e7fa1SMatthew G. Knepley     /* Need to reset boundary values if we transformed */
1099566063dSJacob Faibussowitsch     PetscCall(DMHasBasisTransform(dm, &transform));
1109566063dSJacob Faibussowitsch     if (transform && dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
111ff35dfedSBarry Smith     CHKMEMQ;
1129566063dSJacob Faibussowitsch     PetscCall((*dmlocalsnes->jacobianlocal)(dm, Xloc, A, B, dmlocalsnes->jacobianlocalctx));
113ff35dfedSBarry Smith     CHKMEMQ;
1149566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Xloc));
115ff35dfedSBarry Smith   } else {
116ff35dfedSBarry Smith     MatFDColoring fdcoloring;
1179566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring));
118ff35dfedSBarry Smith     if (!fdcoloring) {
119ff35dfedSBarry Smith       ISColoring coloring;
120ff35dfedSBarry Smith 
1219566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring));
1229566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring));
1239566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&coloring));
124ff35dfedSBarry Smith       switch (dm->coloringtype) {
1259371c9d4SSatish Balay       case IS_COLORING_GLOBAL: PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))SNESComputeFunction_DMLocal, dmlocalsnes)); break;
12698921bdaSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]);
127ff35dfedSBarry Smith       }
1289566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix));
1299566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1309566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring));
1319566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring));
1329566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)fdcoloring));
133ff35dfedSBarry Smith 
134ff35dfedSBarry Smith       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
135ff35dfedSBarry Smith        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
136ff35dfedSBarry Smith        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
137ff35dfedSBarry Smith        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
138140e18c1SBarry Smith        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
139ff35dfedSBarry Smith        */
1409566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dm));
141ff35dfedSBarry Smith     }
1429566063dSJacob Faibussowitsch     PetscCall(MatFDColoringApply(B, fdcoloring, X, snes));
143ff35dfedSBarry Smith   }
144ff35dfedSBarry Smith   /* This will be redundant if the user called both, but it's too common to forget. */
14594ab13aaSBarry Smith   if (A != B) {
1469566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1479566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
148ff35dfedSBarry Smith   }
149ff35dfedSBarry Smith   PetscFunctionReturn(0);
150ff35dfedSBarry Smith }
151ff35dfedSBarry Smith 
152ff35dfedSBarry Smith /*@C
153ff35dfedSBarry Smith    DMSNESSetFunctionLocal - set a local residual evaluation function. This function is called with local vector
154ff35dfedSBarry Smith       containing the local vector information PLUS ghost point information. It should compute a result for all local
155f6dfbefdSBarry Smith       elements and `DMSNES` will automatically accumulate the overlapping values.
156ff35dfedSBarry Smith 
157ff35dfedSBarry Smith    Logically Collective
158ff35dfedSBarry Smith 
1594165533cSJose E. Roman    Input Parameters:
160f6dfbefdSBarry Smith +  dm - `DM` to associate callback with
161ff35dfedSBarry Smith .  func - local residual evaluation
162ff35dfedSBarry Smith -  ctx - optional context for local residual evaluation
163ff35dfedSBarry Smith 
164f6dfbefdSBarry Smith    Level: advanced
165ff35dfedSBarry Smith 
166db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
167ff35dfedSBarry Smith @*/
1689371c9d4SSatish Balay PetscErrorCode DMSNESSetFunctionLocal(DM dm, PetscErrorCode (*func)(DM, Vec, Vec, void *), void *ctx) {
169942e3340SBarry Smith   DMSNES        sdm;
170942e3340SBarry Smith   DMSNES_Local *dmlocalsnes;
171ff35dfedSBarry Smith 
172ff35dfedSBarry Smith   PetscFunctionBegin;
173ff35dfedSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1749566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
1759566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
1761aa26658SKarl Rupp 
177ff35dfedSBarry Smith   dmlocalsnes->residuallocal    = func;
178ff35dfedSBarry Smith   dmlocalsnes->residuallocalctx = ctx;
1791aa26658SKarl Rupp 
1809566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMLocal, dmlocalsnes));
18122c6f798SBarry Smith   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
1829566063dSJacob Faibussowitsch     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMLocal, dmlocalsnes));
183ff35dfedSBarry Smith   }
184ff35dfedSBarry Smith   PetscFunctionReturn(0);
185ff35dfedSBarry Smith }
186ff35dfedSBarry Smith 
187bdd6f66aSToby Isaac /*@C
188bdd6f66aSToby Isaac    DMSNESSetBoundaryLocal - set a local boundary value function. This function is called with local vector
189bdd6f66aSToby Isaac       containing the local vector information PLUS ghost point information. It should insert values into the local
190bdd6f66aSToby Isaac       vector that do not come from the global vector, such as essential boundary condition data.
191bdd6f66aSToby Isaac 
192bdd6f66aSToby Isaac    Logically Collective
193bdd6f66aSToby Isaac 
1944165533cSJose E. Roman    Input Parameters:
195f6dfbefdSBarry Smith +  dm - `DM` to associate callback with
196bdd6f66aSToby Isaac .  func - local boundary value evaluation
197bdd6f66aSToby Isaac -  ctx - optional context for local boundary value evaluation
198bdd6f66aSToby Isaac 
199f6dfbefdSBarry Smith    Level: advanced
200bdd6f66aSToby Isaac 
201db781477SPatrick Sanan .seealso: `DMSNESSetFunctionLocal()`, `DMDASNESSetJacobianLocal()`
202bdd6f66aSToby Isaac @*/
2039371c9d4SSatish Balay PetscErrorCode DMSNESSetBoundaryLocal(DM dm, PetscErrorCode (*func)(DM, Vec, void *), void *ctx) {
204bdd6f66aSToby Isaac   DMSNES        sdm;
205bdd6f66aSToby Isaac   DMSNES_Local *dmlocalsnes;
206bdd6f66aSToby Isaac 
207bdd6f66aSToby Isaac   PetscFunctionBegin;
208bdd6f66aSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2099566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2109566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
211bdd6f66aSToby Isaac 
212bdd6f66aSToby Isaac   dmlocalsnes->boundarylocal    = func;
213bdd6f66aSToby Isaac   dmlocalsnes->boundarylocalctx = ctx;
214bdd6f66aSToby Isaac 
215bdd6f66aSToby Isaac   PetscFunctionReturn(0);
216bdd6f66aSToby Isaac }
217bdd6f66aSToby Isaac 
218ff35dfedSBarry Smith /*@C
219ff35dfedSBarry Smith    DMSNESSetJacobianLocal - set a local Jacobian evaluation function
220ff35dfedSBarry Smith 
221ff35dfedSBarry Smith    Logically Collective
222ff35dfedSBarry Smith 
2234165533cSJose E. Roman    Input Parameters:
224ff35dfedSBarry Smith +  dm - DM to associate callback with
225ff35dfedSBarry Smith .  func - local Jacobian evaluation
226ff35dfedSBarry Smith -  ctx - optional context for local Jacobian evaluation
227ff35dfedSBarry Smith 
228f6dfbefdSBarry Smith    Level: advanced
229ff35dfedSBarry Smith 
230db781477SPatrick Sanan .seealso: `DMSNESSetJacobian()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
231ff35dfedSBarry Smith @*/
2329371c9d4SSatish Balay PetscErrorCode DMSNESSetJacobianLocal(DM dm, PetscErrorCode (*func)(DM, Vec, Mat, Mat, void *), void *ctx) {
233942e3340SBarry Smith   DMSNES        sdm;
234942e3340SBarry Smith   DMSNES_Local *dmlocalsnes;
235ff35dfedSBarry Smith 
236ff35dfedSBarry Smith   PetscFunctionBegin;
237ff35dfedSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2389566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2399566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
2401aa26658SKarl Rupp 
241ff35dfedSBarry Smith   dmlocalsnes->jacobianlocal    = func;
242ff35dfedSBarry Smith   dmlocalsnes->jacobianlocalctx = ctx;
2431aa26658SKarl Rupp 
2449566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMLocal, dmlocalsnes));
24528d58a37SPierre Jolivet   PetscFunctionReturn(0);
24628d58a37SPierre Jolivet }
24728d58a37SPierre Jolivet 
24828d58a37SPierre Jolivet /*@C
249f6dfbefdSBarry Smith    DMSNESGetFunctionLocal - get the local residual evaluation function information set with `DMSNESSetFunctionLocal()`.
25028d58a37SPierre Jolivet 
25128d58a37SPierre Jolivet    Not Collective
25228d58a37SPierre Jolivet 
2534165533cSJose E. Roman    Input Parameter:
254f6dfbefdSBarry Smith .  dm - `DM` with the associated callback
25528d58a37SPierre Jolivet 
2564165533cSJose E. Roman    Output Parameters:
25728d58a37SPierre Jolivet +  func - local residual evaluation
25828d58a37SPierre Jolivet -  ctx - context for local residual evaluation
25928d58a37SPierre Jolivet 
26028d58a37SPierre Jolivet    Level: beginner
26128d58a37SPierre Jolivet 
262db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMSNESSetFunctionLocal()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
26328d58a37SPierre Jolivet @*/
2649371c9d4SSatish Balay PetscErrorCode DMSNESGetFunctionLocal(DM dm, PetscErrorCode (**func)(DM, Vec, Vec, void *), void **ctx) {
26528d58a37SPierre Jolivet   DMSNES        sdm;
26628d58a37SPierre Jolivet   DMSNES_Local *dmlocalsnes;
26728d58a37SPierre Jolivet 
26828d58a37SPierre Jolivet   PetscFunctionBegin;
26928d58a37SPierre Jolivet   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2709566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
2719566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
27228d58a37SPierre Jolivet   if (func) *func = dmlocalsnes->residuallocal;
27328d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->residuallocalctx;
27428d58a37SPierre Jolivet   PetscFunctionReturn(0);
27528d58a37SPierre Jolivet }
27628d58a37SPierre Jolivet 
27728d58a37SPierre Jolivet /*@C
278f6dfbefdSBarry Smith    DMSNESGetBoundaryLocal - get the local boundary value function set with `DMSNESSetBoundaryLocal()`.
27928d58a37SPierre Jolivet 
28028d58a37SPierre Jolivet    Not Collective
28128d58a37SPierre Jolivet 
2824165533cSJose E. Roman    Input Parameter:
283f6dfbefdSBarry Smith .  dm - `DM` with the associated callback
28428d58a37SPierre Jolivet 
2854165533cSJose E. Roman    Output Parameters:
28628d58a37SPierre Jolivet +  func - local boundary value evaluation
28728d58a37SPierre Jolivet -  ctx - context for local boundary value evaluation
28828d58a37SPierre Jolivet 
28928d58a37SPierre Jolivet    Level: intermediate
29028d58a37SPierre Jolivet 
291db781477SPatrick Sanan .seealso: `DMSNESSetFunctionLocal()`, `DMSNESSetBoundaryLocal()`, `DMDASNESSetJacobianLocal()`
29228d58a37SPierre Jolivet @*/
2939371c9d4SSatish Balay PetscErrorCode DMSNESGetBoundaryLocal(DM dm, PetscErrorCode (**func)(DM, Vec, void *), void **ctx) {
29428d58a37SPierre Jolivet   DMSNES        sdm;
29528d58a37SPierre Jolivet   DMSNES_Local *dmlocalsnes;
29628d58a37SPierre Jolivet 
29728d58a37SPierre Jolivet   PetscFunctionBegin;
29828d58a37SPierre Jolivet   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2999566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
3009566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
30128d58a37SPierre Jolivet   if (func) *func = dmlocalsnes->boundarylocal;
30228d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->boundarylocalctx;
30328d58a37SPierre Jolivet   PetscFunctionReturn(0);
30428d58a37SPierre Jolivet }
30528d58a37SPierre Jolivet 
30628d58a37SPierre Jolivet /*@C
307f6dfbefdSBarry Smith    DMSNESGetJacobianLocal - the local Jacobian evaluation function set with `DMSNESSetJacobianLocal()`.
30828d58a37SPierre Jolivet 
30928d58a37SPierre Jolivet    Logically Collective
31028d58a37SPierre Jolivet 
3114165533cSJose E. Roman    Input Parameter:
312f6dfbefdSBarry Smith .  dm - `DM` with the associated callback
31328d58a37SPierre Jolivet 
3144165533cSJose E. Roman    Output Parameters:
31528d58a37SPierre Jolivet +  func - local Jacobian evaluation
31628d58a37SPierre Jolivet -  ctx - context for local Jacobian evaluation
31728d58a37SPierre Jolivet 
31828d58a37SPierre Jolivet    Level: beginner
31928d58a37SPierre Jolivet 
320db781477SPatrick Sanan .seealso: `DMSNESSetJacobianLocal()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
32128d58a37SPierre Jolivet @*/
3229371c9d4SSatish Balay PetscErrorCode DMSNESGetJacobianLocal(DM dm, PetscErrorCode (**func)(DM, Vec, Mat, Mat, void *), void **ctx) {
32328d58a37SPierre Jolivet   DMSNES        sdm;
32428d58a37SPierre Jolivet   DMSNES_Local *dmlocalsnes;
32528d58a37SPierre Jolivet 
32628d58a37SPierre Jolivet   PetscFunctionBegin;
32728d58a37SPierre Jolivet   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3289566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
3299566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
33028d58a37SPierre Jolivet   if (func) *func = dmlocalsnes->jacobianlocal;
33128d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->jacobianlocalctx;
332ff35dfedSBarry Smith   PetscFunctionReturn(0);
333ff35dfedSBarry Smith }
334