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