xref: /petsc/src/snes/utils/dmsnes.c (revision 609bdbee21ea3be08735c64dbe00a9ab27759925)
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 = 0; 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,PETSC_FALSE);CHKERRQ(ierr);
58     ierr = PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION,PETSC_FALSE);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 Argument:
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->data         = kdm->data;
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 Argument:
184 .  dm - DM to be used with SNES
185 
186    Output Argument:
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 
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 Argument:
222 .  dm - DM to be used with SNES
223 
224    Output Argument:
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) sdm->originaldm = dm;
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   }
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 Arguments:
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   PetscErrorCode ierr;
271 
272   PetscFunctionBegin;
273   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
274   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
275   if (!dmdest->dmsnes) {ierr = DMSNESCreate(PetscObjectComm((PetscObject) dmdest), (DMSNES *) &dmdest->dmsnes);CHKERRQ(ierr);}
276   ierr = DMSNESCopy((DMSNES) dmsrc->dmsnes, (DMSNES) dmdest->dmsnes);CHKERRQ(ierr);
277   ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMSNES,NULL,NULL);CHKERRQ(ierr);
278   ierr = DMRefineHookAdd(dmdest,DMRefineHook_DMSNES,NULL,NULL);CHKERRQ(ierr);
279   ierr = DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);CHKERRQ(ierr);
280   PetscFunctionReturn(0);
281 }
282 
283 /*@C
284    DMSNESSetFunction - set SNES residual evaluation function
285 
286    Not Collective
287 
288    Input Arguments:
289 +  dm - DM to be used with SNES
290 .  f - residual evaluation function; see SNESFunction for details
291 -  ctx - context for residual evaluation
292 
293    Level: advanced
294 
295    Note:
296    SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
297    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
298    not. 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   PetscErrorCode ierr;
305   DMSNES         sdm;
306 
307   PetscFunctionBegin;
308   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
309   if (f || ctx) {
310     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
311   }
312   if (f) sdm->ops->computefunction = f;
313   if (ctx) sdm->functionctx = ctx;
314   PetscFunctionReturn(0);
315 }
316 
317 /*@C
318    DMSNESGetFunction - get SNES residual evaluation function
319 
320    Not Collective
321 
322    Input Argument:
323 .  dm - DM to be used with SNES
324 
325    Output Arguments:
326 +  f - residual evaluation function; see SNESFunction for details
327 -  ctx - context for residual evaluation
328 
329    Level: advanced
330 
331    Note:
332    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
333    associated with the DM.
334 
335 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction(), SNESFunction
336 @*/
337 PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
338 {
339   PetscErrorCode ierr;
340   DMSNES         sdm;
341 
342   PetscFunctionBegin;
343   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
344   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
345   if (f) *f = sdm->ops->computefunction;
346   if (ctx) *ctx = sdm->functionctx;
347   PetscFunctionReturn(0);
348 }
349 
350 /*@C
351    DMSNESSetObjective - set SNES objective evaluation function
352 
353    Not Collective
354 
355    Input Arguments:
356 +  dm - DM to be used with SNES
357 .  obj - objective evaluation function; see SNESObjectiveFunction for details
358 -  ctx - context for residual evaluation
359 
360    Level: advanced
361 
362 .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction()
363 @*/
364 PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx)
365 {
366   PetscErrorCode ierr;
367   DMSNES         sdm;
368 
369   PetscFunctionBegin;
370   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
371   if (obj || ctx) {
372     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
373   }
374   if (obj) sdm->ops->computeobjective = obj;
375   if (ctx) sdm->objectivectx = ctx;
376   PetscFunctionReturn(0);
377 }
378 
379 /*@C
380    DMSNESGetObjective - get SNES objective evaluation function
381 
382    Not Collective
383 
384    Input Argument:
385 .  dm - DM to be used with SNES
386 
387    Output Arguments:
388 +  obj- residual evaluation function; see SNESObjectiveFunction for details
389 -  ctx - context for residual evaluation
390 
391    Level: advanced
392 
393    Note:
394    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
395    associated with the DM.
396 
397 .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction()
398 @*/
399 PetscErrorCode DMSNESGetObjective(DM dm,PetscErrorCode (**obj)(SNES,Vec,PetscReal*,void*),void **ctx)
400 {
401   PetscErrorCode ierr;
402   DMSNES         sdm;
403 
404   PetscFunctionBegin;
405   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
406   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
407   if (obj) *obj = sdm->ops->computeobjective;
408   if (ctx) *ctx = sdm->objectivectx;
409   PetscFunctionReturn(0);
410 }
411 
412 /*@C
413    DMSNESSetNGS - set SNES Gauss-Seidel relaxation function
414 
415    Not Collective
416 
417    Input Argument:
418 +  dm - DM to be used with SNES
419 .  f  - relaxation function, see SNESGSFunction
420 -  ctx - context for residual evaluation
421 
422    Level: advanced
423 
424    Note:
425    SNESSetNGS() is normally used, but it calls this function internally because the user context is actually
426    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
427    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
428 
429 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction(), SNESGSFunction
430 @*/
431 PetscErrorCode DMSNESSetNGS(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
432 {
433   PetscErrorCode ierr;
434   DMSNES         sdm;
435 
436   PetscFunctionBegin;
437   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
438   if (f || ctx) {
439     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
440   }
441   if (f) sdm->ops->computegs = f;
442   if (ctx) sdm->gsctx = ctx;
443   PetscFunctionReturn(0);
444 }
445 
446 /*@C
447    DMSNESGetNGS - get SNES Gauss-Seidel relaxation function
448 
449    Not Collective
450 
451    Input Argument:
452 .  dm - DM to be used with SNES
453 
454    Output Arguments:
455 +  f - relaxation function which performs Gauss-Seidel sweeps, see SNESGSFunction
456 -  ctx - context for residual evaluation
457 
458    Level: advanced
459 
460    Note:
461    SNESGetNGS() is normally used, but it calls this function internally because the user context is actually
462    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
463    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
464 
465 .seealso: DMSNESSetContext(), SNESGetNGS(), DMSNESGetJacobian(), DMSNESGetFunction(), SNESNGSFunction
466 @*/
467 PetscErrorCode DMSNESGetNGS(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
468 {
469   PetscErrorCode ierr;
470   DMSNES         sdm;
471 
472   PetscFunctionBegin;
473   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
474   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
475   if (f) *f = sdm->ops->computegs;
476   if (ctx) *ctx = sdm->gsctx;
477   PetscFunctionReturn(0);
478 }
479 
480 /*@C
481    DMSNESSetJacobian - set SNES Jacobian evaluation function
482 
483    Not Collective
484 
485    Input Argument:
486 +  dm - DM to be used with SNES
487 .  J - Jacobian evaluation function
488 -  ctx - context for residual evaluation
489 
490    Level: advanced
491 
492    Note:
493    SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
494    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
495    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
496 
497 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian(), SNESJacobianFunction
498 @*/
499 PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
500 {
501   PetscErrorCode ierr;
502   DMSNES         sdm;
503 
504   PetscFunctionBegin;
505   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
506   if (J || ctx) {
507     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
508   }
509   if (J) sdm->ops->computejacobian = J;
510   if (ctx) sdm->jacobianctx = ctx;
511   PetscFunctionReturn(0);
512 }
513 
514 /*@C
515    DMSNESGetJacobian - get SNES Jacobian evaluation function
516 
517    Not Collective
518 
519    Input Argument:
520 .  dm - DM to be used with SNES
521 
522    Output Arguments:
523 +  J - Jacobian evaluation function; see SNESJacobianFunction for all calling sequence
524 -  ctx - context for residual evaluation
525 
526    Level: advanced
527 
528    Note:
529    SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
530    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
531    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
532 
533 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESJacobianFunction
534 @*/
535 PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
536 {
537   PetscErrorCode ierr;
538   DMSNES         sdm;
539 
540   PetscFunctionBegin;
541   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
542   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
543   if (J) *J = sdm->ops->computejacobian;
544   if (ctx) *ctx = sdm->jacobianctx;
545   PetscFunctionReturn(0);
546 }
547 
548 /*@C
549    DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.
550 
551    Not Collective
552 
553    Input Argument:
554 +  dm - DM to be used with SNES
555 .  b - RHS evaluation function
556 .  J - Picard matrix evaluation function
557 -  ctx - context for residual evaluation
558 
559    Level: advanced
560 
561 .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian()
562 @*/
563 PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*b)(SNES,Vec,Vec,void*),PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
564 {
565   PetscErrorCode ierr;
566   DMSNES         sdm;
567 
568   PetscFunctionBegin;
569   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
570   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
571   if (b) sdm->ops->computepfunction = b;
572   if (J) sdm->ops->computepjacobian = J;
573   if (ctx) sdm->pctx = ctx;
574   PetscFunctionReturn(0);
575 }
576 
577 /*@C
578    DMSNESGetPicard - get SNES Picard iteration evaluation functions
579 
580    Not Collective
581 
582    Input Argument:
583 .  dm - DM to be used with SNES
584 
585    Output Arguments:
586 +  b - RHS evaluation function; see SNESFunction for details
587 .  J  - RHS evaluation function; see SNESJacobianFunction for detailsa
588 -  ctx - context for residual evaluation
589 
590    Level: advanced
591 
592 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
593 @*/
594 PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**b)(SNES,Vec,Vec,void*),PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
595 {
596   PetscErrorCode ierr;
597   DMSNES         sdm;
598 
599   PetscFunctionBegin;
600   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
601   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
602   if (b) *b = sdm->ops->computepfunction;
603   if (J) *J = sdm->ops->computepjacobian;
604   if (ctx) *ctx = sdm->pctx;
605   PetscFunctionReturn(0);
606 }
607