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