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