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