xref: /petsc/src/snes/utils/dmlocalsnes.c (revision 420bcc1b905230dede3c88f397d1a4e60493adde)
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 
13d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDestroy_DMLocal(DMSNES sdm)
14d71ae5a4SJacob Faibussowitsch {
15ff35dfedSBarry Smith   PetscFunctionBegin;
169566063dSJacob Faibussowitsch   PetscCall(PetscFree(sdm->data));
175ed5e208SMatthew G. Knepley   sdm->data = NULL;
183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19ff35dfedSBarry Smith }
20ff35dfedSBarry Smith 
21d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm, DMSNES sdm)
22d71ae5a4SJacob 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   }
293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30ff35dfedSBarry Smith }
31ff35dfedSBarry Smith 
32d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMLocalSNESGetContext(DM dm, DMSNES sdm, DMSNES_Local **dmlocalsnes)
33d71ae5a4SJacob 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;
433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44ff35dfedSBarry Smith }
45ff35dfedSBarry Smith 
46d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes, Vec X, Vec F, void *ctx)
47d71ae5a4SJacob 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   }
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95ff35dfedSBarry Smith }
96ff35dfedSBarry Smith 
97d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes, Vec X, Mat A, Mat B, void *ctx)
98d71ae5a4SJacob 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) {
130d71ae5a4SJacob Faibussowitsch       case IS_COLORING_GLOBAL:
131d71ae5a4SJacob Faibussowitsch         PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))SNESComputeFunction_DMLocal, dmlocalsnes));
132d71ae5a4SJacob Faibussowitsch         break;
133d71ae5a4SJacob Faibussowitsch       default:
134d71ae5a4SJacob 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   }
1573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 
172*420bcc1bSBarry Smith   Calling sequence of `func`:
173*420bcc1bSBarry Smith + dm  - `DM` for the function
174*420bcc1bSBarry Smith . x   - vector to state at which to evaluate residual
175*420bcc1bSBarry Smith . f   - vector to hold the function evaluation
176*420bcc1bSBarry Smith - ctx - optional context passed above
177*420bcc1bSBarry Smith 
178f6dfbefdSBarry Smith   Level: advanced
179ff35dfedSBarry Smith 
180*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
181ff35dfedSBarry Smith @*/
182*420bcc1bSBarry Smith PetscErrorCode DMSNESSetFunctionLocal(DM dm, PetscErrorCode (*func)(DM dm, Vec x, Vec f, void *ctx), void *ctx)
183d71ae5a4SJacob Faibussowitsch {
184942e3340SBarry Smith   DMSNES        sdm;
185942e3340SBarry Smith   DMSNES_Local *dmlocalsnes;
186ff35dfedSBarry Smith 
187ff35dfedSBarry Smith   PetscFunctionBegin;
188ff35dfedSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1899566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
1909566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
1911aa26658SKarl Rupp 
192ff35dfedSBarry Smith   dmlocalsnes->residuallocal    = func;
193ff35dfedSBarry Smith   dmlocalsnes->residuallocalctx = ctx;
1941aa26658SKarl Rupp 
1959566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMLocal, dmlocalsnes));
19622c6f798SBarry Smith   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
1979566063dSJacob Faibussowitsch     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMLocal, dmlocalsnes));
198ff35dfedSBarry Smith   }
1993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
200ff35dfedSBarry Smith }
201ff35dfedSBarry Smith 
202bdd6f66aSToby Isaac /*@C
2030b4db180SJacob Faibussowitsch   DMSNESSetBoundaryLocal - set a function to insert, for example, essential boundary conditions into a ghosted solution vector
204bdd6f66aSToby Isaac 
205bdd6f66aSToby Isaac   Logically Collective
206bdd6f66aSToby Isaac 
2074165533cSJose E. Roman   Input Parameters:
208f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
209bdd6f66aSToby Isaac . func - local boundary value evaluation
210bdd6f66aSToby Isaac - ctx  - optional context for local boundary value evaluation
211bdd6f66aSToby Isaac 
2120b4db180SJacob Faibussowitsch   Calling sequence of `func`:
2130b4db180SJacob Faibussowitsch + dm  - the `DM` context
214baca6076SPierre Jolivet . X   - ghosted solution vector, appropriate locations (such as essential boundary condition nodes) should be filled
215*420bcc1bSBarry Smith - ctx - option context passed in `DMSNESSetBoundaryLocal()`
2160b4db180SJacob Faibussowitsch 
217f6dfbefdSBarry Smith   Level: advanced
218bdd6f66aSToby Isaac 
219*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMSNESSetFunctionLocal()`, `DMDASNESSetJacobianLocal()`
220bdd6f66aSToby Isaac @*/
2210b4db180SJacob Faibussowitsch PetscErrorCode DMSNESSetBoundaryLocal(DM dm, PetscErrorCode (*func)(DM dm, Vec X, void *ctx), void *ctx)
222d71ae5a4SJacob Faibussowitsch {
223bdd6f66aSToby Isaac   DMSNES        sdm;
224bdd6f66aSToby Isaac   DMSNES_Local *dmlocalsnes;
225bdd6f66aSToby Isaac 
226bdd6f66aSToby Isaac   PetscFunctionBegin;
227bdd6f66aSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2289566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2299566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
230bdd6f66aSToby Isaac 
231bdd6f66aSToby Isaac   dmlocalsnes->boundarylocal    = func;
232bdd6f66aSToby Isaac   dmlocalsnes->boundarylocalctx = ctx;
233bdd6f66aSToby Isaac 
2343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
235bdd6f66aSToby Isaac }
236bdd6f66aSToby Isaac 
237ff35dfedSBarry Smith /*@C
238ff35dfedSBarry Smith   DMSNESSetJacobianLocal - set a local Jacobian evaluation function
239ff35dfedSBarry Smith 
240ff35dfedSBarry Smith   Logically Collective
241ff35dfedSBarry Smith 
2424165533cSJose E. Roman   Input Parameters:
243*420bcc1bSBarry Smith + dm   - `DM` to associate callback with
244ff35dfedSBarry Smith . func - local Jacobian evaluation
245ff35dfedSBarry Smith - ctx  - optional context for local Jacobian evaluation
246ff35dfedSBarry Smith 
2470b4db180SJacob Faibussowitsch   Calling sequence of `func`:
2480b4db180SJacob Faibussowitsch + dm  - the `DM` context
2490b4db180SJacob Faibussowitsch . X   - current solution vector (ghosted or not?)
2500b4db180SJacob Faibussowitsch . J   - the Jacobian
2510b4db180SJacob Faibussowitsch . Jp  - approximate Jacobian used to compute the preconditioner, often `J`
2520b4db180SJacob Faibussowitsch - ctx - a user provided context
2530b4db180SJacob Faibussowitsch 
254f6dfbefdSBarry Smith   Level: advanced
255ff35dfedSBarry Smith 
256*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMSNESSetJacobian()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
257ff35dfedSBarry Smith @*/
2580b4db180SJacob Faibussowitsch PetscErrorCode DMSNESSetJacobianLocal(DM dm, PetscErrorCode (*func)(DM dm, Vec X, Mat J, Mat Jp, void *ctx), void *ctx)
259d71ae5a4SJacob Faibussowitsch {
260942e3340SBarry Smith   DMSNES        sdm;
261942e3340SBarry Smith   DMSNES_Local *dmlocalsnes;
262ff35dfedSBarry Smith 
263ff35dfedSBarry Smith   PetscFunctionBegin;
264ff35dfedSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2659566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2669566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
2671aa26658SKarl Rupp 
268ff35dfedSBarry Smith   dmlocalsnes->jacobianlocal    = func;
269ff35dfedSBarry Smith   dmlocalsnes->jacobianlocalctx = ctx;
2701aa26658SKarl Rupp 
2719566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMLocal, dmlocalsnes));
2723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27328d58a37SPierre Jolivet }
27428d58a37SPierre Jolivet 
27528d58a37SPierre Jolivet /*@C
276f6dfbefdSBarry Smith   DMSNESGetFunctionLocal - get the local residual evaluation function information set with `DMSNESSetFunctionLocal()`.
27728d58a37SPierre Jolivet 
27828d58a37SPierre Jolivet   Not Collective
27928d58a37SPierre Jolivet 
2804165533cSJose E. Roman   Input Parameter:
281f6dfbefdSBarry Smith . dm - `DM` with the associated callback
28228d58a37SPierre Jolivet 
2834165533cSJose E. Roman   Output Parameters:
28428d58a37SPierre Jolivet + func - local residual evaluation
28528d58a37SPierre Jolivet - ctx  - context for local residual evaluation
28628d58a37SPierre Jolivet 
28728d58a37SPierre Jolivet   Level: beginner
28828d58a37SPierre Jolivet 
289*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMSNESSetFunction()`, `DMSNESSetFunctionLocal()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
29028d58a37SPierre Jolivet @*/
291d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESGetFunctionLocal(DM dm, PetscErrorCode (**func)(DM, Vec, Vec, void *), void **ctx)
292d71ae5a4SJacob Faibussowitsch {
29328d58a37SPierre Jolivet   DMSNES        sdm;
29428d58a37SPierre Jolivet   DMSNES_Local *dmlocalsnes;
29528d58a37SPierre Jolivet 
29628d58a37SPierre Jolivet   PetscFunctionBegin;
29728d58a37SPierre Jolivet   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2989566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
2999566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
30028d58a37SPierre Jolivet   if (func) *func = dmlocalsnes->residuallocal;
30128d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->residuallocalctx;
3023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30328d58a37SPierre Jolivet }
30428d58a37SPierre Jolivet 
30528d58a37SPierre Jolivet /*@C
306f6dfbefdSBarry Smith   DMSNESGetBoundaryLocal - get the local boundary value function set with `DMSNESSetBoundaryLocal()`.
30728d58a37SPierre Jolivet 
30828d58a37SPierre Jolivet   Not Collective
30928d58a37SPierre Jolivet 
3104165533cSJose E. Roman   Input Parameter:
311f6dfbefdSBarry Smith . dm - `DM` with the associated callback
31228d58a37SPierre Jolivet 
3134165533cSJose E. Roman   Output Parameters:
31428d58a37SPierre Jolivet + func - local boundary value evaluation
31528d58a37SPierre Jolivet - ctx  - context for local boundary value evaluation
31628d58a37SPierre Jolivet 
31728d58a37SPierre Jolivet   Level: intermediate
31828d58a37SPierre Jolivet 
319*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMSNESSetFunctionLocal()`, `DMSNESSetBoundaryLocal()`, `DMDASNESSetJacobianLocal()`
32028d58a37SPierre Jolivet @*/
321d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESGetBoundaryLocal(DM dm, PetscErrorCode (**func)(DM, Vec, void *), void **ctx)
322d71ae5a4SJacob Faibussowitsch {
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->boundarylocal;
33128d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->boundarylocalctx;
3323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33328d58a37SPierre Jolivet }
33428d58a37SPierre Jolivet 
33528d58a37SPierre Jolivet /*@C
336f6dfbefdSBarry Smith   DMSNESGetJacobianLocal - the local Jacobian evaluation function set with `DMSNESSetJacobianLocal()`.
33728d58a37SPierre Jolivet 
33828d58a37SPierre Jolivet   Logically Collective
33928d58a37SPierre Jolivet 
3404165533cSJose E. Roman   Input Parameter:
341f6dfbefdSBarry Smith . dm - `DM` with the associated callback
34228d58a37SPierre Jolivet 
3434165533cSJose E. Roman   Output Parameters:
34428d58a37SPierre Jolivet + func - local Jacobian evaluation
34528d58a37SPierre Jolivet - ctx  - context for local Jacobian evaluation
34628d58a37SPierre Jolivet 
34728d58a37SPierre Jolivet   Level: beginner
34828d58a37SPierre Jolivet 
349*420bcc1bSBarry Smith .seealso: [](ch_snes), `DMSNESSetJacobianLocal()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
35028d58a37SPierre Jolivet @*/
351d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESGetJacobianLocal(DM dm, PetscErrorCode (**func)(DM, Vec, Mat, Mat, void *), void **ctx)
352d71ae5a4SJacob Faibussowitsch {
35328d58a37SPierre Jolivet   DMSNES        sdm;
35428d58a37SPierre Jolivet   DMSNES_Local *dmlocalsnes;
35528d58a37SPierre Jolivet 
35628d58a37SPierre Jolivet   PetscFunctionBegin;
35728d58a37SPierre Jolivet   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3589566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
3599566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
36028d58a37SPierre Jolivet   if (func) *func = dmlocalsnes->jacobianlocal;
36128d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->jacobianlocalctx;
3623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
363ff35dfedSBarry Smith }
364