xref: /petsc/src/snes/utils/dmsnes.c (revision 43208f3fbfb2a0d489265d03b5883f57f56e054f)
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__ "SNESDMComputeFunction"
6 /*@C
7   SNESDMComputeFunction - This is a universal function evaluation routine that
8   may be used with SNESSetFunction() as long as the user context has a DM
9   as its first record and the user has called DMSetLocalFunction().
10 
11   Collective on SNES
12 
13   Input Parameters:
14 + snes - the SNES context
15 . X - input vector
16 . F - function vector
17 - ptr - pointer to a structure that must have a DM as its first entry.
18         This ptr must have been passed into SNESDMComputeFunction() as the context.
19 
20   Level: intermediate
21 
22 .seealso: DMSetLocalFunction(), DMSetLocalJacobian(), SNESSetFunction(), SNESSetJacobian()
23 @*/
24 PetscErrorCode SNESDMComputeFunction(SNES snes, Vec X, Vec F, void *ptr)
25 {
26   DM               dm = *(DM*) ptr;
27   PetscErrorCode (*lf)(DM, Vec, Vec, void *);
28   Vec              localX, localF;
29   PetscInt         N, n;
30   PetscErrorCode   ierr;
31 
32   PetscFunctionBegin;
33   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
34   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
35   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
36   if (!dm) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "Looks like you called SNESSetFromFuntion(snes,SNESDMComputeFunction,) without the DM context");
37   PetscValidHeaderSpecific(dm, DM_CLASSID, 4);
38 
39   /* determine whether X = localX */
40   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
41   ierr = DMGetLocalVector(dm, &localF);CHKERRQ(ierr);
42   ierr = VecGetSize(X, &N);CHKERRQ(ierr);
43   ierr = VecGetSize(localX, &n);CHKERRQ(ierr);
44 
45   if (n != N){ /* X != localX */
46     /* Scatter ghost points to local vector, using the 2-step process
47         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
48     */
49     ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
50     ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
51   } else {
52     ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
53     localX = X;
54   }
55   ierr = DMGetLocalFunction(dm, &lf);CHKERRQ(ierr);
56   ierr = (*lf)(dm, localX, localF, ptr);CHKERRQ(ierr);
57   if (n != N){
58     ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
59   }
60   ierr = DMLocalToGlobalBegin(dm, localF, ADD_VALUES, F);CHKERRQ(ierr);
61   ierr = DMLocalToGlobalEnd(dm, localF, ADD_VALUES, F);CHKERRQ(ierr);
62   ierr = DMRestoreLocalVector(dm, &localF);CHKERRQ(ierr);
63   PetscFunctionReturn(0);
64 }
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "SNESDMComputeJacobian"
68 /*
69   SNESDMComputeJacobian - This is a universal Jacobian evaluation routine that
70   may be used with SNESSetJacobian() as long as the user context has a DM
71   as its first record and the user has called DMSetLocalJacobian().
72 
73   Collective on SNES
74 
75   Input Parameters:
76 + snes - the SNES context
77 . X - input vector
78 . J - Jacobian
79 . B - Jacobian used in preconditioner (usally same as J)
80 . flag - indicates if the matrix changed its structure
81 - ptr - pointer to a structure that must have a DM as its first entry.
82         This ptr must have been passed into SNESDMComputeFunction() as the context.
83 
84   Level: intermediate
85 
86 .seealso: DMSetLocalFunction(), DMSetLocalJacobian(), SNESSetFunction(), SNESSetJacobian()
87 */
88 PetscErrorCode SNESDMComputeJacobian(SNES snes, Vec X, Mat *J, Mat *B, MatStructure *flag, void *ptr)
89 {
90   DM               dm = *(DM*) ptr;
91   PetscErrorCode (*lj)(DM, Vec, Mat, Mat, void *);
92   Vec              localX;
93   PetscErrorCode   ierr;
94 
95   PetscFunctionBegin;
96   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
97   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
98   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
99   ierr = DMGetLocalJacobian(dm, &lj);CHKERRQ(ierr);
100   ierr = (*lj)(dm, localX, *J, *B, ptr);CHKERRQ(ierr);
101   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
102   /* Assemble true Jacobian; if it is different */
103   if (*J != *B) {
104     ierr  = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105     ierr  = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
106   }
107   ierr  = MatSetOption(*B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE);CHKERRQ(ierr);
108   *flag = SAME_NONZERO_PATTERN;
109   PetscFunctionReturn(0);
110 }
111 
112 #undef __FUNCT__
113 #define __FUNCT__ "PetscContainerDestroy_SNESDM"
114 static PetscErrorCode PetscContainerDestroy_SNESDM(void *ctx)
115 {
116   PetscErrorCode ierr;
117   SNESDM sdm = (SNESDM)ctx;
118 
119   PetscFunctionBegin;
120   if (sdm->destroy) {ierr = (*sdm->destroy)(sdm);CHKERRQ(ierr);}
121   ierr = PetscFree(sdm);CHKERRQ(ierr);
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "DMCoarsenHook_SNESDM"
127 /* Attaches the SNESDM to the coarse level.
128  * Under what conditions should we copy versus duplicate?
129  */
130 static PetscErrorCode DMCoarsenHook_SNESDM(DM dm,DM dmc,void *ctx)
131 {
132   PetscErrorCode ierr;
133 
134   PetscFunctionBegin;
135   ierr = DMSNESCopyContext(dm,dmc);CHKERRQ(ierr);
136   PetscFunctionReturn(0);
137 }
138 
139 #undef __FUNCT__
140 #define __FUNCT__ "DMRestrictHook_SNESDM"
141 /* This could restrict auxiliary information to the coarse level.
142  */
143 static PetscErrorCode DMRestrictHook_SNESDM(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
144 {
145 
146   PetscFunctionBegin;
147   PetscFunctionReturn(0);
148 }
149 
150 #undef __FUNCT__
151 #define __FUNCT__ "DMRefineHook_SNESDM"
152 static PetscErrorCode DMRefineHook_SNESDM(DM dm,DM dmf,void *ctx)
153 {
154   PetscErrorCode ierr;
155 
156   PetscFunctionBegin;
157   ierr = DMSNESCopyContext(dm,dmf);CHKERRQ(ierr);
158   PetscFunctionReturn(0);
159 }
160 
161 #undef __FUNCT__
162 #define __FUNCT__ "DMInterpolateHook_SNESDM"
163 /* This could restrict auxiliary information to the coarse level.
164  */
165 static PetscErrorCode DMInterpolateHook_SNESDM(DM dm,Mat Interp,DM dmf,void *ctx)
166 {
167 
168   PetscFunctionBegin;
169   PetscFunctionReturn(0);
170 }
171 
172 #undef __FUNCT__
173 #define __FUNCT__ "DMSNESGetContext"
174 /*@C
175    DMSNESGetContext - get read-only private SNESDM context from a DM
176 
177    Not Collective
178 
179    Input Argument:
180 .  dm - DM to be used with SNES
181 
182    Output Argument:
183 .  snesdm - private SNESDM context
184 
185    Level: developer
186 
187    Notes:
188    Use DMSNESGetContextWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible.
189 
190 .seealso: DMSNESGetContextWrite()
191 @*/
192 PetscErrorCode DMSNESGetContext(DM dm,SNESDM *snesdm)
193 {
194   PetscErrorCode ierr;
195   PetscContainer container;
196   SNESDM         sdm;
197 
198   PetscFunctionBegin;
199   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
200   ierr = PetscObjectQuery((PetscObject)dm,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr);
201   if (container) {
202     ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr);
203   } else {
204     ierr = PetscInfo(dm,"Creating new SNESDM\n");CHKERRQ(ierr);
205     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
206     ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr);
207     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
208     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr);
209     ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
210     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
211     ierr = DMRefineHookAdd(dm,DMRefineHook_SNESDM,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
212     ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr);
213     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
214   }
215   PetscFunctionReturn(0);
216 }
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "DMSNESGetContextWrite"
220 /*@C
221    DMSNESGetContextWrite - get write access to private SNESDM context from a DM
222 
223    Not Collective
224 
225    Input Argument:
226 .  dm - DM to be used with SNES
227 
228    Output Argument:
229 .  snesdm - private SNESDM context
230 
231    Level: developer
232 
233 .seealso: DMSNESGetContext()
234 @*/
235 PetscErrorCode DMSNESGetContextWrite(DM dm,SNESDM *snesdm)
236 {
237   PetscErrorCode ierr;
238   SNESDM         sdm;
239 
240   PetscFunctionBegin;
241   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
242   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
243   if (!sdm->originaldm) sdm->originaldm = dm;
244   if (sdm->originaldm != dm) {  /* Copy on write */
245     PetscContainer container;
246     SNESDM         oldsdm = sdm;
247     ierr = PetscInfo(dm,"Copying SNESDM due to write\n");CHKERRQ(ierr);
248     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
249     ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr);
250     ierr = PetscMemcpy(sdm,oldsdm,sizeof(*sdm));CHKERRQ(ierr);
251     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
252     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr);
253     ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
254     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
255     /* implementation specific copy hooks */
256     if (sdm->duplicate) {ierr = (*sdm->duplicate)(oldsdm,dm);CHKERRQ(ierr);}
257   }
258   *snesdm = sdm;
259   PetscFunctionReturn(0);
260 }
261 
262 #undef __FUNCT__
263 #define __FUNCT__ "DMSNESCopyContext"
264 /*@C
265    DMSNESCopyContext - copies a DM context to a new DM
266 
267    Logically Collective
268 
269    Input Arguments:
270 +  dmsrc - DM to obtain context from
271 -  dmdest - DM to add context to
272 
273    Level: developer
274 
275    Note:
276    The context is copied by reference. This function does not ensure that a context exists.
277 
278 .seealso: DMSNESGetContext(), SNESSetDM()
279 @*/
280 PetscErrorCode DMSNESCopyContext(DM dmsrc,DM dmdest)
281 {
282   PetscErrorCode ierr;
283   PetscContainer container;
284 
285   PetscFunctionBegin;
286   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
287   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
288   ierr = PetscObjectQuery((PetscObject)dmsrc,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr);
289   if (container) {
290     ierr = PetscObjectCompose((PetscObject)dmdest,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
291     ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
292     ierr = DMRefineHookAdd(dmdest,DMRefineHook_SNESDM,DMInterpolateHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
293   }
294   PetscFunctionReturn(0);
295 }
296 
297 #undef __FUNCT__
298 #define __FUNCT__ "DMSNESSetFunction"
299 /*@C
300    DMSNESSetFunction - set SNES residual evaluation function
301 
302    Not Collective
303 
304    Input Arguments:
305 +  dm - DM to be used with SNES
306 .  func - residual evaluation function, see SNESSetFunction() for calling sequence
307 -  ctx - context for residual evaluation
308 
309    Level: advanced
310 
311    Note:
312    SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
313    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
314    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
315 
316 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
317 @*/
318 PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
319 {
320   PetscErrorCode ierr;
321   SNESDM sdm;
322 
323   PetscFunctionBegin;
324   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
325   if (func || ctx) {
326     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
327   }
328   if (func) sdm->computefunction = func;
329   if (ctx)  sdm->functionctx = ctx;
330   PetscFunctionReturn(0);
331 }
332 
333 #undef __FUNCT__
334 #define __FUNCT__ "DMSNESGetFunction"
335 /*@C
336    DMSNESGetFunction - get SNES residual evaluation function
337 
338    Not Collective
339 
340    Input Argument:
341 .  dm - DM to be used with SNES
342 
343    Output Arguments:
344 +  func - residual evaluation function, see SNESSetFunction() for calling sequence
345 -  ctx - context for residual evaluation
346 
347    Level: advanced
348 
349    Note:
350    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
351    associated with the DM.
352 
353 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction()
354 @*/
355 PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
356 {
357   PetscErrorCode ierr;
358   SNESDM sdm;
359 
360   PetscFunctionBegin;
361   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
362   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
363   if (func) *func = sdm->computefunction;
364   if (ctx)  *ctx = sdm->functionctx;
365   PetscFunctionReturn(0);
366 }
367 
368 #undef __FUNCT__
369 #define __FUNCT__ "DMSNESSetObjective"
370 /*@C
371    DMSNESSetObjective - set SNES objective evaluation function
372 
373    Not Collective
374 
375    Input Arguments:
376 +  dm - DM to be used with SNES
377 .  func - residual evaluation function, see SNESSetObjective() for calling sequence
378 -  ctx - context for residual evaluation
379 
380    Level: advanced
381 
382 .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction()
383 @*/
384 PetscErrorCode DMSNESSetObjective(DM dm,SNESObjective func,void *ctx)
385 {
386   PetscErrorCode ierr;
387   SNESDM sdm;
388 
389   PetscFunctionBegin;
390   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
391   if (func || ctx) {
392     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
393   }
394   if (func) sdm->computeobjective = func;
395   if (ctx)  sdm->objectivectx = ctx;
396   PetscFunctionReturn(0);
397 }
398 
399 #undef __FUNCT__
400 #define __FUNCT__ "DMSNESGetObjective"
401 /*@C
402    DMSNESGetObjective - get SNES objective evaluation function
403 
404    Not Collective
405 
406    Input Argument:
407 .  dm - DM to be used with SNES
408 
409    Output Arguments:
410 +  func - residual evaluation function, see SNESSetObjective() for calling sequence
411 -  ctx - context for residual evaluation
412 
413    Level: advanced
414 
415    Note:
416    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
417    associated with the DM.
418 
419 .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction()
420 @*/
421 PetscErrorCode DMSNESGetObjective(DM dm,SNESObjective *func,void **ctx)
422 {
423   PetscErrorCode ierr;
424   SNESDM sdm;
425 
426   PetscFunctionBegin;
427   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
428   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
429   if (func) *func = sdm->computeobjective;
430   if (ctx)  *ctx = sdm->objectivectx;
431   PetscFunctionReturn(0);
432 }
433 
434 #undef __FUNCT__
435 #define __FUNCT__ "DMSNESSetGS"
436 /*@C
437    DMSNESSetGS - set SNES Gauss-Seidel relaxation function
438 
439    Not Collective
440 
441    Input Argument:
442 +  dm - DM to be used with SNES
443 .  func - relaxation function, see SNESSetGS() for calling sequence
444 -  ctx - context for residual evaluation
445 
446    Level: advanced
447 
448    Note:
449    SNESSetGS() is normally used, but it calls this function internally because the user context is actually
450    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
451    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
452 
453 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction()
454 @*/
455 PetscErrorCode DMSNESSetGS(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
456 {
457   PetscErrorCode ierr;
458   SNESDM sdm;
459 
460   PetscFunctionBegin;
461   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
462   if (func || ctx) {
463     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
464   }
465   if (func) sdm->computegs = func;
466   if (ctx)  sdm->gsctx = ctx;
467   PetscFunctionReturn(0);
468 }
469 
470 #undef __FUNCT__
471 #define __FUNCT__ "DMSNESGetGS"
472 /*@C
473    DMSNESGetGS - get SNES Gauss-Seidel relaxation function
474 
475    Not Collective
476 
477    Input Argument:
478 .  dm - DM to be used with SNES
479 
480    Output Arguments:
481 +  func - relaxation function, see SNESSetGS() for calling sequence
482 -  ctx - context for residual evaluation
483 
484    Level: advanced
485 
486    Note:
487    SNESGetGS() is normally used, but it calls this function internally because the user context is actually
488    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
489    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
490 
491 .seealso: DMSNESSetContext(), SNESGetGS(), DMSNESGetJacobian(), DMSNESGetFunction()
492 @*/
493 PetscErrorCode DMSNESGetGS(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
494 {
495   PetscErrorCode ierr;
496   SNESDM sdm;
497 
498   PetscFunctionBegin;
499   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
500   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
501   if (func) *func = sdm->computegs;
502   if (ctx)  *ctx = sdm->gsctx;
503   PetscFunctionReturn(0);
504 }
505 
506 #undef __FUNCT__
507 #define __FUNCT__ "DMSNESSetJacobian"
508 /*@C
509    DMSNESSetJacobian - set SNES Jacobian evaluation function
510 
511    Not Collective
512 
513    Input Argument:
514 +  dm - DM to be used with SNES
515 .  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
516 -  ctx - context for residual evaluation
517 
518    Level: advanced
519 
520    Note:
521    SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
522    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
523    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
524 
525 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian()
526 @*/
527 PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
528 {
529   PetscErrorCode ierr;
530   SNESDM         sdm;
531 
532   PetscFunctionBegin;
533   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
534   if (func || ctx) {
535     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
536   }
537   if (func) sdm->computejacobian = func;
538   if (ctx)  sdm->jacobianctx = ctx;
539   PetscFunctionReturn(0);
540 }
541 
542 #undef __FUNCT__
543 #define __FUNCT__ "DMSNESGetJacobian"
544 /*@C
545    DMSNESGetJacobian - get SNES Jacobian evaluation function
546 
547    Not Collective
548 
549    Input Argument:
550 .  dm - DM to be used with SNES
551 
552    Output Arguments:
553 +  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
554 -  ctx - context for residual evaluation
555 
556    Level: advanced
557 
558    Note:
559    SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
560    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
561    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
562 
563 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
564 @*/
565 PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
566 {
567   PetscErrorCode ierr;
568   SNESDM sdm;
569 
570   PetscFunctionBegin;
571   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
572   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
573   if (func) *func = sdm->computejacobian;
574   if (ctx)  *ctx = sdm->jacobianctx;
575   PetscFunctionReturn(0);
576 }
577 
578 #undef __FUNCT__
579 #define __FUNCT__ "DMSNESSetPicard"
580 /*@C
581    DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.
582 
583    Not Collective
584 
585    Input Argument:
586 +  dm - DM to be used with SNES
587 .  func - RHS evaluation function, see SNESSetFunction() for calling sequence
588 .  pjac - Picard matrix evaluation function, see SNESSetJacobian() for calling sequence
589 -  ctx - context for residual evaluation
590 
591    Level: advanced
592 
593 .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian()
594 @*/
595 PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*pfunc)(SNES,Vec,Vec,void*),PetscErrorCode (*pjac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
596 {
597   PetscErrorCode ierr;
598   SNESDM sdm;
599 
600   PetscFunctionBegin;
601   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
602   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
603   if (pfunc) sdm->computepfunction = pfunc;
604   if (pjac)  sdm->computepjacobian = pjac;
605   if (ctx)   sdm->pctx             = ctx;
606   PetscFunctionReturn(0);
607 }
608 
609 
610 #undef __FUNCT__
611 #define __FUNCT__ "DMSNESGetPicard"
612 /*@C
613    DMSNESGetPicard - get SNES Picard iteration evaluation functions
614 
615    Not Collective
616 
617    Input Argument:
618 .  dm - DM to be used with SNES
619 
620    Output Arguments:
621 +  pfunc - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
622 .  pjac  - RHS evaluation function, see SNESSetFunction() for calling sequence
623 -  ctx - context for residual evaluation
624 
625    Level: advanced
626 
627 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
628 @*/
629 PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**pfunc)(SNES,Vec,Vec,void*),PetscErrorCode (**pjac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
630 {
631   PetscErrorCode ierr;
632   SNESDM sdm;
633 
634   PetscFunctionBegin;
635   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
636   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
637   if (pfunc) *pfunc = sdm->computepfunction;
638   if (pjac) *pjac   = sdm->computepjacobian;
639   if (ctx)  *ctx    = sdm->pctx;
640   PetscFunctionReturn(0);
641 }
642 
643 
644 #undef __FUNCT__
645 #define __FUNCT__ "SNESDefaultComputeFunction_DMLegacy"
646 static PetscErrorCode SNESDefaultComputeFunction_DMLegacy(SNES snes,Vec X,Vec F,void *ctx)
647 {
648   PetscErrorCode ierr;
649   DM             dm;
650 
651   PetscFunctionBegin;
652   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
653   ierr = DMComputeFunction(dm,X,F);CHKERRQ(ierr);
654   PetscFunctionReturn(0);
655 }
656 
657 #undef __FUNCT__
658 #define __FUNCT__ "SNESDefaultComputeJacobian_DMLegacy"
659 static PetscErrorCode SNESDefaultComputeJacobian_DMLegacy(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
660 {
661   PetscErrorCode ierr;
662   DM             dm;
663 
664   PetscFunctionBegin;
665   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
666   ierr = DMComputeJacobian(dm,X,*A,*B,mstr);CHKERRQ(ierr);
667   PetscFunctionReturn(0);
668 }
669 
670 #undef __FUNCT__
671 #define __FUNCT__ "DMSNESSetUpLegacy"
672 /* Sets up calling of legacy DM routines */
673 PetscErrorCode DMSNESSetUpLegacy(DM dm)
674 {
675   PetscErrorCode ierr;
676   SNESDM         sdm;
677 
678   PetscFunctionBegin;
679   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
680   if (!sdm->computefunction) {ierr = DMSNESSetFunction(dm,SNESDefaultComputeFunction_DMLegacy,PETSC_NULL);CHKERRQ(ierr);}
681   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
682   if (!sdm->computejacobian) {
683     if (dm->ops->functionj) {
684       ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobian_DMLegacy,PETSC_NULL);CHKERRQ(ierr);
685     } else {
686       ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobianColor,PETSC_NULL);CHKERRQ(ierr);
687     }
688   }
689   PetscFunctionReturn(0);
690 }
691 
692 /* block functions */
693 
694 #undef __FUNCT__
695 #define __FUNCT__ "DMSNESSetBlockFunction"
696 /*@C
697    DMSNESSetBlockFunction - set SNES residual evaluation function
698 
699    Not Collective
700 
701    Input Arguments:
702 +  dm - DM to be used with SNES
703 .  func - residual evaluation function, see SNESSetFunction() for calling sequence
704 -  ctx - context for residual evaluation
705 
706    Level: developer
707 
708    Note:
709    Mostly for use in DM implementations and transferred to a block function rather than being called from here.
710 
711 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
712 @*/
713 PetscErrorCode DMSNESSetBlockFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
714 {
715   PetscErrorCode ierr;
716   SNESDM sdm;
717 
718   PetscFunctionBegin;
719   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
720   if (func || ctx) {
721     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
722   }
723   if (func) sdm->computeblockfunction = func;
724   if (ctx)  sdm->blockfunctionctx = ctx;
725   PetscFunctionReturn(0);
726 }
727 
728 #undef __FUNCT__
729 #define __FUNCT__ "DMSNESGetBlockFunction"
730 /*@C
731    DMSNESGetBlockFunction - get SNES residual evaluation function
732 
733    Not Collective
734 
735    Input Argument:
736 .  dm - DM to be used with SNES
737 
738    Output Arguments:
739 +  func - residual evaluation function, see SNESSetFunction() for calling sequence
740 -  ctx - context for residual evaluation
741 
742    Level: developer
743 
744 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction()
745 @*/
746 PetscErrorCode DMSNESGetBlockFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
747 {
748   PetscErrorCode ierr;
749   SNESDM sdm;
750 
751   PetscFunctionBegin;
752   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
753   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
754   if (func) *func = sdm->computeblockfunction;
755   if (ctx)  *ctx = sdm->blockfunctionctx;
756   PetscFunctionReturn(0);
757 }
758 
759 
760 #undef __FUNCT__
761 #define __FUNCT__ "DMSNESSetBlockJacobian"
762 /*@C
763    DMSNESSetJacobian - set SNES Jacobian evaluation function
764 
765    Not Collective
766 
767    Input Argument:
768 +  dm - DM to be used with SNES
769 .  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
770 -  ctx - context for residual evaluation
771 
772    Level: advanced
773 
774    Note:
775    Mostly for use in DM implementations and transferred to a block function rather than being called from here.
776 
777 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian()
778 @*/
779 PetscErrorCode DMSNESSetBlockJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
780 {
781   PetscErrorCode ierr;
782   SNESDM sdm;
783 
784   PetscFunctionBegin;
785   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
786   if (func || ctx) {
787     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
788   }
789   if (func) sdm->computeblockjacobian = func;
790   if (ctx)  sdm->blockjacobianctx = ctx;
791   PetscFunctionReturn(0);
792 }
793 
794 #undef __FUNCT__
795 #define __FUNCT__ "DMSNESGetBlockJacobian"
796 /*@C
797    DMSNESGetBlockJacobian - get SNES Jacobian evaluation function
798 
799    Not Collective
800 
801    Input Argument:
802 .  dm - DM to be used with SNES
803 
804    Output Arguments:
805 +  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
806 -  ctx - context for residual evaluation
807 
808    Level: advanced
809 
810 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
811 @*/
812 PetscErrorCode DMSNESGetBlockJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
813 {
814   PetscErrorCode ierr;
815   SNESDM sdm;
816 
817   PetscFunctionBegin;
818   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
819   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
820   if (func) *func = sdm->computeblockjacobian;
821   if (ctx)  *ctx = sdm->blockjacobianctx;
822   PetscFunctionReturn(0);
823 }
824