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