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