1 #include <petsc/private/snesimpl.h>   /*I "petscsnes.h" I*/
2 #include <petsc/private/dmimpl.h>     /*I "petscdm.h" I*/
3 
DMSNESDestroy(DMSNES * kdm)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 
DMSNESLoad(DMSNES kdm,PetscViewer viewer)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 
DMSNESView(DMSNES kdm,PetscViewer viewer)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 
DMSNESCreate(MPI_Comm comm,DMSNES * kdm)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  */
DMCoarsenHook_DMSNES(DM dm,DM dmc,void * ctx)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  */
DMRestrictHook_DMSNES(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void * ctx)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. */
DMSubDomainHook_DMSNES(DM dm,DM subdm,void * ctx)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  */
DMSubDomainRestrictHook_DMSNES(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void * ctx)106 static PetscErrorCode DMSubDomainRestrictHook_DMSNES(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
107 {
108 
109   PetscFunctionBegin;
110   PetscFunctionReturn(0);
111 }
112 
DMRefineHook_DMSNES(DM dm,DM dmf,void * ctx)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  */
DMInterpolateHook_DMSNES(DM dm,Mat Interp,DM dmf,void * ctx)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 @*/
DMSNESCopy(DMSNES kdm,DMSNES nkdm)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 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 @*/
DMGetDMSNES(DM dm,DMSNES * snesdm)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 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 @*/
DMGetDMSNESWrite(DM dm,DMSNES * snesdm)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 Arguments:
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 @*/
DMCopyDMSNES(DM dmsrc,DM dmdest)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 Arguments:
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 @*/
DMSNESSetFunction(DM dm,PetscErrorCode (* f)(SNES,Vec,Vec,void *),void * ctx)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    DMSNESGetFunction - get SNES residual evaluation function
320 
321    Not Collective
322 
323    Input Argument:
324 .  dm - DM to be used with SNES
325 
326    Output Arguments:
327 +  f - residual evaluation function; see SNESFunction for details
328 -  ctx - context for residual evaluation
329 
330    Level: advanced
331 
332    Note:
333    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
334    associated with the DM.
335 
336 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction(), SNESFunction
337 @*/
DMSNESGetFunction(DM dm,PetscErrorCode (** f)(SNES,Vec,Vec,void *),void ** ctx)338 PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
339 {
340   PetscErrorCode ierr;
341   DMSNES         sdm;
342 
343   PetscFunctionBegin;
344   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
345   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
346   if (f) *f = sdm->ops->computefunction;
347   if (ctx) *ctx = sdm->functionctx;
348   PetscFunctionReturn(0);
349 }
350 
351 /*@C
352    DMSNESSetObjective - set SNES objective evaluation function
353 
354    Not Collective
355 
356    Input Arguments:
357 +  dm - DM to be used with SNES
358 .  obj - objective evaluation function; see SNESObjectiveFunction for details
359 -  ctx - context for residual evaluation
360 
361    Level: advanced
362 
363 .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction()
364 @*/
DMSNESSetObjective(DM dm,PetscErrorCode (* obj)(SNES,Vec,PetscReal *,void *),void * ctx)365 PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx)
366 {
367   PetscErrorCode ierr;
368   DMSNES         sdm;
369 
370   PetscFunctionBegin;
371   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
372   if (obj || ctx) {
373     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
374   }
375   if (obj) sdm->ops->computeobjective = obj;
376   if (ctx) sdm->objectivectx = ctx;
377   PetscFunctionReturn(0);
378 }
379 
380 /*@C
381    DMSNESGetObjective - get SNES objective evaluation function
382 
383    Not Collective
384 
385    Input Argument:
386 .  dm - DM to be used with SNES
387 
388    Output Arguments:
389 +  obj- residual evaluation function; see SNESObjectiveFunction for details
390 -  ctx - context for residual evaluation
391 
392    Level: advanced
393 
394    Note:
395    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
396    associated with the DM.
397 
398 .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction()
399 @*/
DMSNESGetObjective(DM dm,PetscErrorCode (** obj)(SNES,Vec,PetscReal *,void *),void ** ctx)400 PetscErrorCode DMSNESGetObjective(DM dm,PetscErrorCode (**obj)(SNES,Vec,PetscReal*,void*),void **ctx)
401 {
402   PetscErrorCode ierr;
403   DMSNES         sdm;
404 
405   PetscFunctionBegin;
406   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
407   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
408   if (obj) *obj = sdm->ops->computeobjective;
409   if (ctx) *ctx = sdm->objectivectx;
410   PetscFunctionReturn(0);
411 }
412 
413 /*@C
414    DMSNESSetNGS - set SNES Gauss-Seidel relaxation function
415 
416    Not Collective
417 
418    Input Argument:
419 +  dm - DM to be used with SNES
420 .  f  - relaxation function, see SNESGSFunction
421 -  ctx - context for residual evaluation
422 
423    Level: advanced
424 
425    Note:
426    SNESSetNGS() is normally used, but it calls this function internally because the user context is actually
427    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
428    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
429 
430 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction(), SNESGSFunction
431 @*/
DMSNESSetNGS(DM dm,PetscErrorCode (* f)(SNES,Vec,Vec,void *),void * ctx)432 PetscErrorCode DMSNESSetNGS(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
433 {
434   PetscErrorCode ierr;
435   DMSNES         sdm;
436 
437   PetscFunctionBegin;
438   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
439   if (f || ctx) {
440     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
441   }
442   if (f) sdm->ops->computegs = f;
443   if (ctx) sdm->gsctx = ctx;
444   PetscFunctionReturn(0);
445 }
446 
447 /*@C
448    DMSNESGetNGS - get SNES Gauss-Seidel relaxation function
449 
450    Not Collective
451 
452    Input Argument:
453 .  dm - DM to be used with SNES
454 
455    Output Arguments:
456 +  f - relaxation function which performs Gauss-Seidel sweeps, see SNESGSFunction
457 -  ctx - context for residual evaluation
458 
459    Level: advanced
460 
461    Note:
462    SNESGetNGS() is normally used, but it calls this function internally because the user context is actually
463    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
464    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
465 
466 .seealso: DMSNESSetContext(), SNESGetNGS(), DMSNESGetJacobian(), DMSNESGetFunction(), SNESNGSFunction
467 @*/
DMSNESGetNGS(DM dm,PetscErrorCode (** f)(SNES,Vec,Vec,void *),void ** ctx)468 PetscErrorCode DMSNESGetNGS(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
469 {
470   PetscErrorCode ierr;
471   DMSNES         sdm;
472 
473   PetscFunctionBegin;
474   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
475   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
476   if (f) *f = sdm->ops->computegs;
477   if (ctx) *ctx = sdm->gsctx;
478   PetscFunctionReturn(0);
479 }
480 
481 /*@C
482    DMSNESSetJacobian - set SNES Jacobian evaluation function
483 
484    Not Collective
485 
486    Input Argument:
487 +  dm - DM to be used with SNES
488 .  J - Jacobian evaluation function
489 -  ctx - context for residual evaluation
490 
491    Level: advanced
492 
493    Note:
494    SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
495    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
496    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
497 
498 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian(), SNESJacobianFunction
499 @*/
DMSNESSetJacobian(DM dm,PetscErrorCode (* J)(SNES,Vec,Mat,Mat,void *),void * ctx)500 PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
501 {
502   PetscErrorCode ierr;
503   DMSNES         sdm;
504 
505   PetscFunctionBegin;
506   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
507   if (J || ctx) {
508     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
509   }
510   if (J) sdm->ops->computejacobian = J;
511   if (ctx) sdm->jacobianctx = ctx;
512   PetscFunctionReturn(0);
513 }
514 
515 /*@C
516    DMSNESGetJacobian - get SNES Jacobian evaluation function
517 
518    Not Collective
519 
520    Input Argument:
521 .  dm - DM to be used with SNES
522 
523    Output Arguments:
524 +  J - Jacobian evaluation function; see SNESJacobianFunction for all calling sequence
525 -  ctx - context for residual evaluation
526 
527    Level: advanced
528 
529    Note:
530    SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
531    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
532    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
533 
534 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESJacobianFunction
535 @*/
DMSNESGetJacobian(DM dm,PetscErrorCode (** J)(SNES,Vec,Mat,Mat,void *),void ** ctx)536 PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
537 {
538   PetscErrorCode ierr;
539   DMSNES         sdm;
540 
541   PetscFunctionBegin;
542   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
543   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
544   if (J) *J = sdm->ops->computejacobian;
545   if (ctx) *ctx = sdm->jacobianctx;
546   PetscFunctionReturn(0);
547 }
548 
549 /*@C
550    DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.
551 
552    Not Collective
553 
554    Input Argument:
555 +  dm - DM to be used with SNES
556 .  b - RHS evaluation function
557 .  J - Picard matrix evaluation function
558 -  ctx - context for residual evaluation
559 
560    Level: advanced
561 
562 .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian()
563 @*/
DMSNESSetPicard(DM dm,PetscErrorCode (* b)(SNES,Vec,Vec,void *),PetscErrorCode (* J)(SNES,Vec,Mat,Mat,void *),void * ctx)564 PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*b)(SNES,Vec,Vec,void*),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 (b) sdm->ops->computepfunction = b;
573   if (J) sdm->ops->computepjacobian = J;
574   if (ctx) sdm->pctx = ctx;
575   PetscFunctionReturn(0);
576 }
577 
578 /*@C
579    DMSNESGetPicard - get SNES Picard iteration evaluation functions
580 
581    Not Collective
582 
583    Input Argument:
584 .  dm - DM to be used with SNES
585 
586    Output Arguments:
587 +  b - RHS evaluation function; see SNESFunction for details
588 .  J  - RHS evaluation function; see SNESJacobianFunction for detailsa
589 -  ctx - context for residual evaluation
590 
591    Level: advanced
592 
593 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
594 @*/
DMSNESGetPicard(DM dm,PetscErrorCode (** b)(SNES,Vec,Vec,void *),PetscErrorCode (** J)(SNES,Vec,Mat,Mat,void *),void ** ctx)595 PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**b)(SNES,Vec,Vec,void*),PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
596 {
597   PetscErrorCode ierr;
598   DMSNES         sdm;
599 
600   PetscFunctionBegin;
601   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
602   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
603   if (b) *b = sdm->ops->computepfunction;
604   if (J) *J = sdm->ops->computepjacobian;
605   if (ctx) *ctx = sdm->pctx;
606   PetscFunctionReturn(0);
607 }
608