xref: /petsc/src/ts/utils/dmts.c (revision 8cbdbec6f8317ddf7886f91eb9c6bd083b543c50)
1 #include <petsc-private/tsimpl.h>     /*I "petscts.h" I*/
2 #include <petsc-private/dmimpl.h>     /*I "petscdm.h" I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMTSDestroy"
6 static PetscErrorCode DMTSDestroy(DMTS *kdm)
7 {
8   PetscErrorCode ierr;
9 
10   PetscFunctionBegin;
11  if (!*kdm) PetscFunctionReturn(0);
12   PetscValidHeaderSpecific((*kdm),DMTS_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__ "DMTSLoad"
21 PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer)
22 {
23   PetscErrorCode ierr;
24 
25   PetscFunctionBegin;
26   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,PETSC_FUNCTION);CHKERRQ(ierr);
27   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,PETSC_FUNCTION);CHKERRQ(ierr);
28   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,PETSC_FUNCTION);CHKERRQ(ierr);
29   if (kdm->ops->ifunctionload) {
30     ierr = (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);CHKERRQ(ierr);
31   }
32   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,PETSC_FUNCTION);CHKERRQ(ierr);
33   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,PETSC_FUNCTION);CHKERRQ(ierr);
34   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,PETSC_FUNCTION);CHKERRQ(ierr);
35   if (kdm->ops->ijacobianload) {
36     ierr = (*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer);CHKERRQ(ierr);
37   }
38   PetscFunctionReturn(0);
39 }
40 
41 #undef __FUNCT__
42 #define __FUNCT__ "DMTSView"
43 PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer)
44 {
45   PetscErrorCode ierr;
46   PetscBool      isascii,isbinary;
47 
48   PetscFunctionBegin;
49   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
50   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
51   if (isascii) {
52 #if defined(PETSC_SERIALIZE_FUNCTIONS)
53     const char *fname;
54 
55     ierr = PetscFPTFind(kdm->ops->ifunction,&fname);CHKERRQ(ierr);
56     if (fname) {
57       ierr = PetscViewerASCIIPrintf(viewer,"  IFunction used by TS: %s\n",fname);CHKERRQ(ierr);
58     }
59     ierr = PetscFPTFind(kdm->ops->ijacobian,&fname);CHKERRQ(ierr);
60     if (fname) {
61       ierr = PetscViewerASCIIPrintf(viewer,"  IJacobian function used by TS: %s\n",fname);CHKERRQ(ierr);
62     }
63 #endif
64   } else if (isbinary) {
65     ierr = PetscViewerBinaryWrite(viewer,(void*)kdm->ops->ifunction,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
66     ierr = PetscViewerBinaryWrite(viewer,(void*)kdm->ops->ifunctionview,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
67     ierr = PetscViewerBinaryWrite(viewer,(void*)kdm->ops->ifunctionload,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
68     if (kdm->ops->ifunctionview) {
69       ierr = (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);CHKERRQ(ierr);
70     }
71     ierr = PetscViewerBinaryWrite(viewer,(void*)kdm->ops->ijacobian,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
72     ierr = PetscViewerBinaryWrite(viewer,(void*)kdm->ops->ijacobianview,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
73     ierr = PetscViewerBinaryWrite(viewer,(void*)kdm->ops->ijacobianload,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
74     if (kdm->ops->ijacobianview) {
75       ierr = (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);CHKERRQ(ierr);
76     }
77   }
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "DMTSCreate"
83 static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
84 {
85   PetscErrorCode ierr;
86 
87   PetscFunctionBegin;
88 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
89   ierr = TSInitializePackage(PETSC_NULL);CHKERRQ(ierr);
90 #endif
91   ierr = PetscHeaderCreate(*kdm, _p_DMTS, struct _DMTSOps, DMTS_CLASSID, -1, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);CHKERRQ(ierr);
92   ierr = PetscMemzero((*kdm)->ops, sizeof(struct _DMTSOps));CHKERRQ(ierr);
93   PetscFunctionReturn(0);
94 }
95 
96 #undef __FUNCT__
97 #define __FUNCT__ "DMCoarsenHook_DMTS"
98 /* Attaches the DMTS to the coarse level.
99  * Under what conditions should we copy versus duplicate?
100  */
101 static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
102 {
103   PetscErrorCode ierr;
104 
105   PetscFunctionBegin;
106   ierr = DMCopyDMTS(dm,dmc);CHKERRQ(ierr);
107   PetscFunctionReturn(0);
108 }
109 
110 #undef __FUNCT__
111 #define __FUNCT__ "DMRestrictHook_DMTS"
112 /* This could restrict auxiliary information to the coarse level.
113  */
114 static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
115 {
116 
117   PetscFunctionBegin;
118   PetscFunctionReturn(0);
119 }
120 
121 #undef __FUNCT__
122 #define __FUNCT__ "DMSubDomainHook_DMTS"
123 static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
124 {
125   PetscErrorCode ierr;
126 
127   PetscFunctionBegin;
128   ierr = DMCopyDMTS(dm,subdm);CHKERRQ(ierr);
129   PetscFunctionReturn(0);
130 }
131 
132 #undef __FUNCT__
133 #define __FUNCT__ "DMSubDomainRestrictHook_DMTS"
134 /* This could restrict auxiliary information to the coarse level.
135  */
136 static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
137 {
138   PetscFunctionBegin;
139   PetscFunctionReturn(0);
140 }
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "DMTSCopy"
144 /*@C
145    DMTSCopy - copies the information in a DMTS to another DMTS
146 
147    Not Collective
148 
149    Input Argument:
150 +  kdm - Original DMTS
151 -  nkdm - DMTS to receive the data, should have been created with DMTSCreate()
152 
153    Level: developer
154 
155 .seealso: DMTSCreate(), DMTSDestroy()
156 @*/
157 PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
158 {
159   PetscErrorCode ierr;
160 
161   PetscFunctionBegin;
162   PetscValidHeaderSpecific(kdm,DMTS_CLASSID,1);
163   PetscValidHeaderSpecific(nkdm,DMTS_CLASSID,2);
164   nkdm->ops->rhsfunction    = kdm->ops->rhsfunction;
165   nkdm->ops->rhsjacobian    = kdm->ops->rhsjacobian;
166   nkdm->ops->ifunction      = kdm->ops->ifunction;
167   nkdm->ops->ijacobian      = kdm->ops->ijacobian;
168   nkdm->ops->solution       = kdm->ops->solution;
169   nkdm->ops->destroy        = kdm->ops->destroy;
170   nkdm->ops->duplicate      = kdm->ops->duplicate;
171 
172   nkdm->rhsfunctionctx      = kdm->rhsfunctionctx;
173   nkdm->rhsjacobianctx      = kdm->rhsjacobianctx;
174   nkdm->ifunctionctx        = kdm->ifunctionctx;
175   nkdm->ijacobianctx        = kdm->ijacobianctx;
176   nkdm->solutionctx         = kdm->solutionctx;
177 
178   nkdm->data             = kdm->data;
179 
180   /*
181   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
182   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
183   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
184   */
185 
186   /* implementation specific copy hooks */
187   if (kdm->ops->duplicate) {ierr = (*kdm->ops->duplicate)(kdm,nkdm);CHKERRQ(ierr);}
188   PetscFunctionReturn(0);
189 }
190 
191 #undef __FUNCT__
192 #define __FUNCT__ "DMGetDMTS"
193 /*@C
194    DMGetDMTS - get read-only private DMTS context from a DM
195 
196    Not Collective
197 
198    Input Argument:
199 .  dm - DM to be used with TS
200 
201    Output Argument:
202 .  tsdm - private DMTS context
203 
204    Level: developer
205 
206    Notes:
207    Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.
208 
209 .seealso: DMGetDMTSWrite()
210 @*/
211 PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
212 {
213   PetscErrorCode ierr;
214 
215   PetscFunctionBegin;
216   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
217   *tsdm = (DMTS) dm->dmts;
218   if (!*tsdm) {
219     ierr = PetscInfo(dm,"Creating new DMTS\n");CHKERRQ(ierr);
220     ierr = DMTSCreate(((PetscObject)dm)->comm,tsdm);CHKERRQ(ierr);
221     dm->dmts = (PetscObject) *tsdm;
222     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr);
223     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr);
224   }
225   PetscFunctionReturn(0);
226 }
227 
228 #undef __FUNCT__
229 #define __FUNCT__ "DMGetDMTSWrite"
230 /*@C
231    DMGetDMTSWrite - get write access to private DMTS context from a DM
232 
233    Not Collective
234 
235    Input Argument:
236 .  dm - DM to be used with TS
237 
238    Output Argument:
239 .  tsdm - private DMTS context
240 
241    Level: developer
242 
243 .seealso: DMGetDMTS()
244 @*/
245 PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
246 {
247   PetscErrorCode ierr;
248   DMTS           sdm;
249 
250   PetscFunctionBegin;
251   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
252   ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr);
253   if (!sdm->originaldm) sdm->originaldm = dm;
254   if (sdm->originaldm != dm) {  /* Copy on write */
255     DMTS          oldsdm = sdm;
256     ierr = PetscInfo(dm,"Copying DMTS due to write\n");CHKERRQ(ierr);
257     ierr = DMTSCreate(((PetscObject)dm)->comm,&sdm);CHKERRQ(ierr);
258     ierr = DMTSCopy(oldsdm,sdm);CHKERRQ(ierr);
259     ierr = DMTSDestroy((DMTS*)&dm->dmts);CHKERRQ(ierr);
260     dm->dmts = (PetscObject) sdm;
261   }
262   *tsdm = sdm;
263   PetscFunctionReturn(0);
264 }
265 
266 #undef __FUNCT__
267 #define __FUNCT__ "DMCopyDMTS"
268 /*@C
269    DMCopyDMTS - copies a DM context to a new DM
270 
271    Logically Collective
272 
273    Input Arguments:
274 +  dmsrc - DM to obtain context from
275 -  dmdest - DM to add context to
276 
277    Level: developer
278 
279    Note:
280    The context is copied by reference. This function does not ensure that a context exists.
281 
282 .seealso: DMGetDMTS(), TSSetDM()
283 @*/
284 PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
285 {
286   PetscErrorCode ierr;
287 
288   PetscFunctionBegin;
289   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
290   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
291   ierr = DMTSDestroy((DMTS*)&dmdest->dmts);CHKERRQ(ierr);
292   dmdest->dmts = dmsrc->dmts;
293   ierr = PetscObjectReference(dmdest->dmts);CHKERRQ(ierr);
294   ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr);
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "DMTSSetIFunction"
300 /*@C
301    DMTSSetIFunction - set TS implicit function evaluation function
302 
303    Not Collective
304 
305    Input Arguments:
306 +  dm - DM to be used with TS
307 .  func - function evaluation function, see TSSetIFunction() for calling sequence
308 -  ctx - context for residual evaluation
309 
310    Level: advanced
311 
312    Note:
313    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
314    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
315    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
316 
317 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
318 @*/
319 PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
320 {
321   PetscErrorCode ierr;
322   DMTS           tsdm;
323 
324   PetscFunctionBegin;
325   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
326   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
327   if (func) tsdm->ops->ifunction = func;
328   if (ctx)  tsdm->ifunctionctx = ctx;
329   PetscFunctionReturn(0);
330 }
331 
332 #undef __FUNCT__
333 #define __FUNCT__ "DMTSGetIFunction"
334 /*@C
335    DMTSGetIFunction - get TS implicit residual evaluation function
336 
337    Not Collective
338 
339    Input Argument:
340 .  dm - DM to be used with TS
341 
342    Output Arguments:
343 +  func - function evaluation function, see TSSetIFunction() for calling sequence
344 -  ctx - context for residual evaluation
345 
346    Level: advanced
347 
348    Note:
349    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
350    associated with the DM.
351 
352 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
353 @*/
354 PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
355 {
356   PetscErrorCode ierr;
357   DMTS           tsdm;
358 
359   PetscFunctionBegin;
360   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
361   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
362   if (func) *func = tsdm->ops->ifunction;
363   if (ctx)  *ctx = tsdm->ifunctionctx;
364   PetscFunctionReturn(0);
365 }
366 
367 
368 #undef __FUNCT__
369 #define __FUNCT__ "DMTSSetRHSFunction"
370 /*@C
371    DMTSSetRHSFunction - set TS explicit residual evaluation function
372 
373    Not Collective
374 
375    Input Arguments:
376 +  dm - DM to be used with TS
377 .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
378 -  ctx - context for residual evaluation
379 
380    Level: advanced
381 
382    Note:
383    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
384    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
385    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
386 
387 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
388 @*/
389 PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
390 {
391   PetscErrorCode ierr;
392   DMTS           tsdm;
393 
394   PetscFunctionBegin;
395   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
396   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
397   if (func) tsdm->ops->rhsfunction = func;
398   if (ctx)  tsdm->rhsfunctionctx = ctx;
399   PetscFunctionReturn(0);
400 }
401 
402 #undef __FUNCT__
403 #define __FUNCT__ "DMTSGetSolutionFunction"
404 /*@C
405    DMTSGetSolutionFunction - gets the TS solution evaluation function
406 
407    Not Collective
408 
409    Input Arguments:
410 .  dm - DM to be used with TS
411 
412    Output Parameters:
413 +  func - solution function evaluation function, see TSSetSolution() for calling sequence
414 -  ctx - context for solution evaluation
415 
416    Level: advanced
417 
418 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
419 @*/
420 PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
421 {
422   PetscErrorCode ierr;
423   DMTS           tsdm;
424 
425   PetscFunctionBegin;
426   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
427   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
428   if (func) *func = tsdm->ops->solution;
429   if (ctx)  *ctx  = tsdm->solutionctx;
430   PetscFunctionReturn(0);
431 }
432 
433 #undef __FUNCT__
434 #define __FUNCT__ "DMTSSetSolutionFunction"
435 /*@C
436    DMTSSetSolutionFunction - set TS solution evaluation function
437 
438    Not Collective
439 
440    Input Arguments:
441 +  dm - DM to be used with TS
442 .  func - solution function evaluation function, see TSSetSolution() for calling sequence
443 -  ctx - context for solution evaluation
444 
445    Level: advanced
446 
447    Note:
448    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
449    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
450    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
451 
452 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
453 @*/
454 PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
455 {
456   PetscErrorCode ierr;
457   DMTS           tsdm;
458 
459   PetscFunctionBegin;
460   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
461   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
462   if (func) tsdm->ops->solution = func;
463   if (ctx)  tsdm->solutionctx   = ctx;
464   PetscFunctionReturn(0);
465 }
466 
467 #undef __FUNCT__
468 #define __FUNCT__ "DMTSGetRHSFunction"
469 /*@C
470    DMTSGetRHSFunction - get TS explicit residual evaluation function
471 
472    Not Collective
473 
474    Input Argument:
475 .  dm - DM to be used with TS
476 
477    Output Arguments:
478 +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
479 -  ctx - context for residual evaluation
480 
481    Level: advanced
482 
483    Note:
484    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
485    associated with the DM.
486 
487 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
488 @*/
489 PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
490 {
491   PetscErrorCode ierr;
492   DMTS           tsdm;
493 
494   PetscFunctionBegin;
495   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
496   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
497   if (func) *func = tsdm->ops->rhsfunction;
498   if (ctx)  *ctx = tsdm->rhsfunctionctx;
499   PetscFunctionReturn(0);
500 }
501 
502 #undef __FUNCT__
503 #define __FUNCT__ "DMTSSetIJacobian"
504 /*@C
505    DMTSSetIJacobian - set TS Jacobian evaluation function
506 
507    Not Collective
508 
509    Input Argument:
510 +  dm - DM to be used with TS
511 .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
512 -  ctx - context for residual evaluation
513 
514    Level: advanced
515 
516    Note:
517    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
518    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
519    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
520 
521 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
522 @*/
523 PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
524 {
525   PetscErrorCode ierr;
526   DMTS           sdm;
527 
528   PetscFunctionBegin;
529   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
530   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
531   if (func) sdm->ops->ijacobian = func;
532   if (ctx)  sdm->ijacobianctx   = ctx;
533   PetscFunctionReturn(0);
534 }
535 
536 #undef __FUNCT__
537 #define __FUNCT__ "DMTSGetIJacobian"
538 /*@C
539    DMTSGetIJacobian - get TS Jacobian evaluation function
540 
541    Not Collective
542 
543    Input Argument:
544 .  dm - DM to be used with TS
545 
546    Output Arguments:
547 +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
548 -  ctx - context for residual evaluation
549 
550    Level: advanced
551 
552    Note:
553    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
554    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
555    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
556 
557 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
558 @*/
559 PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
560 {
561   PetscErrorCode ierr;
562   DMTS           tsdm;
563 
564   PetscFunctionBegin;
565   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
566   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
567   if (func) *func = tsdm->ops->ijacobian;
568   if (ctx)  *ctx = tsdm->ijacobianctx;
569   PetscFunctionReturn(0);
570 }
571 
572 
573 #undef __FUNCT__
574 #define __FUNCT__ "DMTSSetRHSJacobian"
575 /*@C
576    DMTSSetRHSJacobian - set TS Jacobian evaluation function
577 
578    Not Collective
579 
580    Input Argument:
581 +  dm - DM to be used with TS
582 .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
583 -  ctx - context for residual evaluation
584 
585    Level: advanced
586 
587    Note:
588    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
589    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
590    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
591 
592 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
593 @*/
594 PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
595 {
596   PetscErrorCode ierr;
597   DMTS           tsdm;
598 
599   PetscFunctionBegin;
600   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
601   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
602   if (func) tsdm->ops->rhsjacobian = func;
603   if (ctx)  tsdm->rhsjacobianctx = ctx;
604   PetscFunctionReturn(0);
605 }
606 
607 #undef __FUNCT__
608 #define __FUNCT__ "DMTSGetRHSJacobian"
609 /*@C
610    DMTSGetRHSJacobian - get TS Jacobian evaluation function
611 
612    Not Collective
613 
614    Input Argument:
615 .  dm - DM to be used with TS
616 
617    Output Arguments:
618 +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
619 -  ctx - context for residual evaluation
620 
621    Level: advanced
622 
623    Note:
624    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
625    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
626    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
627 
628 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
629 @*/
630 PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
631 {
632   PetscErrorCode ierr;
633   DMTS           tsdm;
634 
635   PetscFunctionBegin;
636   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
637   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
638   if (func) *func = tsdm->ops->rhsjacobian;
639   if (ctx)  *ctx = tsdm->rhsjacobianctx;
640   PetscFunctionReturn(0);
641 }
642 
643 #undef __FUNCT__
644 #define __FUNCT__ "DMTSSetIFunctionSerialize"
645 /*@C
646    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
647 
648    Not Collective
649 
650    Input Arguments:
651 +  dm - DM to be used with TS
652 .  view - viewer function
653 -  load - loading function
654 
655    Level: advanced
656 
657 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
658 @*/
659 PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
660 {
661   PetscErrorCode ierr;
662   DMTS           tsdm;
663 
664   PetscFunctionBegin;
665   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
666   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
667   tsdm->ops->ifunctionview = view;
668   tsdm->ops->ifunctionload = load;
669   PetscFunctionReturn(0);
670 }
671 
672 #undef __FUNCT__
673 #define __FUNCT__ "DMTSSetIJacobianSerialize"
674 /*@C
675    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
676 
677    Not Collective
678 
679    Input Arguments:
680 +  dm - DM to be used with TS
681 .  view - viewer function
682 -  load - loading function
683 
684    Level: advanced
685 
686 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
687 @*/
688 PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
689 {
690   PetscErrorCode ierr;
691   DMTS           tsdm;
692 
693   PetscFunctionBegin;
694   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
695   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
696   tsdm->ops->ijacobianview = view;
697   tsdm->ops->ijacobianload = load;
698   PetscFunctionReturn(0);
699 }
700