xref: /petsc/src/snes/utils/dmlocalsnes.c (revision 6493148fa132b91628b1fcea3e6c935adfa8483c)
1af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
3ff35dfedSBarry Smith 
4ff35dfedSBarry Smith typedef struct {
5*6493148fSStefano Zampini   PetscErrorCode (*objectivelocal)(DM, Vec, PetscReal *, void *);
6ff35dfedSBarry Smith   PetscErrorCode (*residuallocal)(DM, Vec, Vec, void *);
7d1e9a80fSBarry Smith   PetscErrorCode (*jacobianlocal)(DM, Vec, Mat, Mat, void *);
8bdd6f66aSToby Isaac   PetscErrorCode (*boundarylocal)(DM, Vec, void *);
9*6493148fSStefano Zampini   void *objectivelocalctx;
10ff35dfedSBarry Smith   void *residuallocalctx;
11ff35dfedSBarry Smith   void *jacobianlocalctx;
12bdd6f66aSToby Isaac   void *boundarylocalctx;
13942e3340SBarry Smith } DMSNES_Local;
14ff35dfedSBarry Smith 
15d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDestroy_DMLocal(DMSNES sdm)
16d71ae5a4SJacob Faibussowitsch {
17ff35dfedSBarry Smith   PetscFunctionBegin;
189566063dSJacob Faibussowitsch   PetscCall(PetscFree(sdm->data));
195ed5e208SMatthew G. Knepley   sdm->data = NULL;
203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21ff35dfedSBarry Smith }
22ff35dfedSBarry Smith 
23d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm, DMSNES sdm)
24d71ae5a4SJacob Faibussowitsch {
25ff35dfedSBarry Smith   PetscFunctionBegin;
265ed5e208SMatthew G. Knepley   if (sdm->data != oldsdm->data) {
279566063dSJacob Faibussowitsch     PetscCall(PetscFree(sdm->data));
284dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMSNES_Local **)&sdm->data));
299566063dSJacob Faibussowitsch     if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data, oldsdm->data, sizeof(DMSNES_Local)));
30ff35dfedSBarry Smith   }
313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32ff35dfedSBarry Smith }
33ff35dfedSBarry Smith 
34d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMLocalSNESGetContext(DM dm, DMSNES sdm, DMSNES_Local **dmlocalsnes)
35d71ae5a4SJacob Faibussowitsch {
36ff35dfedSBarry Smith   PetscFunctionBegin;
370298fd71SBarry Smith   *dmlocalsnes = NULL;
38ff35dfedSBarry Smith   if (!sdm->data) {
394dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMSNES_Local **)&sdm->data));
401aa26658SKarl Rupp 
4122c6f798SBarry Smith     sdm->ops->destroy   = DMSNESDestroy_DMLocal;
4222c6f798SBarry Smith     sdm->ops->duplicate = DMSNESDuplicate_DMLocal;
43ff35dfedSBarry Smith   }
44942e3340SBarry Smith   *dmlocalsnes = (DMSNES_Local *)sdm->data;
453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46ff35dfedSBarry Smith }
47ff35dfedSBarry Smith 
48*6493148fSStefano Zampini static PetscErrorCode SNESComputeObjective_DMLocal(SNES snes, Vec X, PetscReal *obj, void *ctx)
49*6493148fSStefano Zampini {
50*6493148fSStefano Zampini   DMSNES_Local *dmlocalsnes = (DMSNES_Local *)ctx;
51*6493148fSStefano Zampini   DM            dm;
52*6493148fSStefano Zampini   Vec           Xloc;
53*6493148fSStefano Zampini   PetscBool     transform;
54*6493148fSStefano Zampini 
55*6493148fSStefano Zampini   PetscFunctionBegin;
56*6493148fSStefano Zampini   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
57*6493148fSStefano Zampini   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
58*6493148fSStefano Zampini   PetscCall(SNESGetDM(snes, &dm));
59*6493148fSStefano Zampini   PetscCall(DMGetLocalVector(dm, &Xloc));
60*6493148fSStefano Zampini   PetscCall(VecZeroEntries(Xloc));
61*6493148fSStefano Zampini   /* Non-conforming routines needs boundary values before G2L */
62*6493148fSStefano Zampini   if (dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
63*6493148fSStefano Zampini   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
64*6493148fSStefano Zampini   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
65*6493148fSStefano Zampini   /* Need to reset boundary values if we transformed */
66*6493148fSStefano Zampini   PetscCall(DMHasBasisTransform(dm, &transform));
67*6493148fSStefano Zampini   if (transform && dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
68*6493148fSStefano Zampini   CHKMEMQ;
69*6493148fSStefano Zampini   PetscCall((*dmlocalsnes->objectivelocal)(dm, Xloc, obj, dmlocalsnes->objectivelocalctx));
70*6493148fSStefano Zampini   CHKMEMQ;
71*6493148fSStefano Zampini   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, obj, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
72*6493148fSStefano Zampini   PetscCall(DMRestoreLocalVector(dm, &Xloc));
73*6493148fSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
74*6493148fSStefano Zampini }
75*6493148fSStefano Zampini 
76d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes, Vec X, Vec F, void *ctx)
77d71ae5a4SJacob Faibussowitsch {
78942e3340SBarry Smith   DMSNES_Local *dmlocalsnes = (DMSNES_Local *)ctx;
797a1e7fa1SMatthew G. Knepley   DM            dm;
80ff35dfedSBarry Smith   Vec           Xloc, Floc;
817a1e7fa1SMatthew G. Knepley   PetscBool     transform;
82ff35dfedSBarry Smith 
83ff35dfedSBarry Smith   PetscFunctionBegin;
84ff35dfedSBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
85ff35dfedSBarry Smith   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
86ff35dfedSBarry Smith   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
879566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
889566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
899566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Floc));
909566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Xloc));
919566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Floc));
927a1e7fa1SMatthew G. Knepley   /* Non-conforming routines needs boundary values before G2L */
939566063dSJacob Faibussowitsch   if (dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
949566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
959566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
967a1e7fa1SMatthew G. Knepley   /* Need to reset boundary values if we transformed */
979566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dm, &transform));
989566063dSJacob Faibussowitsch   if (transform && dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
99ff35dfedSBarry Smith   CHKMEMQ;
1009566063dSJacob Faibussowitsch   PetscCall((*dmlocalsnes->residuallocal)(dm, Xloc, Floc, dmlocalsnes->residuallocalctx));
101ff35dfedSBarry Smith   CHKMEMQ;
1029566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
1039566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
1049566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
1059566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Floc));
1069566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
107afa9dec6SMatthew G. Knepley   {
108afa9dec6SMatthew G. Knepley     char        name[PETSC_MAX_PATH_LEN];
109f2f0fc17SMatthew G. Knepley     char        oldname[PETSC_MAX_PATH_LEN];
110f2f0fc17SMatthew G. Knepley     const char *tmp;
111afa9dec6SMatthew G. Knepley     PetscInt    it;
112afa9dec6SMatthew G. Knepley 
1139566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(snes, &it));
1149566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Solution, Iterate %d", (int)it));
1159566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)X, &tmp));
1169566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(oldname, tmp, PETSC_MAX_PATH_LEN - 1));
1179566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)X, name));
1189566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(X, (PetscObject)snes, "-dmsnes_solution_vec_view"));
1199566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)X, oldname));
1209566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Residual, Iterate %d", (int)it));
1219566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)F, name));
1229566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(F, (PetscObject)snes, "-dmsnes_residual_vec_view"));
123afa9dec6SMatthew G. Knepley   }
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125ff35dfedSBarry Smith }
126ff35dfedSBarry Smith 
127d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes, Vec X, Mat A, Mat B, void *ctx)
128d71ae5a4SJacob Faibussowitsch {
129942e3340SBarry Smith   DMSNES_Local *dmlocalsnes = (DMSNES_Local *)ctx;
1307a1e7fa1SMatthew G. Knepley   DM            dm;
131ff35dfedSBarry Smith   Vec           Xloc;
1327a1e7fa1SMatthew G. Knepley   PetscBool     transform;
133ff35dfedSBarry Smith 
134ff35dfedSBarry Smith   PetscFunctionBegin;
1359566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
136ff35dfedSBarry Smith   if (dmlocalsnes->jacobianlocal) {
1379566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Xloc));
1389566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Xloc));
1397a1e7fa1SMatthew G. Knepley     /* Non-conforming routines needs boundary values before G2L */
1409566063dSJacob Faibussowitsch     if (dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
1419566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1429566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1437a1e7fa1SMatthew G. Knepley     /* Need to reset boundary values if we transformed */
1449566063dSJacob Faibussowitsch     PetscCall(DMHasBasisTransform(dm, &transform));
1459566063dSJacob Faibussowitsch     if (transform && dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm, Xloc, dmlocalsnes->boundarylocalctx));
146ff35dfedSBarry Smith     CHKMEMQ;
1479566063dSJacob Faibussowitsch     PetscCall((*dmlocalsnes->jacobianlocal)(dm, Xloc, A, B, dmlocalsnes->jacobianlocalctx));
148ff35dfedSBarry Smith     CHKMEMQ;
1499566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Xloc));
150ff35dfedSBarry Smith   } else {
151ff35dfedSBarry Smith     MatFDColoring fdcoloring;
1529566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring));
153ff35dfedSBarry Smith     if (!fdcoloring) {
154ff35dfedSBarry Smith       ISColoring coloring;
155ff35dfedSBarry Smith 
1569566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring));
1579566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring));
1589566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&coloring));
159ff35dfedSBarry Smith       switch (dm->coloringtype) {
160d71ae5a4SJacob Faibussowitsch       case IS_COLORING_GLOBAL:
161d71ae5a4SJacob Faibussowitsch         PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))SNESComputeFunction_DMLocal, dmlocalsnes));
162d71ae5a4SJacob Faibussowitsch         break;
163d71ae5a4SJacob Faibussowitsch       default:
164d71ae5a4SJacob Faibussowitsch         SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]);
165ff35dfedSBarry Smith       }
1669566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix));
1679566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1689566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring));
1699566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring));
1709566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)fdcoloring));
171ff35dfedSBarry Smith 
172ff35dfedSBarry Smith       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
173ff35dfedSBarry Smith        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
174ff35dfedSBarry Smith        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
175ff35dfedSBarry Smith        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
176140e18c1SBarry Smith        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
177ff35dfedSBarry Smith        */
1789566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dm));
179ff35dfedSBarry Smith     }
1809566063dSJacob Faibussowitsch     PetscCall(MatFDColoringApply(B, fdcoloring, X, snes));
181ff35dfedSBarry Smith   }
182ff35dfedSBarry Smith   /* This will be redundant if the user called both, but it's too common to forget. */
18394ab13aaSBarry Smith   if (A != B) {
1849566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1859566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
186ff35dfedSBarry Smith   }
1873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
188ff35dfedSBarry Smith }
189ff35dfedSBarry Smith 
190ff35dfedSBarry Smith /*@C
191*6493148fSStefano Zampini   DMSNESSetObjectiveLocal - set a local objective evaluation function. This function is called with local vector
192*6493148fSStefano Zampini   containing the local vector information PLUS ghost point information. It should compute a result for all local
193*6493148fSStefano Zampini   elements and `DMSNES` will automatically accumulate the overlapping values.
194*6493148fSStefano Zampini 
195*6493148fSStefano Zampini   Logically Collective
196*6493148fSStefano Zampini 
197*6493148fSStefano Zampini   Input Parameters:
198*6493148fSStefano Zampini + dm   - `DM` to associate callback with
199*6493148fSStefano Zampini . func - local objective evaluation
200*6493148fSStefano Zampini - ctx  - optional context for local residual evaluation
201*6493148fSStefano Zampini 
202*6493148fSStefano Zampini   Level: advanced
203*6493148fSStefano Zampini 
204*6493148fSStefano Zampini .seealso: `DMSNESSetFunctionLocal()`, `DMSNESSetJacobianLocal()`
205*6493148fSStefano Zampini @*/
206*6493148fSStefano Zampini PetscErrorCode DMSNESSetObjectiveLocal(DM dm, PetscErrorCode (*func)(DM, Vec, PetscReal *, void *), void *ctx)
207*6493148fSStefano Zampini {
208*6493148fSStefano Zampini   DMSNES        sdm;
209*6493148fSStefano Zampini   DMSNES_Local *dmlocalsnes;
210*6493148fSStefano Zampini 
211*6493148fSStefano Zampini   PetscFunctionBegin;
212*6493148fSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
213*6493148fSStefano Zampini   PetscCall(DMGetDMSNESWrite(dm, &sdm));
214*6493148fSStefano Zampini   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
215*6493148fSStefano Zampini 
216*6493148fSStefano Zampini   dmlocalsnes->objectivelocal    = func;
217*6493148fSStefano Zampini   dmlocalsnes->objectivelocalctx = ctx;
218*6493148fSStefano Zampini 
219*6493148fSStefano Zampini   PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMLocal, dmlocalsnes));
220*6493148fSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
221*6493148fSStefano Zampini }
222*6493148fSStefano Zampini 
223*6493148fSStefano Zampini /*@C
224ff35dfedSBarry Smith   DMSNESSetFunctionLocal - set a local residual evaluation function. This function is called with local vector
225ff35dfedSBarry Smith   containing the local vector information PLUS ghost point information. It should compute a result for all local
226f6dfbefdSBarry Smith   elements and `DMSNES` will automatically accumulate the overlapping values.
227ff35dfedSBarry Smith 
228ff35dfedSBarry Smith   Logically Collective
229ff35dfedSBarry Smith 
2304165533cSJose E. Roman   Input Parameters:
231f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
232ff35dfedSBarry Smith . func - local residual evaluation
233ff35dfedSBarry Smith - ctx  - optional context for local residual evaluation
234ff35dfedSBarry Smith 
235420bcc1bSBarry Smith   Calling sequence of `func`:
236420bcc1bSBarry Smith + dm  - `DM` for the function
237420bcc1bSBarry Smith . x   - vector to state at which to evaluate residual
238420bcc1bSBarry Smith . f   - vector to hold the function evaluation
239420bcc1bSBarry Smith - ctx - optional context passed above
240420bcc1bSBarry Smith 
241f6dfbefdSBarry Smith   Level: advanced
242ff35dfedSBarry Smith 
243*6493148fSStefano Zampini .seealso: [](ch_snes), `DMSNESSetFunction()`, `DMSNESSetJacobianLocal()`
244ff35dfedSBarry Smith @*/
245420bcc1bSBarry Smith PetscErrorCode DMSNESSetFunctionLocal(DM dm, PetscErrorCode (*func)(DM dm, Vec x, Vec f, void *ctx), void *ctx)
246d71ae5a4SJacob Faibussowitsch {
247942e3340SBarry Smith   DMSNES        sdm;
248942e3340SBarry Smith   DMSNES_Local *dmlocalsnes;
249ff35dfedSBarry Smith 
250ff35dfedSBarry Smith   PetscFunctionBegin;
251ff35dfedSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2529566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2539566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
2541aa26658SKarl Rupp 
255ff35dfedSBarry Smith   dmlocalsnes->residuallocal    = func;
256ff35dfedSBarry Smith   dmlocalsnes->residuallocalctx = ctx;
2571aa26658SKarl Rupp 
2589566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMLocal, dmlocalsnes));
25922c6f798SBarry Smith   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
2609566063dSJacob Faibussowitsch     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMLocal, dmlocalsnes));
261ff35dfedSBarry Smith   }
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
263ff35dfedSBarry Smith }
264ff35dfedSBarry Smith 
265bdd6f66aSToby Isaac /*@C
2660b4db180SJacob Faibussowitsch   DMSNESSetBoundaryLocal - set a function to insert, for example, essential boundary conditions into a ghosted solution vector
267bdd6f66aSToby Isaac 
268bdd6f66aSToby Isaac   Logically Collective
269bdd6f66aSToby Isaac 
2704165533cSJose E. Roman   Input Parameters:
271f6dfbefdSBarry Smith + dm   - `DM` to associate callback with
272bdd6f66aSToby Isaac . func - local boundary value evaluation
273bdd6f66aSToby Isaac - ctx  - optional context for local boundary value evaluation
274bdd6f66aSToby Isaac 
2750b4db180SJacob Faibussowitsch   Calling sequence of `func`:
2760b4db180SJacob Faibussowitsch + dm  - the `DM` context
277baca6076SPierre Jolivet . X   - ghosted solution vector, appropriate locations (such as essential boundary condition nodes) should be filled
278420bcc1bSBarry Smith - ctx - option context passed in `DMSNESSetBoundaryLocal()`
2790b4db180SJacob Faibussowitsch 
280f6dfbefdSBarry Smith   Level: advanced
281bdd6f66aSToby Isaac 
282*6493148fSStefano Zampini .seealso: [](ch_snes), `DMSNESSetObjectiveLocal()`, `DMSNESSetFunctionLocal()`, `DMSNESSetJacobianLocal()`
283bdd6f66aSToby Isaac @*/
2840b4db180SJacob Faibussowitsch PetscErrorCode DMSNESSetBoundaryLocal(DM dm, PetscErrorCode (*func)(DM dm, Vec X, void *ctx), void *ctx)
285d71ae5a4SJacob Faibussowitsch {
286bdd6f66aSToby Isaac   DMSNES        sdm;
287bdd6f66aSToby Isaac   DMSNES_Local *dmlocalsnes;
288bdd6f66aSToby Isaac 
289bdd6f66aSToby Isaac   PetscFunctionBegin;
290bdd6f66aSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2919566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2929566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
293bdd6f66aSToby Isaac 
294bdd6f66aSToby Isaac   dmlocalsnes->boundarylocal    = func;
295bdd6f66aSToby Isaac   dmlocalsnes->boundarylocalctx = ctx;
296bdd6f66aSToby Isaac 
2973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
298bdd6f66aSToby Isaac }
299bdd6f66aSToby Isaac 
300ff35dfedSBarry Smith /*@C
301ff35dfedSBarry Smith   DMSNESSetJacobianLocal - set a local Jacobian evaluation function
302ff35dfedSBarry Smith 
303ff35dfedSBarry Smith   Logically Collective
304ff35dfedSBarry Smith 
3054165533cSJose E. Roman   Input Parameters:
306420bcc1bSBarry Smith + dm   - `DM` to associate callback with
307ff35dfedSBarry Smith . func - local Jacobian evaluation
308ff35dfedSBarry Smith - ctx  - optional context for local Jacobian evaluation
309ff35dfedSBarry Smith 
3100b4db180SJacob Faibussowitsch   Calling sequence of `func`:
3110b4db180SJacob Faibussowitsch + dm  - the `DM` context
3120b4db180SJacob Faibussowitsch . X   - current solution vector (ghosted or not?)
3130b4db180SJacob Faibussowitsch . J   - the Jacobian
3140b4db180SJacob Faibussowitsch . Jp  - approximate Jacobian used to compute the preconditioner, often `J`
3150b4db180SJacob Faibussowitsch - ctx - a user provided context
3160b4db180SJacob Faibussowitsch 
317f6dfbefdSBarry Smith   Level: advanced
318ff35dfedSBarry Smith 
319*6493148fSStefano Zampini .seealso: [](ch_snes), `DMSNESSetObjectiveLocal()`, `DMSNESSetFunctionLocal()`, `DMSNESSetBoundaryLocal()`
320ff35dfedSBarry Smith @*/
3210b4db180SJacob Faibussowitsch PetscErrorCode DMSNESSetJacobianLocal(DM dm, PetscErrorCode (*func)(DM dm, Vec X, Mat J, Mat Jp, void *ctx), void *ctx)
322d71ae5a4SJacob Faibussowitsch {
323942e3340SBarry Smith   DMSNES        sdm;
324942e3340SBarry Smith   DMSNES_Local *dmlocalsnes;
325ff35dfedSBarry Smith 
326ff35dfedSBarry Smith   PetscFunctionBegin;
327ff35dfedSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3289566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3299566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
3301aa26658SKarl Rupp 
331ff35dfedSBarry Smith   dmlocalsnes->jacobianlocal    = func;
332ff35dfedSBarry Smith   dmlocalsnes->jacobianlocalctx = ctx;
3331aa26658SKarl Rupp 
3349566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMLocal, dmlocalsnes));
3353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33628d58a37SPierre Jolivet }
33728d58a37SPierre Jolivet 
33828d58a37SPierre Jolivet /*@C
339*6493148fSStefano Zampini   DMSNESGetObjectiveLocal - get the local objective evaluation function information set with `DMSNESSetObjectiveLocal()`.
340*6493148fSStefano Zampini 
341*6493148fSStefano Zampini   Not Collective
342*6493148fSStefano Zampini 
343*6493148fSStefano Zampini   Input Parameter:
344*6493148fSStefano Zampini . dm - `DM` with the associated callback
345*6493148fSStefano Zampini 
346*6493148fSStefano Zampini   Output Parameters:
347*6493148fSStefano Zampini + func - local objective evaluation
348*6493148fSStefano Zampini - ctx  - context for local residual evaluation
349*6493148fSStefano Zampini 
350*6493148fSStefano Zampini   Level: beginner
351*6493148fSStefano Zampini 
352*6493148fSStefano Zampini .seealso: `DMSNESSetObjective()`, `DMSNESSetObjectiveLocal()`, `DMSNESSetFunctionLocal()`
353*6493148fSStefano Zampini @*/
354*6493148fSStefano Zampini PetscErrorCode DMSNESGetObjectiveLocal(DM dm, PetscErrorCode (**func)(DM, Vec, PetscReal *, void *), void **ctx)
355*6493148fSStefano Zampini {
356*6493148fSStefano Zampini   DMSNES        sdm;
357*6493148fSStefano Zampini   DMSNES_Local *dmlocalsnes;
358*6493148fSStefano Zampini 
359*6493148fSStefano Zampini   PetscFunctionBegin;
360*6493148fSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
361*6493148fSStefano Zampini   PetscCall(DMGetDMSNES(dm, &sdm));
362*6493148fSStefano Zampini   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
363*6493148fSStefano Zampini   if (func) *func = dmlocalsnes->objectivelocal;
364*6493148fSStefano Zampini   if (ctx) *ctx = dmlocalsnes->objectivelocalctx;
365*6493148fSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
366*6493148fSStefano Zampini }
367*6493148fSStefano Zampini 
368*6493148fSStefano Zampini /*@C
369f6dfbefdSBarry Smith   DMSNESGetFunctionLocal - get the local residual evaluation function information set with `DMSNESSetFunctionLocal()`.
37028d58a37SPierre Jolivet 
37128d58a37SPierre Jolivet   Not Collective
37228d58a37SPierre Jolivet 
3734165533cSJose E. Roman   Input Parameter:
374f6dfbefdSBarry Smith . dm - `DM` with the associated callback
37528d58a37SPierre Jolivet 
3764165533cSJose E. Roman   Output Parameters:
37728d58a37SPierre Jolivet + func - local residual evaluation
37828d58a37SPierre Jolivet - ctx  - context for local residual evaluation
37928d58a37SPierre Jolivet 
38028d58a37SPierre Jolivet   Level: beginner
38128d58a37SPierre Jolivet 
382*6493148fSStefano Zampini .seealso: [](ch_snes), `DMSNESSetFunction()`, `DMSNESSetFunctionLocal()`, `DMSNESSetJacobianLocal()`
38328d58a37SPierre Jolivet @*/
384d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESGetFunctionLocal(DM dm, PetscErrorCode (**func)(DM, Vec, Vec, void *), void **ctx)
385d71ae5a4SJacob Faibussowitsch {
38628d58a37SPierre Jolivet   DMSNES        sdm;
38728d58a37SPierre Jolivet   DMSNES_Local *dmlocalsnes;
38828d58a37SPierre Jolivet 
38928d58a37SPierre Jolivet   PetscFunctionBegin;
39028d58a37SPierre Jolivet   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3919566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
3929566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
39328d58a37SPierre Jolivet   if (func) *func = dmlocalsnes->residuallocal;
39428d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->residuallocalctx;
3953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39628d58a37SPierre Jolivet }
39728d58a37SPierre Jolivet 
39828d58a37SPierre Jolivet /*@C
399f6dfbefdSBarry Smith   DMSNESGetBoundaryLocal - get the local boundary value function set with `DMSNESSetBoundaryLocal()`.
40028d58a37SPierre Jolivet 
40128d58a37SPierre Jolivet   Not Collective
40228d58a37SPierre Jolivet 
4034165533cSJose E. Roman   Input Parameter:
404f6dfbefdSBarry Smith . dm - `DM` with the associated callback
40528d58a37SPierre Jolivet 
4064165533cSJose E. Roman   Output Parameters:
40728d58a37SPierre Jolivet + func - local boundary value evaluation
40828d58a37SPierre Jolivet - ctx  - context for local boundary value evaluation
40928d58a37SPierre Jolivet 
41028d58a37SPierre Jolivet   Level: intermediate
41128d58a37SPierre Jolivet 
412*6493148fSStefano Zampini .seealso: [](ch_snes), `DMSNESSetFunctionLocal()`, `DMSNESSetBoundaryLocal()`, `DMSNESSetJacobianLocal()`
41328d58a37SPierre Jolivet @*/
414d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESGetBoundaryLocal(DM dm, PetscErrorCode (**func)(DM, Vec, void *), void **ctx)
415d71ae5a4SJacob Faibussowitsch {
41628d58a37SPierre Jolivet   DMSNES        sdm;
41728d58a37SPierre Jolivet   DMSNES_Local *dmlocalsnes;
41828d58a37SPierre Jolivet 
41928d58a37SPierre Jolivet   PetscFunctionBegin;
42028d58a37SPierre Jolivet   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4219566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
4229566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
42328d58a37SPierre Jolivet   if (func) *func = dmlocalsnes->boundarylocal;
42428d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->boundarylocalctx;
4253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42628d58a37SPierre Jolivet }
42728d58a37SPierre Jolivet 
42828d58a37SPierre Jolivet /*@C
429f6dfbefdSBarry Smith   DMSNESGetJacobianLocal - the local Jacobian evaluation function set with `DMSNESSetJacobianLocal()`.
43028d58a37SPierre Jolivet 
43128d58a37SPierre Jolivet   Logically Collective
43228d58a37SPierre Jolivet 
4334165533cSJose E. Roman   Input Parameter:
434f6dfbefdSBarry Smith . dm - `DM` with the associated callback
43528d58a37SPierre Jolivet 
4364165533cSJose E. Roman   Output Parameters:
43728d58a37SPierre Jolivet + func - local Jacobian evaluation
43828d58a37SPierre Jolivet - ctx  - context for local Jacobian evaluation
43928d58a37SPierre Jolivet 
44028d58a37SPierre Jolivet   Level: beginner
44128d58a37SPierre Jolivet 
442*6493148fSStefano Zampini .seealso: [](ch_snes), `DMSNESSetJacobianLocal()`, `DMSNESSetJacobian()`
44328d58a37SPierre Jolivet @*/
444d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESGetJacobianLocal(DM dm, PetscErrorCode (**func)(DM, Vec, Mat, Mat, void *), void **ctx)
445d71ae5a4SJacob Faibussowitsch {
44628d58a37SPierre Jolivet   DMSNES        sdm;
44728d58a37SPierre Jolivet   DMSNES_Local *dmlocalsnes;
44828d58a37SPierre Jolivet 
44928d58a37SPierre Jolivet   PetscFunctionBegin;
45028d58a37SPierre Jolivet   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4519566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
4529566063dSJacob Faibussowitsch   PetscCall(DMLocalSNESGetContext(dm, sdm, &dmlocalsnes));
45328d58a37SPierre Jolivet   if (func) *func = dmlocalsnes->jacobianlocal;
45428d58a37SPierre Jolivet   if (ctx) *ctx = dmlocalsnes->jacobianlocalctx;
4553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
456ff35dfedSBarry Smith }
457