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