xref: /petsc/src/snes/utils/dmsnes.c (revision bfcb38ea38335faa6e7f8d97f6bc6ce9aa2a1dd1)
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 
166   /*
167   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
168   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
169   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
170   */
171 
172   /* implementation specific copy hooks */
173   if (kdm->ops->duplicate) {ierr = (*kdm->ops->duplicate)(kdm,nkdm);CHKERRQ(ierr);}
174   PetscFunctionReturn(0);
175 }
176 
177 /*@C
178    DMGetDMSNES - get read-only private DMSNES context from a DM
179 
180    Not Collective
181 
182    Input Argument:
183 .  dm - DM to be used with SNES
184 
185    Output Argument:
186 .  snesdm - private DMSNES context
187 
188    Level: developer
189 
190    Notes:
191    Use DMGetDMSNESWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible.
192 
193 .seealso: DMGetDMSNESWrite()
194 @*/
195 PetscErrorCode DMGetDMSNES(DM dm,DMSNES *snesdm)
196 {
197   PetscErrorCode ierr;
198 
199   PetscFunctionBegin;
200   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
201   *snesdm = (DMSNES) dm->dmsnes;
202   if (!*snesdm) {
203     ierr = PetscInfo(dm,"Creating new DMSNES\n");CHKERRQ(ierr);
204     ierr = DMSNESCreate(PetscObjectComm((PetscObject)dm),snesdm);CHKERRQ(ierr);
205 
206     dm->dmsnes = (PetscObject) *snesdm;
207 
208     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMSNES,DMRestrictHook_DMSNES,NULL);CHKERRQ(ierr);
209     ierr = DMRefineHookAdd(dm,DMRefineHook_DMSNES,DMInterpolateHook_DMSNES,NULL);CHKERRQ(ierr);
210     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);CHKERRQ(ierr);
211   }
212   PetscFunctionReturn(0);
213 }
214 
215 /*@C
216    DMGetDMSNESWrite - get write access to private DMSNES context from a DM
217 
218    Not Collective
219 
220    Input Argument:
221 .  dm - DM to be used with SNES
222 
223    Output Argument:
224 .  snesdm - private DMSNES context
225 
226    Level: developer
227 
228 .seealso: DMGetDMSNES()
229 @*/
230 PetscErrorCode DMGetDMSNESWrite(DM dm,DMSNES *snesdm)
231 {
232   PetscErrorCode ierr;
233   DMSNES         sdm;
234 
235   PetscFunctionBegin;
236   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
237   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
238   if (!sdm->originaldm) sdm->originaldm = dm;
239   if (sdm->originaldm != dm) {  /* Copy on write */
240     DMSNES oldsdm = sdm;
241     ierr       = PetscInfo(dm,"Copying DMSNES due to write\n");CHKERRQ(ierr);
242     ierr       = DMSNESCreate(PetscObjectComm((PetscObject)dm),&sdm);CHKERRQ(ierr);
243     ierr       = DMSNESCopy(oldsdm,sdm);CHKERRQ(ierr);
244     ierr       = DMSNESDestroy((DMSNES*)&dm->dmsnes);CHKERRQ(ierr);
245     dm->dmsnes = (PetscObject)sdm;
246   }
247   *snesdm = sdm;
248   PetscFunctionReturn(0);
249 }
250 
251 /*@C
252    DMCopyDMSNES - copies a DM context to a new DM
253 
254    Logically Collective
255 
256    Input Arguments:
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   PetscErrorCode ierr;
270 
271   PetscFunctionBegin;
272   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
273   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
274   if (!dmdest->dmsnes) {ierr = DMSNESCreate(PetscObjectComm((PetscObject) dmdest), (DMSNES *) &dmdest->dmsnes);CHKERRQ(ierr);}
275   ierr = DMSNESCopy((DMSNES) dmsrc->dmsnes, (DMSNES) dmdest->dmsnes);CHKERRQ(ierr);
276   ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMSNES,NULL,NULL);CHKERRQ(ierr);
277   ierr = DMRefineHookAdd(dmdest,DMRefineHook_DMSNES,NULL,NULL);CHKERRQ(ierr);
278   ierr = DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);CHKERRQ(ierr);
279   PetscFunctionReturn(0);
280 }
281 
282 /*@C
283    DMSNESSetFunction - set SNES residual evaluation function
284 
285    Not Collective
286 
287    Input Arguments:
288 +  dm - DM to be used with SNES
289 .  f - residual evaluation function; see SNESFunction for details
290 -  ctx - context for residual evaluation
291 
292    Level: advanced
293 
294    Note:
295    SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
296    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
297    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
298 
299 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESFunction
300 @*/
301 PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
302 {
303   PetscErrorCode ierr;
304   DMSNES         sdm;
305 
306   PetscFunctionBegin;
307   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
308   if (f || ctx) {
309     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
310   }
311   if (f) sdm->ops->computefunction = f;
312   if (ctx) sdm->functionctx = ctx;
313   PetscFunctionReturn(0);
314 }
315 
316 /*@C
317    DMSNESGetFunction - get SNES residual evaluation function
318 
319    Not Collective
320 
321    Input Argument:
322 .  dm - DM to be used with SNES
323 
324    Output Arguments:
325 +  f - residual evaluation function; see SNESFunction for details
326 -  ctx - context for residual evaluation
327 
328    Level: advanced
329 
330    Note:
331    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
332    associated with the DM.
333 
334 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction(), SNESFunction
335 @*/
336 PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
337 {
338   PetscErrorCode ierr;
339   DMSNES         sdm;
340 
341   PetscFunctionBegin;
342   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
343   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
344   if (f) *f = sdm->ops->computefunction;
345   if (ctx) *ctx = sdm->functionctx;
346   PetscFunctionReturn(0);
347 }
348 
349 /*@C
350    DMSNESSetObjective - set SNES objective evaluation function
351 
352    Not Collective
353 
354    Input Arguments:
355 +  dm - DM to be used with SNES
356 .  obj - objective evaluation function; see SNESObjectiveFunction for details
357 -  ctx - context for residual evaluation
358 
359    Level: advanced
360 
361 .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction()
362 @*/
363 PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx)
364 {
365   PetscErrorCode ierr;
366   DMSNES         sdm;
367 
368   PetscFunctionBegin;
369   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
370   if (obj || ctx) {
371     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
372   }
373   if (obj) sdm->ops->computeobjective = obj;
374   if (ctx) sdm->objectivectx = ctx;
375   PetscFunctionReturn(0);
376 }
377 
378 /*@C
379    DMSNESGetObjective - get SNES objective evaluation function
380 
381    Not Collective
382 
383    Input Argument:
384 .  dm - DM to be used with SNES
385 
386    Output Arguments:
387 +  obj- residual evaluation function; see SNESObjectiveFunction for details
388 -  ctx - context for residual evaluation
389 
390    Level: advanced
391 
392    Note:
393    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
394    associated with the DM.
395 
396 .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction()
397 @*/
398 PetscErrorCode DMSNESGetObjective(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   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
406   if (obj) *obj = sdm->ops->computeobjective;
407   if (ctx) *ctx = sdm->objectivectx;
408   PetscFunctionReturn(0);
409 }
410 
411 /*@C
412    DMSNESSetNGS - set SNES Gauss-Seidel relaxation function
413 
414    Not Collective
415 
416    Input Argument:
417 +  dm - DM to be used with SNES
418 .  f  - relaxation function, see SNESGSFunction
419 -  ctx - context for residual evaluation
420 
421    Level: advanced
422 
423    Note:
424    SNESSetNGS() is normally used, but it calls this function internally because the user context is actually
425    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
426    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
427 
428 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction(), SNESGSFunction
429 @*/
430 PetscErrorCode DMSNESSetNGS(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
431 {
432   PetscErrorCode ierr;
433   DMSNES         sdm;
434 
435   PetscFunctionBegin;
436   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
437   if (f || ctx) {
438     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
439   }
440   if (f) sdm->ops->computegs = f;
441   if (ctx) sdm->gsctx = ctx;
442   PetscFunctionReturn(0);
443 }
444 
445 /*@C
446    DMSNESGetNGS - get SNES Gauss-Seidel relaxation function
447 
448    Not Collective
449 
450    Input Argument:
451 .  dm - DM to be used with SNES
452 
453    Output Arguments:
454 +  f - relaxation function which performs Gauss-Seidel sweeps, see SNESGSFunction
455 -  ctx - context for residual evaluation
456 
457    Level: advanced
458 
459    Note:
460    SNESGetNGS() is normally used, but it calls this function internally because the user context is actually
461    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
462    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
463 
464 .seealso: DMSNESSetContext(), SNESGetNGS(), DMSNESGetJacobian(), DMSNESGetFunction(), SNESNGSFunction
465 @*/
466 PetscErrorCode DMSNESGetNGS(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
467 {
468   PetscErrorCode ierr;
469   DMSNES         sdm;
470 
471   PetscFunctionBegin;
472   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
473   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
474   if (f) *f = sdm->ops->computegs;
475   if (ctx) *ctx = sdm->gsctx;
476   PetscFunctionReturn(0);
477 }
478 
479 /*@C
480    DMSNESSetJacobian - set SNES Jacobian evaluation function
481 
482    Not Collective
483 
484    Input Argument:
485 +  dm - DM to be used with SNES
486 .  J - Jacobian evaluation function
487 -  ctx - context for residual evaluation
488 
489    Level: advanced
490 
491    Note:
492    SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
493    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
494    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
495 
496 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian(), SNESJacobianFunction
497 @*/
498 PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
499 {
500   PetscErrorCode ierr;
501   DMSNES         sdm;
502 
503   PetscFunctionBegin;
504   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
505   if (J || ctx) {
506     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
507   }
508   if (J) sdm->ops->computejacobian = J;
509   if (ctx) sdm->jacobianctx = ctx;
510   PetscFunctionReturn(0);
511 }
512 
513 /*@C
514    DMSNESGetJacobian - get SNES Jacobian evaluation function
515 
516    Not Collective
517 
518    Input Argument:
519 .  dm - DM to be used with SNES
520 
521    Output Arguments:
522 +  J - Jacobian evaluation function; see SNESJacobianFunction for all calling sequence
523 -  ctx - context for residual evaluation
524 
525    Level: advanced
526 
527    Note:
528    SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
529    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
530    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
531 
532 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESJacobianFunction
533 @*/
534 PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
535 {
536   PetscErrorCode ierr;
537   DMSNES         sdm;
538 
539   PetscFunctionBegin;
540   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
541   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
542   if (J) *J = sdm->ops->computejacobian;
543   if (ctx) *ctx = sdm->jacobianctx;
544   PetscFunctionReturn(0);
545 }
546 
547 /*@C
548    DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.
549 
550    Not Collective
551 
552    Input Argument:
553 +  dm - DM to be used with SNES
554 .  b - RHS evaluation function
555 .  J - Picard matrix evaluation function
556 -  ctx - context for residual evaluation
557 
558    Level: advanced
559 
560 .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian()
561 @*/
562 PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*b)(SNES,Vec,Vec,void*),PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
563 {
564   PetscErrorCode ierr;
565   DMSNES         sdm;
566 
567   PetscFunctionBegin;
568   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
569   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
570   if (b) sdm->ops->computepfunction = b;
571   if (J) sdm->ops->computepjacobian = J;
572   if (ctx) sdm->pctx = ctx;
573   PetscFunctionReturn(0);
574 }
575 
576 /*@C
577    DMSNESGetPicard - get SNES Picard iteration evaluation functions
578 
579    Not Collective
580 
581    Input Argument:
582 .  dm - DM to be used with SNES
583 
584    Output Arguments:
585 +  b - RHS evaluation function; see SNESFunction for details
586 .  J  - RHS evaluation function; see SNESJacobianFunction for detailsa
587 -  ctx - context for residual evaluation
588 
589    Level: advanced
590 
591 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
592 @*/
593 PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**b)(SNES,Vec,Vec,void*),PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
594 {
595   PetscErrorCode ierr;
596   DMSNES         sdm;
597 
598   PetscFunctionBegin;
599   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
600   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
601   if (b) *b = sdm->ops->computepfunction;
602   if (J) *J = sdm->ops->computepjacobian;
603   if (ctx) *ctx = sdm->pctx;
604   PetscFunctionReturn(0);
605 }
606