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