xref: /petsc/src/snes/utils/dmsnes.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2 #include <petsc/private/dmimpl.h>   /*I "petscdm.h" I*/
3 
4 static PetscErrorCode DMSNESUnsetFunctionContext_DMSNES(DMSNES sdm) {
5   PetscFunctionBegin;
6   PetscCall(PetscObjectCompose((PetscObject)sdm, "function ctx", NULL));
7   sdm->functionctxcontainer = NULL;
8   PetscFunctionReturn(0);
9 }
10 
11 static PetscErrorCode DMSNESUnsetJacobianContext_DMSNES(DMSNES sdm) {
12   PetscFunctionBegin;
13   PetscCall(PetscObjectCompose((PetscObject)sdm, "jacobian ctx", NULL));
14   sdm->jacobianctxcontainer = NULL;
15   PetscFunctionReturn(0);
16 }
17 
18 static PetscErrorCode DMSNESDestroy(DMSNES *kdm) {
19   PetscFunctionBegin;
20   if (!*kdm) PetscFunctionReturn(0);
21   PetscValidHeaderSpecific((*kdm), DMSNES_CLASSID, 1);
22   if (--((PetscObject)(*kdm))->refct > 0) {
23     *kdm = NULL;
24     PetscFunctionReturn(0);
25   }
26   PetscCall(DMSNESUnsetFunctionContext_DMSNES(*kdm));
27   PetscCall(DMSNESUnsetJacobianContext_DMSNES(*kdm));
28   if ((*kdm)->ops->destroy) PetscCall(((*kdm)->ops->destroy)(*kdm));
29   PetscCall(PetscHeaderDestroy(kdm));
30   PetscFunctionReturn(0);
31 }
32 
33 PetscErrorCode DMSNESLoad(DMSNES kdm, PetscViewer viewer) {
34   PetscFunctionBegin;
35   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->computefunction, 1, NULL, PETSC_FUNCTION));
36   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->computejacobian, 1, NULL, PETSC_FUNCTION));
37   PetscFunctionReturn(0);
38 }
39 
40 PetscErrorCode DMSNESView(DMSNES kdm, PetscViewer viewer) {
41   PetscBool isascii, isbinary;
42 
43   PetscFunctionBegin;
44   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
45   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
46   if (isascii) {
47 #if defined(PETSC_SERIALIZE_FUNCTIONS) && defined(PETSC_SERIALIZE_FUNCTIONS_VIEW)
48     const char *fname;
49 
50     PetscCall(PetscFPTFind(kdm->ops->computefunction, &fname));
51     if (fname) { PetscCall(PetscViewerASCIIPrintf(viewer, "Function used by SNES: %s\n", fname)); }
52     PetscCall(PetscFPTFind(kdm->ops->computejacobian, &fname));
53     if (fname) { PetscCall(PetscViewerASCIIPrintf(viewer, "Jacobian function used by SNES: %s\n", fname)); }
54 #endif
55   } else if (isbinary) {
56     struct {
57       PetscErrorCode (*func)(SNES, Vec, Vec, void *);
58     } funcstruct;
59     struct {
60       PetscErrorCode (*jac)(SNES, Vec, Mat, Mat, void *);
61     } jacstruct;
62     funcstruct.func = kdm->ops->computefunction;
63     jacstruct.jac   = kdm->ops->computejacobian;
64     PetscCall(PetscViewerBinaryWrite(viewer, &funcstruct, 1, PETSC_FUNCTION));
65     PetscCall(PetscViewerBinaryWrite(viewer, &jacstruct, 1, PETSC_FUNCTION));
66   }
67   PetscFunctionReturn(0);
68 }
69 
70 static PetscErrorCode DMSNESCreate(MPI_Comm comm, DMSNES *kdm) {
71   PetscFunctionBegin;
72   PetscCall(SNESInitializePackage());
73   PetscCall(PetscHeaderCreate(*kdm, DMSNES_CLASSID, "DMSNES", "DMSNES", "DMSNES", comm, DMSNESDestroy, DMSNESView));
74   PetscFunctionReturn(0);
75 }
76 
77 /* Attaches the DMSNES to the coarse level.
78  * Under what conditions should we copy versus duplicate?
79  */
80 static PetscErrorCode DMCoarsenHook_DMSNES(DM dm, DM dmc, void *ctx) {
81   PetscFunctionBegin;
82   PetscCall(DMCopyDMSNES(dm, dmc));
83   PetscFunctionReturn(0);
84 }
85 
86 /* This could restrict auxiliary information to the coarse level.
87  */
88 static PetscErrorCode DMRestrictHook_DMSNES(DM dm, Mat Restrict, Vec rscale, Mat Inject, DM dmc, void *ctx) {
89   PetscFunctionBegin;
90   PetscFunctionReturn(0);
91 }
92 
93 /* Attaches the DMSNES to the subdomain. */
94 static PetscErrorCode DMSubDomainHook_DMSNES(DM dm, DM subdm, void *ctx) {
95   PetscFunctionBegin;
96   PetscCall(DMCopyDMSNES(dm, subdm));
97   PetscFunctionReturn(0);
98 }
99 
100 /* This could restrict auxiliary information to the coarse level.
101  */
102 static PetscErrorCode DMSubDomainRestrictHook_DMSNES(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) {
103   PetscFunctionBegin;
104   PetscFunctionReturn(0);
105 }
106 
107 static PetscErrorCode DMRefineHook_DMSNES(DM dm, DM dmf, void *ctx) {
108   PetscFunctionBegin;
109   PetscCall(DMCopyDMSNES(dm, dmf));
110   PetscFunctionReturn(0);
111 }
112 
113 /* This could restrict auxiliary information to the coarse level.
114  */
115 static PetscErrorCode DMInterpolateHook_DMSNES(DM dm, Mat Interp, DM dmf, void *ctx) {
116   PetscFunctionBegin;
117   PetscFunctionReturn(0);
118 }
119 
120 /*@C
121    DMSNESCopy - copies the information in a DMSNES to another DMSNES
122 
123    Not Collective
124 
125    Input Parameters:
126 +  kdm - Original DMSNES
127 -  nkdm - DMSNES to receive the data, should have been created with DMSNESCreate()
128 
129    Level: developer
130 
131 .seealso: `DMSNESCreate()`, `DMSNESDestroy()`
132 @*/
133 PetscErrorCode DMSNESCopy(DMSNES kdm, DMSNES nkdm) {
134   PetscFunctionBegin;
135   PetscValidHeaderSpecific(kdm, DMSNES_CLASSID, 1);
136   PetscValidHeaderSpecific(nkdm, DMSNES_CLASSID, 2);
137   nkdm->ops->computefunction  = kdm->ops->computefunction;
138   nkdm->ops->computejacobian  = kdm->ops->computejacobian;
139   nkdm->ops->computegs        = kdm->ops->computegs;
140   nkdm->ops->computeobjective = kdm->ops->computeobjective;
141   nkdm->ops->computepjacobian = kdm->ops->computepjacobian;
142   nkdm->ops->computepfunction = kdm->ops->computepfunction;
143   nkdm->ops->destroy          = kdm->ops->destroy;
144   nkdm->ops->duplicate        = kdm->ops->duplicate;
145 
146   nkdm->gsctx                = kdm->gsctx;
147   nkdm->pctx                 = kdm->pctx;
148   nkdm->objectivectx         = kdm->objectivectx;
149   nkdm->originaldm           = kdm->originaldm;
150   nkdm->functionctxcontainer = kdm->functionctxcontainer;
151   nkdm->jacobianctxcontainer = kdm->jacobianctxcontainer;
152   if (nkdm->functionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "function ctx", (PetscObject)nkdm->functionctxcontainer));
153   if (nkdm->jacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "jacobian ctx", (PetscObject)nkdm->jacobianctxcontainer));
154 
155   /*
156   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
157   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
158   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
159   */
160 
161   /* implementation specific copy hooks */
162   PetscTryTypeMethod(kdm, duplicate, nkdm);
163   PetscFunctionReturn(0);
164 }
165 
166 /*@C
167    DMGetDMSNES - get read-only private DMSNES context from a DM
168 
169    Not Collective
170 
171    Input Parameter:
172 .  dm - DM to be used with SNES
173 
174    Output Parameter:
175 .  snesdm - private DMSNES context
176 
177    Level: developer
178 
179    Notes:
180    Use DMGetDMSNESWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible.
181 
182 .seealso: `DMGetDMSNESWrite()`
183 @*/
184 PetscErrorCode DMGetDMSNES(DM dm, DMSNES *snesdm) {
185   PetscFunctionBegin;
186   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
187   *snesdm = (DMSNES)dm->dmsnes;
188   if (!*snesdm) {
189     PetscCall(PetscInfo(dm, "Creating new DMSNES\n"));
190     PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dm), snesdm));
191 
192     dm->dmsnes            = (PetscObject)*snesdm;
193     (*snesdm)->originaldm = dm;
194     PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_DMSNES, DMRestrictHook_DMSNES, NULL));
195     PetscCall(DMRefineHookAdd(dm, DMRefineHook_DMSNES, DMInterpolateHook_DMSNES, NULL));
196     PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_DMSNES, DMSubDomainRestrictHook_DMSNES, NULL));
197   }
198   PetscFunctionReturn(0);
199 }
200 
201 /*@C
202    DMGetDMSNESWrite - get write access to private DMSNES context from a DM
203 
204    Not Collective
205 
206    Input Parameter:
207 .  dm - DM to be used with SNES
208 
209    Output Parameter:
210 .  snesdm - private DMSNES context
211 
212    Level: developer
213 
214 .seealso: `DMGetDMSNES()`
215 @*/
216 PetscErrorCode DMGetDMSNESWrite(DM dm, DMSNES *snesdm) {
217   DMSNES sdm;
218 
219   PetscFunctionBegin;
220   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
221   PetscCall(DMGetDMSNES(dm, &sdm));
222   PetscCheck(sdm->originaldm, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DMSNES has a NULL originaldm");
223   if (sdm->originaldm != dm) { /* Copy on write */
224     DMSNES oldsdm = sdm;
225     PetscCall(PetscInfo(dm, "Copying DMSNES due to write\n"));
226     PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dm), &sdm));
227     PetscCall(DMSNESCopy(oldsdm, sdm));
228     PetscCall(DMSNESDestroy((DMSNES *)&dm->dmsnes));
229     dm->dmsnes      = (PetscObject)sdm;
230     sdm->originaldm = dm;
231   }
232   *snesdm = sdm;
233   PetscFunctionReturn(0);
234 }
235 
236 /*@C
237    DMCopyDMSNES - copies a DM context to a new DM
238 
239    Logically Collective
240 
241    Input Parameters:
242 +  dmsrc - DM to obtain context from
243 -  dmdest - DM to add context to
244 
245    Level: developer
246 
247    Note:
248    The context is copied by reference. This function does not ensure that a context exists.
249 
250 .seealso: `DMGetDMSNES()`, `SNESSetDM()`
251 @*/
252 PetscErrorCode DMCopyDMSNES(DM dmsrc, DM dmdest) {
253   PetscFunctionBegin;
254   PetscValidHeaderSpecific(dmsrc, DM_CLASSID, 1);
255   PetscValidHeaderSpecific(dmdest, DM_CLASSID, 2);
256   if (!dmdest->dmsnes) PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dmdest), (DMSNES *)&dmdest->dmsnes));
257   PetscCall(DMSNESCopy((DMSNES)dmsrc->dmsnes, (DMSNES)dmdest->dmsnes));
258   PetscCall(DMCoarsenHookAdd(dmdest, DMCoarsenHook_DMSNES, NULL, NULL));
259   PetscCall(DMRefineHookAdd(dmdest, DMRefineHook_DMSNES, NULL, NULL));
260   PetscCall(DMSubDomainHookAdd(dmdest, DMSubDomainHook_DMSNES, DMSubDomainRestrictHook_DMSNES, NULL));
261   PetscFunctionReturn(0);
262 }
263 
264 /*@C
265    DMSNESSetFunction - set SNES residual evaluation function
266 
267    Not Collective
268 
269    Input Parameters:
270 +  dm - DM to be used with SNES
271 .  f - residual evaluation function; see SNESFunction for details
272 -  ctx - context for residual evaluation
273 
274    Level: advanced
275 
276    Note:
277    SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
278    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
279    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
280 
281 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESFunction`
282 @*/
283 PetscErrorCode DMSNESSetFunction(DM dm, PetscErrorCode (*f)(SNES, Vec, Vec, void *), void *ctx) {
284   DMSNES sdm;
285 
286   PetscFunctionBegin;
287   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
288   PetscCall(DMGetDMSNESWrite(dm, &sdm));
289   if (f) sdm->ops->computefunction = f;
290   if (ctx) {
291     PetscContainer ctxcontainer;
292     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)sdm), &ctxcontainer));
293     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
294     PetscCall(PetscObjectCompose((PetscObject)sdm, "function ctx", (PetscObject)ctxcontainer));
295     sdm->functionctxcontainer = ctxcontainer;
296     PetscCall(PetscContainerDestroy(&ctxcontainer));
297   }
298   PetscFunctionReturn(0);
299 }
300 
301 /*@C
302    DMSNESSetFunctionContextDestroy - set `SNES` residual evaluation context destroy function
303 
304    Not Collective
305 
306    Input Parameters:
307 +  dm - `DM` to be used with `SNES`
308 -  f - residual evaluation context destroy function
309 
310    Level: advanced
311 
312 .seealso: `DMSNESSetFunction()`, `SNESSetFunction()`
313 @*/
314 PetscErrorCode DMSNESSetFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *)) {
315   DMSNES sdm;
316 
317   PetscFunctionBegin;
318   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
319   PetscCall(DMGetDMSNESWrite(dm, &sdm));
320   if (sdm->functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(sdm->functionctxcontainer, f));
321   PetscFunctionReturn(0);
322 }
323 
324 PetscErrorCode DMSNESUnsetFunctionContext_Internal(DM dm) {
325   DMSNES sdm;
326 
327   PetscFunctionBegin;
328   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
329   PetscCall(DMGetDMSNESWrite(dm, &sdm));
330   PetscCall(DMSNESUnsetFunctionContext_DMSNES(sdm));
331   PetscFunctionReturn(0);
332 }
333 
334 /*@C
335    DMSNESSetMFFunction - set SNES residual evaluation function used in applying the matrix-free Jacobian with -snes_mf_operator
336 
337    Logically Collective on dm
338 
339    Input Parameters:
340 +  dm - DM to be used with SNES
341 -  f - residual evaluation function; see SNESFunction for details
342 
343    Level: advanced
344 
345 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESFunction`, `DMSNESSetFunction()`
346 @*/
347 PetscErrorCode DMSNESSetMFFunction(DM dm, PetscErrorCode (*f)(SNES, Vec, Vec, void *), void *ctx) {
348   DMSNES sdm;
349 
350   PetscFunctionBegin;
351   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
352   if (f || ctx) { PetscCall(DMGetDMSNESWrite(dm, &sdm)); }
353   if (f) sdm->ops->computemffunction = f;
354   if (ctx) sdm->mffunctionctx = ctx;
355   PetscFunctionReturn(0);
356 }
357 
358 /*@C
359    DMSNESGetFunction - get SNES residual evaluation function
360 
361    Not Collective
362 
363    Input Parameter:
364 .  dm - DM to be used with SNES
365 
366    Output Parameters:
367 +  f - residual evaluation function; see SNESFunction for details
368 -  ctx - context for residual evaluation
369 
370    Level: advanced
371 
372    Note:
373    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
374    associated with the DM.
375 
376 .seealso: `DMSNESSetContext()`, `DMSNESSetFunction()`, `SNESSetFunction()`, `SNESFunction`
377 @*/
378 PetscErrorCode DMSNESGetFunction(DM dm, PetscErrorCode (**f)(SNES, Vec, Vec, void *), void **ctx) {
379   DMSNES sdm;
380 
381   PetscFunctionBegin;
382   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
383   PetscCall(DMGetDMSNES(dm, &sdm));
384   if (f) *f = sdm->ops->computefunction;
385   if (ctx) {
386     if (sdm->functionctxcontainer) PetscCall(PetscContainerGetPointer(sdm->functionctxcontainer, ctx));
387     else *ctx = NULL;
388   }
389   PetscFunctionReturn(0);
390 }
391 
392 /*@C
393    DMSNESSetObjective - set SNES objective evaluation function
394 
395    Not Collective
396 
397    Input Parameters:
398 +  dm - DM to be used with SNES
399 .  obj - objective evaluation function; see SNESObjectiveFunction for details
400 -  ctx - context for residual evaluation
401 
402    Level: advanced
403 
404 .seealso: `DMSNESSetContext()`, `SNESGetObjective()`, `DMSNESSetFunction()`
405 @*/
406 PetscErrorCode DMSNESSetObjective(DM dm, PetscErrorCode (*obj)(SNES, Vec, PetscReal *, void *), void *ctx) {
407   DMSNES sdm;
408 
409   PetscFunctionBegin;
410   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
411   if (obj || ctx) { PetscCall(DMGetDMSNESWrite(dm, &sdm)); }
412   if (obj) sdm->ops->computeobjective = obj;
413   if (ctx) sdm->objectivectx = ctx;
414   PetscFunctionReturn(0);
415 }
416 
417 /*@C
418    DMSNESGetObjective - get SNES objective evaluation function
419 
420    Not Collective
421 
422    Input Parameter:
423 .  dm - DM to be used with SNES
424 
425    Output Parameters:
426 +  obj- residual evaluation function; see SNESObjectiveFunction for details
427 -  ctx - context for residual evaluation
428 
429    Level: advanced
430 
431    Note:
432    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
433    associated with the DM.
434 
435 .seealso: `DMSNESSetContext()`, `DMSNESSetObjective()`, `SNESSetFunction()`
436 @*/
437 PetscErrorCode DMSNESGetObjective(DM dm, PetscErrorCode (**obj)(SNES, Vec, PetscReal *, void *), void **ctx) {
438   DMSNES sdm;
439 
440   PetscFunctionBegin;
441   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
442   PetscCall(DMGetDMSNES(dm, &sdm));
443   if (obj) *obj = sdm->ops->computeobjective;
444   if (ctx) *ctx = sdm->objectivectx;
445   PetscFunctionReturn(0);
446 }
447 
448 /*@C
449    DMSNESSetNGS - set SNES Gauss-Seidel relaxation function
450 
451    Not Collective
452 
453    Input Parameters:
454 +  dm - DM to be used with SNES
455 .  f  - relaxation function, see SNESGSFunction
456 -  ctx - context for residual evaluation
457 
458    Level: advanced
459 
460    Note:
461    SNESSetNGS() is normally used, but it calls this function internally because the user context is actually
462    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
463    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
464 
465 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `DMSNESSetFunction()`, `SNESGSFunction`
466 @*/
467 PetscErrorCode DMSNESSetNGS(DM dm, PetscErrorCode (*f)(SNES, Vec, Vec, void *), void *ctx) {
468   DMSNES sdm;
469 
470   PetscFunctionBegin;
471   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
472   if (f || ctx) { PetscCall(DMGetDMSNESWrite(dm, &sdm)); }
473   if (f) sdm->ops->computegs = f;
474   if (ctx) sdm->gsctx = ctx;
475   PetscFunctionReturn(0);
476 }
477 
478 /*@C
479    DMSNESGetNGS - get SNES Gauss-Seidel relaxation function
480 
481    Not Collective
482 
483    Input Parameter:
484 .  dm - DM to be used with SNES
485 
486    Output Parameters:
487 +  f - relaxation function which performs Gauss-Seidel sweeps, see SNESGSFunction
488 -  ctx - context for residual evaluation
489 
490    Level: advanced
491 
492    Note:
493    SNESGetNGS() is normally used, but it calls this function internally because the user context is actually
494    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
495    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
496 
497 .seealso: `DMSNESSetContext()`, `SNESGetNGS()`, `DMSNESGetJacobian()`, `DMSNESGetFunction()`, `SNESNGSFunction`
498 @*/
499 PetscErrorCode DMSNESGetNGS(DM dm, PetscErrorCode (**f)(SNES, Vec, Vec, void *), void **ctx) {
500   DMSNES sdm;
501 
502   PetscFunctionBegin;
503   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
504   PetscCall(DMGetDMSNES(dm, &sdm));
505   if (f) *f = sdm->ops->computegs;
506   if (ctx) *ctx = sdm->gsctx;
507   PetscFunctionReturn(0);
508 }
509 
510 /*@C
511    DMSNESSetJacobian - set SNES Jacobian evaluation function
512 
513    Not Collective
514 
515    Input Parameters:
516 +  dm - DM to be used with SNES
517 .  J - Jacobian evaluation function
518 -  ctx - context for residual evaluation
519 
520    Level: advanced
521 
522    Note:
523    SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
524    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
525    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
526 
527 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESGetJacobian()`, `SNESSetJacobian()`, `SNESJacobianFunction`
528 @*/
529 PetscErrorCode DMSNESSetJacobian(DM dm, PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *), void *ctx) {
530   DMSNES sdm;
531 
532   PetscFunctionBegin;
533   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
534   if (J || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm));
535   if (J) sdm->ops->computejacobian = J;
536   if (ctx) {
537     PetscContainer ctxcontainer;
538     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)sdm), &ctxcontainer));
539     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
540     PetscCall(PetscObjectCompose((PetscObject)sdm, "jacobian ctx", (PetscObject)ctxcontainer));
541     sdm->jacobianctxcontainer = ctxcontainer;
542     PetscCall(PetscContainerDestroy(&ctxcontainer));
543   }
544   PetscFunctionReturn(0);
545 }
546 
547 /*@C
548    DMSNESSetJacobianContextDestroy - set `SNES` Jacobian evaluation context destroy function
549 
550    Not Collective
551 
552    Input Parameters:
553 +  dm - `DM` to be used with `SNES`
554 -  f - Jacobian evaluation contex destroy function
555 
556    Level: advanced
557 
558 .seealso: `DMSNESSetJacobian()`
559 @*/
560 PetscErrorCode DMSNESSetJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *)) {
561   DMSNES sdm;
562 
563   PetscFunctionBegin;
564   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
565   PetscCall(DMGetDMSNESWrite(dm, &sdm));
566   if (sdm->jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(sdm->jacobianctxcontainer, f));
567   PetscFunctionReturn(0);
568 }
569 
570 PetscErrorCode DMSNESUnsetJacobianContext_Internal(DM dm) {
571   DMSNES sdm;
572 
573   PetscFunctionBegin;
574   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
575   PetscCall(DMGetDMSNESWrite(dm, &sdm));
576   PetscCall(DMSNESUnsetJacobianContext_DMSNES(sdm));
577   PetscFunctionReturn(0);
578 }
579 
580 /*@C
581    DMSNESGetJacobian - get SNES Jacobian evaluation function
582 
583    Not Collective
584 
585    Input Parameter:
586 .  dm - DM to be used with SNES
587 
588    Output Parameters:
589 +  J - Jacobian evaluation function; see SNESJacobianFunction for all calling sequence
590 -  ctx - context for residual evaluation
591 
592    Level: advanced
593 
594    Note:
595    SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
596    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
597    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
598 
599 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESJacobianFunction`
600 @*/
601 PetscErrorCode DMSNESGetJacobian(DM dm, PetscErrorCode (**J)(SNES, Vec, Mat, Mat, void *), void **ctx) {
602   DMSNES sdm;
603 
604   PetscFunctionBegin;
605   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
606   PetscCall(DMGetDMSNES(dm, &sdm));
607   if (J) *J = sdm->ops->computejacobian;
608   if (ctx) {
609     if (sdm->jacobianctxcontainer) PetscCall(PetscContainerGetPointer(sdm->jacobianctxcontainer, ctx));
610     else *ctx = NULL;
611   }
612   PetscFunctionReturn(0);
613 }
614 
615 /*@C
616    DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.
617 
618    Not Collective
619 
620    Input Parameters:
621 +  dm - DM to be used with SNES
622 .  b - RHS evaluation function
623 .  J - Picard matrix evaluation function
624 -  ctx - context for residual evaluation
625 
626    Level: advanced
627 
628 .seealso: `SNESSetPicard()`, `DMSNESSetFunction()`, `DMSNESSetJacobian()`
629 @*/
630 PetscErrorCode DMSNESSetPicard(DM dm, PetscErrorCode (*b)(SNES, Vec, Vec, void *), PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *), void *ctx) {
631   DMSNES sdm;
632 
633   PetscFunctionBegin;
634   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
635   PetscCall(DMGetDMSNES(dm, &sdm));
636   if (b) sdm->ops->computepfunction = b;
637   if (J) sdm->ops->computepjacobian = J;
638   if (ctx) sdm->pctx = ctx;
639   PetscFunctionReturn(0);
640 }
641 
642 /*@C
643    DMSNESGetPicard - get SNES Picard iteration evaluation functions
644 
645    Not Collective
646 
647    Input Parameter:
648 .  dm - DM to be used with SNES
649 
650    Output Parameters:
651 +  b - RHS evaluation function; see SNESFunction for details
652 .  J  - RHS evaluation function; see SNESJacobianFunction for detailsa
653 -  ctx - context for residual evaluation
654 
655    Level: advanced
656 
657 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`
658 @*/
659 PetscErrorCode DMSNESGetPicard(DM dm, PetscErrorCode (**b)(SNES, Vec, Vec, void *), PetscErrorCode (**J)(SNES, Vec, Mat, Mat, void *), void **ctx) {
660   DMSNES sdm;
661 
662   PetscFunctionBegin;
663   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
664   PetscCall(DMGetDMSNES(dm, &sdm));
665   if (b) *b = sdm->ops->computepfunction;
666   if (J) *J = sdm->ops->computepjacobian;
667   if (ctx) *ctx = sdm->pctx;
668   PetscFunctionReturn(0);
669 }
670