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