xref: /petsc/src/snes/utils/dmsnes.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2af0996ceSBarry Smith #include <petsc/private/dmimpl.h>   /*I "petscdm.h" I*/
36cab3a1bSJed Brown 
49371c9d4SSatish Balay static PetscErrorCode DMSNESUnsetFunctionContext_DMSNES(DMSNES sdm) {
5800f99ffSJeremy L Thompson   PetscFunctionBegin;
6800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)sdm, "function ctx", NULL));
7800f99ffSJeremy L Thompson   sdm->functionctxcontainer = NULL;
8800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
9800f99ffSJeremy L Thompson }
10800f99ffSJeremy L Thompson 
119371c9d4SSatish Balay static PetscErrorCode DMSNESUnsetJacobianContext_DMSNES(DMSNES sdm) {
12800f99ffSJeremy L Thompson   PetscFunctionBegin;
13800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)sdm, "jacobian ctx", NULL));
14800f99ffSJeremy L Thompson   sdm->jacobianctxcontainer = NULL;
15800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
16800f99ffSJeremy L Thompson }
17800f99ffSJeremy L Thompson 
189371c9d4SSatish Balay static PetscErrorCode DMSNESDestroy(DMSNES *kdm) {
196cab3a1bSJed Brown   PetscFunctionBegin;
2022c6f798SBarry Smith   if (!*kdm) PetscFunctionReturn(0);
2122c6f798SBarry Smith   PetscValidHeaderSpecific((*kdm), DMSNES_CLASSID, 1);
229371c9d4SSatish Balay   if (--((PetscObject)(*kdm))->refct > 0) {
239371c9d4SSatish Balay     *kdm = NULL;
249371c9d4SSatish Balay     PetscFunctionReturn(0);
259371c9d4SSatish Balay   }
26800f99ffSJeremy L Thompson   PetscCall(DMSNESUnsetFunctionContext_DMSNES(*kdm));
27800f99ffSJeremy L Thompson   PetscCall(DMSNESUnsetJacobianContext_DMSNES(*kdm));
289566063dSJacob Faibussowitsch   if ((*kdm)->ops->destroy) PetscCall(((*kdm)->ops->destroy)(*kdm));
299566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(kdm));
3022c6f798SBarry Smith   PetscFunctionReturn(0);
3122c6f798SBarry Smith }
3222c6f798SBarry Smith 
339371c9d4SSatish Balay PetscErrorCode DMSNESLoad(DMSNES kdm, PetscViewer viewer) {
342d53ad75SBarry Smith   PetscFunctionBegin;
359566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->computefunction, 1, NULL, PETSC_FUNCTION));
369566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->computejacobian, 1, NULL, PETSC_FUNCTION));
372d53ad75SBarry Smith   PetscFunctionReturn(0);
382d53ad75SBarry Smith }
392d53ad75SBarry Smith 
409371c9d4SSatish Balay PetscErrorCode DMSNESView(DMSNES kdm, PetscViewer viewer) {
412d53ad75SBarry Smith   PetscBool isascii, isbinary;
422d53ad75SBarry Smith 
432d53ad75SBarry Smith   PetscFunctionBegin;
449566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
459566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
462d53ad75SBarry Smith   if (isascii) {
474b5d3663SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) && defined(PETSC_SERIALIZE_FUNCTIONS_VIEW)
482d53ad75SBarry Smith     const char *fname;
492d53ad75SBarry Smith 
509566063dSJacob Faibussowitsch     PetscCall(PetscFPTFind(kdm->ops->computefunction, &fname));
5148a46eb9SPierre Jolivet     if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, "Function used by SNES: %s\n", fname));
529566063dSJacob Faibussowitsch     PetscCall(PetscFPTFind(kdm->ops->computejacobian, &fname));
5348a46eb9SPierre Jolivet     if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, "Jacobian function used by SNES: %s\n", fname));
54c7a10e08SBarry Smith #endif
552d53ad75SBarry Smith   } else if (isbinary) {
563964eb88SJed Brown     struct {
573964eb88SJed Brown       PetscErrorCode (*func)(SNES, Vec, Vec, void *);
589200755eSBarry Smith     } funcstruct;
599200755eSBarry Smith     struct {
60d1e9a80fSBarry Smith       PetscErrorCode (*jac)(SNES, Vec, Mat, Mat, void *);
619200755eSBarry Smith     } jacstruct;
629200755eSBarry Smith     funcstruct.func = kdm->ops->computefunction;
639200755eSBarry Smith     jacstruct.jac   = kdm->ops->computejacobian;
649566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer, &funcstruct, 1, PETSC_FUNCTION));
659566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer, &jacstruct, 1, PETSC_FUNCTION));
662d53ad75SBarry Smith   }
672d53ad75SBarry Smith   PetscFunctionReturn(0);
682d53ad75SBarry Smith }
692d53ad75SBarry Smith 
709371c9d4SSatish Balay static PetscErrorCode DMSNESCreate(MPI_Comm comm, DMSNES *kdm) {
7122c6f798SBarry Smith   PetscFunctionBegin;
729566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
739566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*kdm, DMSNES_CLASSID, "DMSNES", "DMSNES", "DMSNES", comm, DMSNESDestroy, DMSNESView));
746cab3a1bSJed Brown   PetscFunctionReturn(0);
756cab3a1bSJed Brown }
766cab3a1bSJed Brown 
77942e3340SBarry Smith /* Attaches the DMSNES to the coarse level.
786cab3a1bSJed Brown  * Under what conditions should we copy versus duplicate?
796cab3a1bSJed Brown  */
809371c9d4SSatish Balay static PetscErrorCode DMCoarsenHook_DMSNES(DM dm, DM dmc, void *ctx) {
816cab3a1bSJed Brown   PetscFunctionBegin;
829566063dSJacob Faibussowitsch   PetscCall(DMCopyDMSNES(dm, dmc));
836cab3a1bSJed Brown   PetscFunctionReturn(0);
846cab3a1bSJed Brown }
856cab3a1bSJed Brown 
86dfe15315SJed Brown /* This could restrict auxiliary information to the coarse level.
87caa4e7f2SJed Brown  */
889371c9d4SSatish Balay static PetscErrorCode DMRestrictHook_DMSNES(DM dm, Mat Restrict, Vec rscale, Mat Inject, DM dmc, void *ctx) {
89caa4e7f2SJed Brown   PetscFunctionBegin;
90caa4e7f2SJed Brown   PetscFunctionReturn(0);
91caa4e7f2SJed Brown }
92caa4e7f2SJed Brown 
93be081cd6SPeter Brune /* Attaches the DMSNES to the subdomain. */
949371c9d4SSatish Balay static PetscErrorCode DMSubDomainHook_DMSNES(DM dm, DM subdm, void *ctx) {
95be081cd6SPeter Brune   PetscFunctionBegin;
969566063dSJacob Faibussowitsch   PetscCall(DMCopyDMSNES(dm, subdm));
97be081cd6SPeter Brune   PetscFunctionReturn(0);
98be081cd6SPeter Brune }
99be081cd6SPeter Brune 
100be081cd6SPeter Brune /* This could restrict auxiliary information to the coarse level.
101be081cd6SPeter Brune  */
1029371c9d4SSatish Balay static PetscErrorCode DMSubDomainRestrictHook_DMSNES(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) {
103be081cd6SPeter Brune   PetscFunctionBegin;
104be081cd6SPeter Brune   PetscFunctionReturn(0);
105be081cd6SPeter Brune }
106be081cd6SPeter Brune 
1079371c9d4SSatish Balay static PetscErrorCode DMRefineHook_DMSNES(DM dm, DM dmf, void *ctx) {
10803a0fabfSPeter Brune   PetscFunctionBegin;
1099566063dSJacob Faibussowitsch   PetscCall(DMCopyDMSNES(dm, dmf));
11003a0fabfSPeter Brune   PetscFunctionReturn(0);
11103a0fabfSPeter Brune }
11203a0fabfSPeter Brune 
11303a0fabfSPeter Brune /* This could restrict auxiliary information to the coarse level.
11403a0fabfSPeter Brune  */
1159371c9d4SSatish Balay static PetscErrorCode DMInterpolateHook_DMSNES(DM dm, Mat Interp, DM dmf, void *ctx) {
11603a0fabfSPeter Brune   PetscFunctionBegin;
11703a0fabfSPeter Brune   PetscFunctionReturn(0);
11803a0fabfSPeter Brune }
11903a0fabfSPeter Brune 
12022c6f798SBarry Smith /*@C
121*f6dfbefdSBarry Smith    DMSNESCopy - copies the information in a `DMSNES` to another `DMSNES`
12222c6f798SBarry Smith 
12322c6f798SBarry Smith    Not Collective
12422c6f798SBarry Smith 
1254165533cSJose E. Roman    Input Parameters:
126*f6dfbefdSBarry Smith +  kdm - Original `DMSNES`
127*f6dfbefdSBarry Smith -  nkdm - `DMSNES` to receive the data, should have been created with `DMSNESCreate()`
12822c6f798SBarry Smith 
12922c6f798SBarry Smith    Level: developer
13022c6f798SBarry Smith 
131*f6dfbefdSBarry Smith .seealso: `DMSNES`, `DMSNESCreate()`, `DMSNESDestroy()`
13222c6f798SBarry Smith @*/
1339371c9d4SSatish Balay PetscErrorCode DMSNESCopy(DMSNES kdm, DMSNES nkdm) {
13422c6f798SBarry Smith   PetscFunctionBegin;
13522c6f798SBarry Smith   PetscValidHeaderSpecific(kdm, DMSNES_CLASSID, 1);
13622c6f798SBarry Smith   PetscValidHeaderSpecific(nkdm, DMSNES_CLASSID, 2);
13722c6f798SBarry Smith   nkdm->ops->computefunction  = kdm->ops->computefunction;
1382bc4d0c4SPeter Brune   nkdm->ops->computejacobian  = kdm->ops->computejacobian;
13922c6f798SBarry Smith   nkdm->ops->computegs        = kdm->ops->computegs;
14022c6f798SBarry Smith   nkdm->ops->computeobjective = kdm->ops->computeobjective;
14122c6f798SBarry Smith   nkdm->ops->computepjacobian = kdm->ops->computepjacobian;
14222c6f798SBarry Smith   nkdm->ops->computepfunction = kdm->ops->computepfunction;
14322c6f798SBarry Smith   nkdm->ops->destroy          = kdm->ops->destroy;
14422c6f798SBarry Smith   nkdm->ops->duplicate        = kdm->ops->duplicate;
14522c6f798SBarry Smith 
14622c6f798SBarry Smith   nkdm->gsctx                = kdm->gsctx;
14722c6f798SBarry Smith   nkdm->pctx                 = kdm->pctx;
14822c6f798SBarry Smith   nkdm->objectivectx         = kdm->objectivectx;
149af903c1dSJunchao Zhang   nkdm->originaldm           = kdm->originaldm;
150800f99ffSJeremy L Thompson   nkdm->functionctxcontainer = kdm->functionctxcontainer;
151800f99ffSJeremy L Thompson   nkdm->jacobianctxcontainer = kdm->jacobianctxcontainer;
152800f99ffSJeremy L Thompson   if (nkdm->functionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "function ctx", (PetscObject)nkdm->functionctxcontainer));
153800f99ffSJeremy L Thompson   if (nkdm->jacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "jacobian ctx", (PetscObject)nkdm->jacobianctxcontainer));
15422c6f798SBarry Smith 
15522c6f798SBarry Smith   /*
15622c6f798SBarry Smith   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
15722c6f798SBarry Smith   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
15822c6f798SBarry Smith   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
15922c6f798SBarry Smith   */
16022c6f798SBarry Smith 
16122c6f798SBarry Smith   /* implementation specific copy hooks */
162dbbe0bcdSBarry Smith   PetscTryTypeMethod(kdm, duplicate, nkdm);
16322c6f798SBarry Smith   PetscFunctionReturn(0);
16422c6f798SBarry Smith }
16522c6f798SBarry Smith 
1666cab3a1bSJed Brown /*@C
167*f6dfbefdSBarry Smith    DMGetDMSNES - get read-only private `DMSNES` context from a `DM`
1686cab3a1bSJed Brown 
1696cab3a1bSJed Brown    Not Collective
1706cab3a1bSJed Brown 
1714165533cSJose E. Roman    Input Parameter:
172*f6dfbefdSBarry Smith .  dm - `DM` to be used with `SNES`
1736cab3a1bSJed Brown 
1744165533cSJose E. Roman    Output Parameter:
175*f6dfbefdSBarry Smith .  snesdm - private `DMSNES` context
1766cab3a1bSJed Brown 
1776cab3a1bSJed Brown    Level: developer
1786cab3a1bSJed Brown 
179*f6dfbefdSBarry Smith    Note:
180*f6dfbefdSBarry Smith    Use `DMGetDMSNESWrite()` if write access is needed. The DMSNESSetXXX API should be used wherever possible.
1816cab3a1bSJed Brown 
182db781477SPatrick Sanan .seealso: `DMGetDMSNESWrite()`
1836cab3a1bSJed Brown @*/
1849371c9d4SSatish Balay PetscErrorCode DMGetDMSNES(DM dm, DMSNES *snesdm) {
1856cab3a1bSJed Brown   PetscFunctionBegin;
1866cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
187b4615a05SBarry Smith   *snesdm = (DMSNES)dm->dmsnes;
18822c6f798SBarry Smith   if (!*snesdm) {
1899566063dSJacob Faibussowitsch     PetscCall(PetscInfo(dm, "Creating new DMSNES\n"));
1909566063dSJacob Faibussowitsch     PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dm), snesdm));
1911aa26658SKarl Rupp 
192b4615a05SBarry Smith     dm->dmsnes            = (PetscObject)*snesdm;
193af903c1dSJunchao Zhang     (*snesdm)->originaldm = dm;
1949566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_DMSNES, DMRestrictHook_DMSNES, NULL));
1959566063dSJacob Faibussowitsch     PetscCall(DMRefineHookAdd(dm, DMRefineHook_DMSNES, DMInterpolateHook_DMSNES, NULL));
1969566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_DMSNES, DMSubDomainRestrictHook_DMSNES, NULL));
1976cab3a1bSJed Brown   }
1986cab3a1bSJed Brown   PetscFunctionReturn(0);
1996cab3a1bSJed Brown }
2006cab3a1bSJed Brown 
2016cab3a1bSJed Brown /*@C
202*f6dfbefdSBarry Smith    DMGetDMSNESWrite - get write access to private `DMSNES` context from a `DM`
2036cab3a1bSJed Brown 
2046cab3a1bSJed Brown    Not Collective
2056cab3a1bSJed Brown 
2064165533cSJose E. Roman    Input Parameter:
207*f6dfbefdSBarry Smith .  dm - `DM` to be used with `SNES`
2086cab3a1bSJed Brown 
2094165533cSJose E. Roman    Output Parameter:
210*f6dfbefdSBarry Smith .  snesdm - private `DMSNES` context
2116cab3a1bSJed Brown 
2126cab3a1bSJed Brown    Level: developer
2136cab3a1bSJed Brown 
214db781477SPatrick Sanan .seealso: `DMGetDMSNES()`
2156cab3a1bSJed Brown @*/
2169371c9d4SSatish Balay PetscErrorCode DMGetDMSNESWrite(DM dm, DMSNES *snesdm) {
217942e3340SBarry Smith   DMSNES sdm;
2186cab3a1bSJed Brown 
2196cab3a1bSJed Brown   PetscFunctionBegin;
2206cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2219566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
22228b400f6SJacob Faibussowitsch   PetscCheck(sdm->originaldm, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DMSNES has a NULL originaldm");
2236cab3a1bSJed Brown   if (sdm->originaldm != dm) { /* Copy on write */
224b4615a05SBarry Smith     DMSNES oldsdm = sdm;
2259566063dSJacob Faibussowitsch     PetscCall(PetscInfo(dm, "Copying DMSNES due to write\n"));
2269566063dSJacob Faibussowitsch     PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dm), &sdm));
2279566063dSJacob Faibussowitsch     PetscCall(DMSNESCopy(oldsdm, sdm));
2289566063dSJacob Faibussowitsch     PetscCall(DMSNESDestroy((DMSNES *)&dm->dmsnes));
229b4615a05SBarry Smith     dm->dmsnes      = (PetscObject)sdm;
230af903c1dSJunchao Zhang     sdm->originaldm = dm;
2316cab3a1bSJed Brown   }
2326cab3a1bSJed Brown   *snesdm = sdm;
2336cab3a1bSJed Brown   PetscFunctionReturn(0);
2346cab3a1bSJed Brown }
2356cab3a1bSJed Brown 
2366cab3a1bSJed Brown /*@C
237*f6dfbefdSBarry Smith    DMCopyDMSNES - copies a `DMSNES` context to a new `DM`
2386cab3a1bSJed Brown 
2396cab3a1bSJed Brown    Logically Collective
2406cab3a1bSJed Brown 
2414165533cSJose E. Roman    Input Parameters:
242*f6dfbefdSBarry Smith +  dmsrc - `DM` to obtain context from
243*f6dfbefdSBarry Smith -  dmdest - `DM` to add context to
2446cab3a1bSJed Brown 
2456cab3a1bSJed Brown    Level: developer
2466cab3a1bSJed Brown 
2476cab3a1bSJed Brown    Note:
2486cab3a1bSJed Brown    The context is copied by reference. This function does not ensure that a context exists.
2496cab3a1bSJed Brown 
250db781477SPatrick Sanan .seealso: `DMGetDMSNES()`, `SNESSetDM()`
2516cab3a1bSJed Brown @*/
2529371c9d4SSatish Balay PetscErrorCode DMCopyDMSNES(DM dmsrc, DM dmdest) {
2536cab3a1bSJed Brown   PetscFunctionBegin;
2546cab3a1bSJed Brown   PetscValidHeaderSpecific(dmsrc, DM_CLASSID, 1);
2556cab3a1bSJed Brown   PetscValidHeaderSpecific(dmdest, DM_CLASSID, 2);
2569566063dSJacob Faibussowitsch   if (!dmdest->dmsnes) PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dmdest), (DMSNES *)&dmdest->dmsnes));
2579566063dSJacob Faibussowitsch   PetscCall(DMSNESCopy((DMSNES)dmsrc->dmsnes, (DMSNES)dmdest->dmsnes));
2589566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dmdest, DMCoarsenHook_DMSNES, NULL, NULL));
2599566063dSJacob Faibussowitsch   PetscCall(DMRefineHookAdd(dmdest, DMRefineHook_DMSNES, NULL, NULL));
2609566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dmdest, DMSubDomainHook_DMSNES, DMSubDomainRestrictHook_DMSNES, NULL));
2616cab3a1bSJed Brown   PetscFunctionReturn(0);
2626cab3a1bSJed Brown }
2636cab3a1bSJed Brown 
2646cab3a1bSJed Brown /*@C
265*f6dfbefdSBarry Smith    DMSNESSetFunction - set `SNES` residual evaluation function
2666cab3a1bSJed Brown 
2676cab3a1bSJed Brown    Not Collective
2686cab3a1bSJed Brown 
2694165533cSJose E. Roman    Input Parameters:
270*f6dfbefdSBarry Smith +  dm - DM to be used with `SNES`
271*f6dfbefdSBarry Smith .  f - residual evaluation function; see `SNESFunction` for details
2726cab3a1bSJed Brown -  ctx - context for residual evaluation
2736cab3a1bSJed Brown 
2746cab3a1bSJed Brown    Level: advanced
2756cab3a1bSJed Brown 
2766cab3a1bSJed Brown    Note:
277*f6dfbefdSBarry Smith    `SNESSetFunction()` is normally used, but it calls this function internally because the user context is actually
278*f6dfbefdSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
279*f6dfbefdSBarry Smith    not.
280*f6dfbefdSBarry Smith 
281*f6dfbefdSBarry Smith    Developer Note:
282*f6dfbefdSBarry Smith    If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
2836cab3a1bSJed Brown 
284db781477SPatrick Sanan .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESFunction`
2856cab3a1bSJed Brown @*/
2869371c9d4SSatish Balay PetscErrorCode DMSNESSetFunction(DM dm, PetscErrorCode (*f)(SNES, Vec, Vec, void *), void *ctx) {
287942e3340SBarry Smith   DMSNES sdm;
2886cab3a1bSJed Brown 
2896cab3a1bSJed Brown   PetscFunctionBegin;
2906cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2919566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
292f8b49ee9SBarry Smith   if (f) sdm->ops->computefunction = f;
293800f99ffSJeremy L Thompson   if (ctx) {
294800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
295800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)sdm), &ctxcontainer));
296800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
297800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)sdm, "function ctx", (PetscObject)ctxcontainer));
298800f99ffSJeremy L Thompson     sdm->functionctxcontainer = ctxcontainer;
299800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
300800f99ffSJeremy L Thompson   }
301800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
302800f99ffSJeremy L Thompson }
303800f99ffSJeremy L Thompson 
304800f99ffSJeremy L Thompson /*@C
3055cb80ecdSBarry Smith    DMSNESSetFunctionContextDestroy - set `SNES` residual evaluation context destroy function
306800f99ffSJeremy L Thompson 
307800f99ffSJeremy L Thompson    Not Collective
308800f99ffSJeremy L Thompson 
309800f99ffSJeremy L Thompson    Input Parameters:
3105cb80ecdSBarry Smith +  dm - `DM` to be used with `SNES`
3115cb80ecdSBarry Smith -  f - residual evaluation context destroy function
312800f99ffSJeremy L Thompson 
313800f99ffSJeremy L Thompson    Level: advanced
314800f99ffSJeremy L Thompson 
315800f99ffSJeremy L Thompson .seealso: `DMSNESSetFunction()`, `SNESSetFunction()`
316800f99ffSJeremy L Thompson @*/
3179371c9d4SSatish Balay PetscErrorCode DMSNESSetFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *)) {
318800f99ffSJeremy L Thompson   DMSNES sdm;
319800f99ffSJeremy L Thompson 
320800f99ffSJeremy L Thompson   PetscFunctionBegin;
321800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
322800f99ffSJeremy L Thompson   PetscCall(DMGetDMSNESWrite(dm, &sdm));
323800f99ffSJeremy L Thompson   if (sdm->functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(sdm->functionctxcontainer, f));
324800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
325800f99ffSJeremy L Thompson }
326800f99ffSJeremy L Thompson 
3279371c9d4SSatish Balay PetscErrorCode DMSNESUnsetFunctionContext_Internal(DM dm) {
328800f99ffSJeremy L Thompson   DMSNES sdm;
329800f99ffSJeremy L Thompson 
330800f99ffSJeremy L Thompson   PetscFunctionBegin;
331800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
332800f99ffSJeremy L Thompson   PetscCall(DMGetDMSNESWrite(dm, &sdm));
333800f99ffSJeremy L Thompson   PetscCall(DMSNESUnsetFunctionContext_DMSNES(sdm));
3346cab3a1bSJed Brown   PetscFunctionReturn(0);
3356cab3a1bSJed Brown }
3366cab3a1bSJed Brown 
3376cab3a1bSJed Brown /*@C
338*f6dfbefdSBarry Smith    DMSNESSetMFFunction - set `SNES` residual evaluation function used in applying the matrix-free Jacobian with -snes_mf_operator
339bbc1464cSBarry Smith 
340bbc1464cSBarry Smith    Logically Collective on dm
341bbc1464cSBarry Smith 
3424165533cSJose E. Roman    Input Parameters:
343*f6dfbefdSBarry Smith +  dm - `DM` to be used with `SNES`
344*f6dfbefdSBarry Smith -  f - residual evaluation function; see `SNESFunction` for details
345bbc1464cSBarry Smith 
346bbc1464cSBarry Smith    Level: advanced
347bbc1464cSBarry Smith 
348db781477SPatrick Sanan .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESFunction`, `DMSNESSetFunction()`
349bbc1464cSBarry Smith @*/
3509371c9d4SSatish Balay PetscErrorCode DMSNESSetMFFunction(DM dm, PetscErrorCode (*f)(SNES, Vec, Vec, void *), void *ctx) {
351bbc1464cSBarry Smith   DMSNES sdm;
352bbc1464cSBarry Smith 
353bbc1464cSBarry Smith   PetscFunctionBegin;
354bbc1464cSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
35548a46eb9SPierre Jolivet   if (f || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm));
356bbc1464cSBarry Smith   if (f) sdm->ops->computemffunction = f;
357bbc1464cSBarry Smith   if (ctx) sdm->mffunctionctx = ctx;
358bbc1464cSBarry Smith   PetscFunctionReturn(0);
359bbc1464cSBarry Smith }
360bbc1464cSBarry Smith 
361bbc1464cSBarry Smith /*@C
362*f6dfbefdSBarry Smith    DMSNESGetFunction - get `SNES` residual evaluation function
3636cab3a1bSJed Brown 
3646cab3a1bSJed Brown    Not Collective
3656cab3a1bSJed Brown 
3664165533cSJose E. Roman    Input Parameter:
367*f6dfbefdSBarry Smith .  dm - `DM` to be used with `SNES`
3686cab3a1bSJed Brown 
3694165533cSJose E. Roman    Output Parameters:
370*f6dfbefdSBarry Smith +  f - residual evaluation function; see `SNESFunction` for details
3716cab3a1bSJed Brown -  ctx - context for residual evaluation
3726cab3a1bSJed Brown 
3736cab3a1bSJed Brown    Level: advanced
3746cab3a1bSJed Brown 
3756cab3a1bSJed Brown    Note:
376*f6dfbefdSBarry Smith    `SNESGetFunction()` is normally used, but it calls this function internally because the user context is actually
377*f6dfbefdSBarry Smith    associated with the `DM`.
3786cab3a1bSJed Brown 
379db781477SPatrick Sanan .seealso: `DMSNESSetContext()`, `DMSNESSetFunction()`, `SNESSetFunction()`, `SNESFunction`
3806cab3a1bSJed Brown @*/
3819371c9d4SSatish Balay PetscErrorCode DMSNESGetFunction(DM dm, PetscErrorCode (**f)(SNES, Vec, Vec, void *), void **ctx) {
382942e3340SBarry Smith   DMSNES sdm;
3836cab3a1bSJed Brown 
3846cab3a1bSJed Brown   PetscFunctionBegin;
3856cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3869566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
387f8b49ee9SBarry Smith   if (f) *f = sdm->ops->computefunction;
388800f99ffSJeremy L Thompson   if (ctx) {
389800f99ffSJeremy L Thompson     if (sdm->functionctxcontainer) PetscCall(PetscContainerGetPointer(sdm->functionctxcontainer, ctx));
390800f99ffSJeremy L Thompson     else *ctx = NULL;
391800f99ffSJeremy L Thompson   }
3926cab3a1bSJed Brown   PetscFunctionReturn(0);
3936cab3a1bSJed Brown }
3946cab3a1bSJed Brown 
3952a4ee8f2SPeter Brune /*@C
396*f6dfbefdSBarry Smith    DMSNESSetObjective - set `SNES` objective evaluation function
3972a4ee8f2SPeter Brune 
3982a4ee8f2SPeter Brune    Not Collective
3992a4ee8f2SPeter Brune 
4004165533cSJose E. Roman    Input Parameters:
401*f6dfbefdSBarry Smith +  dm - `DM` to be used with `SNES`
402*f6dfbefdSBarry Smith .  obj - objective evaluation function; see `SNESObjectiveFunction` for details
4032a4ee8f2SPeter Brune -  ctx - context for residual evaluation
4042a4ee8f2SPeter Brune 
4052a4ee8f2SPeter Brune    Level: advanced
4062a4ee8f2SPeter Brune 
407db781477SPatrick Sanan .seealso: `DMSNESSetContext()`, `SNESGetObjective()`, `DMSNESSetFunction()`
4082a4ee8f2SPeter Brune @*/
4099371c9d4SSatish Balay PetscErrorCode DMSNESSetObjective(DM dm, PetscErrorCode (*obj)(SNES, Vec, PetscReal *, void *), void *ctx) {
410942e3340SBarry Smith   DMSNES sdm;
4112a4ee8f2SPeter Brune 
4122a4ee8f2SPeter Brune   PetscFunctionBegin;
4132a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
41448a46eb9SPierre Jolivet   if (obj || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm));
415f8b49ee9SBarry Smith   if (obj) sdm->ops->computeobjective = obj;
4162a4ee8f2SPeter Brune   if (ctx) sdm->objectivectx = ctx;
4172a4ee8f2SPeter Brune   PetscFunctionReturn(0);
4182a4ee8f2SPeter Brune }
4192a4ee8f2SPeter Brune 
4202a4ee8f2SPeter Brune /*@C
421*f6dfbefdSBarry Smith    DMSNESGetObjective - get `SNES` objective evaluation function
4222a4ee8f2SPeter Brune 
4232a4ee8f2SPeter Brune    Not Collective
4242a4ee8f2SPeter Brune 
4254165533cSJose E. Roman    Input Parameter:
426*f6dfbefdSBarry Smith .  dm - `DM` to be used with `SNES`
4272a4ee8f2SPeter Brune 
4284165533cSJose E. Roman    Output Parameters:
429*f6dfbefdSBarry Smith +  obj- residual evaluation function; see `SNESObjectiveFunction` for details
4302a4ee8f2SPeter Brune -  ctx - context for residual evaluation
4312a4ee8f2SPeter Brune 
4322a4ee8f2SPeter Brune    Level: advanced
4332a4ee8f2SPeter Brune 
4342a4ee8f2SPeter Brune    Note:
435*f6dfbefdSBarry Smith    `SNESGetFunction()` is normally used, but it calls this function internally because the user context is actually
436*f6dfbefdSBarry Smith    associated with the `DM`.
4372a4ee8f2SPeter Brune 
438db781477SPatrick Sanan .seealso: `DMSNESSetContext()`, `DMSNESSetObjective()`, `SNESSetFunction()`
4392a4ee8f2SPeter Brune @*/
4409371c9d4SSatish Balay PetscErrorCode DMSNESGetObjective(DM dm, PetscErrorCode (**obj)(SNES, Vec, PetscReal *, void *), void **ctx) {
441942e3340SBarry Smith   DMSNES sdm;
4422a4ee8f2SPeter Brune 
4432a4ee8f2SPeter Brune   PetscFunctionBegin;
4442a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4459566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
446f8b49ee9SBarry Smith   if (obj) *obj = sdm->ops->computeobjective;
4472a4ee8f2SPeter Brune   if (ctx) *ctx = sdm->objectivectx;
4482a4ee8f2SPeter Brune   PetscFunctionReturn(0);
4492a4ee8f2SPeter Brune }
4502a4ee8f2SPeter Brune 
4516cab3a1bSJed Brown /*@C
452*f6dfbefdSBarry Smith    DMSNESSetNGS - set `SNES` Gauss-Seidel relaxation function
4536cab3a1bSJed Brown 
4546cab3a1bSJed Brown    Not Collective
4556cab3a1bSJed Brown 
4564165533cSJose E. Roman    Input Parameters:
457*f6dfbefdSBarry Smith +  dm - `DM` to be used with `SNES`
458*f6dfbefdSBarry Smith .  f  - relaxation function, see `SNESGSFunction`
4596cab3a1bSJed Brown -  ctx - context for residual evaluation
4606cab3a1bSJed Brown 
4616cab3a1bSJed Brown    Level: advanced
4626cab3a1bSJed Brown 
4636cab3a1bSJed Brown    Note:
464*f6dfbefdSBarry Smith    `SNESSetNGS()` is normally used, but it calls this function internally because the user context is actually
465*f6dfbefdSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
466*f6dfbefdSBarry Smith    not.
467*f6dfbefdSBarry Smith 
468*f6dfbefdSBarry Smith    Dveloper Note:
469*f6dfbefdSBarry Smith    If `DM `took a more central role at some later date, this could become the primary method of supplying the smoother
4706cab3a1bSJed Brown 
471db781477SPatrick Sanan .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `DMSNESSetFunction()`, `SNESGSFunction`
4726cab3a1bSJed Brown @*/
4739371c9d4SSatish Balay PetscErrorCode DMSNESSetNGS(DM dm, PetscErrorCode (*f)(SNES, Vec, Vec, void *), void *ctx) {
474942e3340SBarry Smith   DMSNES sdm;
4756cab3a1bSJed Brown 
4766cab3a1bSJed Brown   PetscFunctionBegin;
4776cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
47848a46eb9SPierre Jolivet   if (f || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm));
479be95d8f1SBarry Smith   if (f) sdm->ops->computegs = f;
4806cab3a1bSJed Brown   if (ctx) sdm->gsctx = ctx;
4816cab3a1bSJed Brown   PetscFunctionReturn(0);
4826cab3a1bSJed Brown }
4836cab3a1bSJed Brown 
4846cab3a1bSJed Brown /*@C
485*f6dfbefdSBarry Smith    DMSNESGetNGS - get `SNES` Gauss-Seidel relaxation function
4866cab3a1bSJed Brown 
4876cab3a1bSJed Brown    Not Collective
4886cab3a1bSJed Brown 
4894165533cSJose E. Roman    Input Parameter:
490*f6dfbefdSBarry Smith .  dm - `DM` to be used with `SNES`
4916cab3a1bSJed Brown 
4924165533cSJose E. Roman    Output Parameters:
493*f6dfbefdSBarry Smith +  f - relaxation function which performs Gauss-Seidel sweeps, see `SNESGSFunction`
4946cab3a1bSJed Brown -  ctx - context for residual evaluation
4956cab3a1bSJed Brown 
4966cab3a1bSJed Brown    Level: advanced
4976cab3a1bSJed Brown 
4986cab3a1bSJed Brown    Note:
499*f6dfbefdSBarry Smith    `SNESGetNGS()` is normally used, but it calls this function internally because the user context is actually
500*f6dfbefdSBarry Smith    associated with the `DM`.
501*f6dfbefdSBarry Smith 
502*f6dfbefdSBarry Smith    Developer Note:
503*f6dfbefdSBarry Smith    This makes the interface consistent regardless of whether the user interacts with a `DM` or
504*f6dfbefdSBarry Smith    not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
5056cab3a1bSJed Brown 
506db781477SPatrick Sanan .seealso: `DMSNESSetContext()`, `SNESGetNGS()`, `DMSNESGetJacobian()`, `DMSNESGetFunction()`, `SNESNGSFunction`
5076cab3a1bSJed Brown @*/
5089371c9d4SSatish Balay PetscErrorCode DMSNESGetNGS(DM dm, PetscErrorCode (**f)(SNES, Vec, Vec, void *), void **ctx) {
509942e3340SBarry Smith   DMSNES sdm;
5106cab3a1bSJed Brown 
5116cab3a1bSJed Brown   PetscFunctionBegin;
5126cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5139566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
514be95d8f1SBarry Smith   if (f) *f = sdm->ops->computegs;
5156cab3a1bSJed Brown   if (ctx) *ctx = sdm->gsctx;
5166cab3a1bSJed Brown   PetscFunctionReturn(0);
5176cab3a1bSJed Brown }
5186cab3a1bSJed Brown 
5196cab3a1bSJed Brown /*@C
520*f6dfbefdSBarry Smith    DMSNESSetJacobian - set `SNES` Jacobian evaluation function
5216cab3a1bSJed Brown 
5226cab3a1bSJed Brown    Not Collective
5236cab3a1bSJed Brown 
5244165533cSJose E. Roman    Input Parameters:
525*f6dfbefdSBarry Smith +  dm - `DM` to be used with `SNES`
526f8b49ee9SBarry Smith .  J - Jacobian evaluation function
5276cab3a1bSJed Brown -  ctx - context for residual evaluation
5286cab3a1bSJed Brown 
5296cab3a1bSJed Brown    Level: advanced
5306cab3a1bSJed Brown 
5316cab3a1bSJed Brown    Note:
532*f6dfbefdSBarry Smith    `SNESSetJacobian()` is normally used, but it calls this function internally because the user context is actually
533*f6dfbefdSBarry Smith    associated with the `DM`.
534*f6dfbefdSBarry Smith 
535*f6dfbefdSBarry Smith    Developer Note:
536*f6dfbefdSBarry Smith    This makes the interface consistent regardless of whether the user interacts with a `DM` or
537*f6dfbefdSBarry Smith    not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
5386cab3a1bSJed Brown 
539db781477SPatrick Sanan .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESGetJacobian()`, `SNESSetJacobian()`, `SNESJacobianFunction`
5406cab3a1bSJed Brown @*/
5419371c9d4SSatish Balay PetscErrorCode DMSNESSetJacobian(DM dm, PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *), void *ctx) {
542942e3340SBarry Smith   DMSNES sdm;
5436cab3a1bSJed Brown 
5446cab3a1bSJed Brown   PetscFunctionBegin;
5456cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
546800f99ffSJeremy L Thompson   if (J || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm));
547f8b49ee9SBarry Smith   if (J) sdm->ops->computejacobian = J;
548800f99ffSJeremy L Thompson   if (ctx) {
549800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
550800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)sdm), &ctxcontainer));
551800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
552800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)sdm, "jacobian ctx", (PetscObject)ctxcontainer));
553800f99ffSJeremy L Thompson     sdm->jacobianctxcontainer = ctxcontainer;
554800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
555800f99ffSJeremy L Thompson   }
556800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
557800f99ffSJeremy L Thompson }
558800f99ffSJeremy L Thompson 
559800f99ffSJeremy L Thompson /*@C
5605cb80ecdSBarry Smith    DMSNESSetJacobianContextDestroy - set `SNES` Jacobian evaluation context destroy function
561800f99ffSJeremy L Thompson 
562800f99ffSJeremy L Thompson    Not Collective
563800f99ffSJeremy L Thompson 
564800f99ffSJeremy L Thompson    Input Parameters:
5655cb80ecdSBarry Smith +  dm - `DM` to be used with `SNES`
5665cb80ecdSBarry Smith -  f - Jacobian evaluation contex destroy function
567800f99ffSJeremy L Thompson 
568800f99ffSJeremy L Thompson    Level: advanced
569800f99ffSJeremy L Thompson 
57090f54644SBarry Smith .seealso: `DMSNESSetJacobian()`
571800f99ffSJeremy L Thompson @*/
5729371c9d4SSatish Balay PetscErrorCode DMSNESSetJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *)) {
573800f99ffSJeremy L Thompson   DMSNES sdm;
574800f99ffSJeremy L Thompson 
575800f99ffSJeremy L Thompson   PetscFunctionBegin;
576800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
577800f99ffSJeremy L Thompson   PetscCall(DMGetDMSNESWrite(dm, &sdm));
578800f99ffSJeremy L Thompson   if (sdm->jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(sdm->jacobianctxcontainer, f));
579800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
580800f99ffSJeremy L Thompson }
581800f99ffSJeremy L Thompson 
5829371c9d4SSatish Balay PetscErrorCode DMSNESUnsetJacobianContext_Internal(DM dm) {
583800f99ffSJeremy L Thompson   DMSNES sdm;
584800f99ffSJeremy L Thompson 
585800f99ffSJeremy L Thompson   PetscFunctionBegin;
586800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
587800f99ffSJeremy L Thompson   PetscCall(DMGetDMSNESWrite(dm, &sdm));
588800f99ffSJeremy L Thompson   PetscCall(DMSNESUnsetJacobianContext_DMSNES(sdm));
5896cab3a1bSJed Brown   PetscFunctionReturn(0);
5906cab3a1bSJed Brown }
5916cab3a1bSJed Brown 
5926cab3a1bSJed Brown /*@C
593*f6dfbefdSBarry Smith    DMSNESGetJacobian - get `SNES` Jacobian evaluation function
5946cab3a1bSJed Brown 
5956cab3a1bSJed Brown    Not Collective
5966cab3a1bSJed Brown 
5974165533cSJose E. Roman    Input Parameter:
598*f6dfbefdSBarry Smith .  dm - `DM` to be used with `SNES`
5996cab3a1bSJed Brown 
6004165533cSJose E. Roman    Output Parameters:
601*f6dfbefdSBarry Smith +  J - Jacobian evaluation function; see `SNESJacobianFunction` for all calling sequence
6026cab3a1bSJed Brown -  ctx - context for residual evaluation
6036cab3a1bSJed Brown 
6046cab3a1bSJed Brown    Level: advanced
6056cab3a1bSJed Brown 
6066cab3a1bSJed Brown    Note:
607*f6dfbefdSBarry Smith    `SNESGetJacobian()` is normally used, but it calls this function internally because the user context is actually
608*f6dfbefdSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
609*f6dfbefdSBarry Smith    not.
610*f6dfbefdSBarry Smith 
611*f6dfbefdSBarry Smith    Developer Note:
612*f6dfbefdSBarry Smith    If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
6136cab3a1bSJed Brown 
614db781477SPatrick Sanan .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESJacobianFunction`
6156cab3a1bSJed Brown @*/
6169371c9d4SSatish Balay PetscErrorCode DMSNESGetJacobian(DM dm, PetscErrorCode (**J)(SNES, Vec, Mat, Mat, void *), void **ctx) {
617942e3340SBarry Smith   DMSNES sdm;
6186cab3a1bSJed Brown 
6196cab3a1bSJed Brown   PetscFunctionBegin;
6206cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6219566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
622f8b49ee9SBarry Smith   if (J) *J = sdm->ops->computejacobian;
623800f99ffSJeremy L Thompson   if (ctx) {
624800f99ffSJeremy L Thompson     if (sdm->jacobianctxcontainer) PetscCall(PetscContainerGetPointer(sdm->jacobianctxcontainer, ctx));
625800f99ffSJeremy L Thompson     else *ctx = NULL;
626800f99ffSJeremy L Thompson   }
6276cab3a1bSJed Brown   PetscFunctionReturn(0);
6286cab3a1bSJed Brown }
6296cab3a1bSJed Brown 
630e03ab78fSPeter Brune /*@C
631e03ab78fSPeter Brune    DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.
632e03ab78fSPeter Brune 
633e03ab78fSPeter Brune    Not Collective
634e03ab78fSPeter Brune 
6354165533cSJose E. Roman    Input Parameters:
636*f6dfbefdSBarry Smith +  dm - `DM` to be used with `SNES`
637f8b49ee9SBarry Smith .  b - RHS evaluation function
638f8b49ee9SBarry Smith .  J - Picard matrix evaluation function
639e03ab78fSPeter Brune -  ctx - context for residual evaluation
640e03ab78fSPeter Brune 
641e03ab78fSPeter Brune    Level: advanced
642e03ab78fSPeter Brune 
643db781477SPatrick Sanan .seealso: `SNESSetPicard()`, `DMSNESSetFunction()`, `DMSNESSetJacobian()`
644e03ab78fSPeter Brune @*/
6459371c9d4SSatish Balay PetscErrorCode DMSNESSetPicard(DM dm, PetscErrorCode (*b)(SNES, Vec, Vec, void *), PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *), void *ctx) {
646942e3340SBarry Smith   DMSNES sdm;
647e03ab78fSPeter Brune 
648e03ab78fSPeter Brune   PetscFunctionBegin;
649e03ab78fSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6509566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
651f8b49ee9SBarry Smith   if (b) sdm->ops->computepfunction = b;
652f8b49ee9SBarry Smith   if (J) sdm->ops->computepjacobian = J;
653e03ab78fSPeter Brune   if (ctx) sdm->pctx = ctx;
654e03ab78fSPeter Brune   PetscFunctionReturn(0);
655e03ab78fSPeter Brune }
656e03ab78fSPeter Brune 
6577971a8bfSPeter Brune /*@C
658*f6dfbefdSBarry Smith    DMSNESGetPicard - get `SNES` Picard iteration evaluation functions
6597971a8bfSPeter Brune 
6607971a8bfSPeter Brune    Not Collective
6617971a8bfSPeter Brune 
6624165533cSJose E. Roman    Input Parameter:
663*f6dfbefdSBarry Smith .  dm - `DM` to be used with `SNES`
6647971a8bfSPeter Brune 
6654165533cSJose E. Roman    Output Parameters:
666*f6dfbefdSBarry Smith +  b - RHS evaluation function; see `SNESFunction` for details
667*f6dfbefdSBarry Smith .  J  - RHS evaluation function; see `SNESJacobianFunction` for detailsa
6687971a8bfSPeter Brune -  ctx - context for residual evaluation
6697971a8bfSPeter Brune 
6707971a8bfSPeter Brune    Level: advanced
6717971a8bfSPeter Brune 
672db781477SPatrick Sanan .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`
6737971a8bfSPeter Brune @*/
6749371c9d4SSatish Balay PetscErrorCode DMSNESGetPicard(DM dm, PetscErrorCode (**b)(SNES, Vec, Vec, void *), PetscErrorCode (**J)(SNES, Vec, Mat, Mat, void *), void **ctx) {
675942e3340SBarry Smith   DMSNES sdm;
6767971a8bfSPeter Brune 
6777971a8bfSPeter Brune   PetscFunctionBegin;
6787971a8bfSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6799566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
680f8b49ee9SBarry Smith   if (b) *b = sdm->ops->computepfunction;
681f8b49ee9SBarry Smith   if (J) *J = sdm->ops->computepjacobian;
6827971a8bfSPeter Brune   if (ctx) *ctx = sdm->pctx;
6837971a8bfSPeter Brune   PetscFunctionReturn(0);
6847971a8bfSPeter Brune }
685