xref: /petsc/src/ts/utils/dmts.c (revision 89669be4d29968dc8d4c19ce1b69194a6a561ea4)
1 #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
2 #include <petsc/private/dmimpl.h>
3 
4 static PetscErrorCode DMTSDestroy(DMTS *kdm)
5 {
6   PetscFunctionBegin;
7   if (!*kdm) PetscFunctionReturn(0);
8   PetscValidHeaderSpecific((*kdm),DMTS_CLASSID,1);
9   if (--((PetscObject)(*kdm))->refct > 0) {*kdm = NULL; PetscFunctionReturn(0);}
10   if ((*kdm)->ops->destroy) PetscCall(((*kdm)->ops->destroy)(*kdm));
11   PetscCall(PetscHeaderDestroy(kdm));
12   PetscFunctionReturn(0);
13 }
14 
15 PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer)
16 {
17   PetscFunctionBegin;
18   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,NULL,PETSC_FUNCTION));
19   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,NULL,PETSC_FUNCTION));
20   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,NULL,PETSC_FUNCTION));
21   if (kdm->ops->ifunctionload) {
22     PetscCall((*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer));
23   }
24   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,NULL,PETSC_FUNCTION));
25   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,NULL,PETSC_FUNCTION));
26   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,NULL,PETSC_FUNCTION));
27   if (kdm->ops->ijacobianload) {
28     PetscCall((*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer));
29   }
30   PetscFunctionReturn(0);
31 }
32 
33 PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer)
34 {
35   PetscBool      isascii,isbinary;
36 
37   PetscFunctionBegin;
38   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
39   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
40   if (isascii) {
41 #if defined(PETSC_SERIALIZE_FUNCTIONS)
42     const char *fname;
43 
44     PetscCall(PetscFPTFind(kdm->ops->ifunction,&fname));
45     if (fname) {
46       PetscCall(PetscViewerASCIIPrintf(viewer,"  IFunction used by TS: %s\n",fname));
47     }
48     PetscCall(PetscFPTFind(kdm->ops->ijacobian,&fname));
49     if (fname) {
50       PetscCall(PetscViewerASCIIPrintf(viewer,"  IJacobian function used by TS: %s\n",fname));
51     }
52 #endif
53   } else if (isbinary) {
54     struct {
55       TSIFunction ifunction;
56     } funcstruct;
57     struct {
58       PetscErrorCode (*ifunctionview)(void*,PetscViewer);
59     } funcviewstruct;
60     struct {
61       PetscErrorCode (*ifunctionload)(void**,PetscViewer);
62     } funcloadstruct;
63     struct {
64       TSIJacobian ijacobian;
65     } jacstruct;
66     struct {
67       PetscErrorCode (*ijacobianview)(void*,PetscViewer);
68     } jacviewstruct;
69     struct {
70       PetscErrorCode (*ijacobianload)(void**,PetscViewer);
71     } jacloadstruct;
72 
73     funcstruct.ifunction         = kdm->ops->ifunction;
74     funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
75     funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
76     PetscCall(PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION));
77     PetscCall(PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION));
78     PetscCall(PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION));
79     if (kdm->ops->ifunctionview) {
80       PetscCall((*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer));
81     }
82     jacstruct.ijacobian = kdm->ops->ijacobian;
83     jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
84     jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
85     PetscCall(PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION));
86     PetscCall(PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION));
87     PetscCall(PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION));
88     if (kdm->ops->ijacobianview) {
89       PetscCall((*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer));
90     }
91   }
92   PetscFunctionReturn(0);
93 }
94 
95 static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
96 {
97   PetscFunctionBegin;
98   PetscCall(TSInitializePackage());
99   PetscCall(PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView));
100   PetscFunctionReturn(0);
101 }
102 
103 /* Attaches the DMTS to the coarse level.
104  * Under what conditions should we copy versus duplicate?
105  */
106 static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
107 {
108   PetscFunctionBegin;
109   PetscCall(DMCopyDMTS(dm,dmc));
110   PetscFunctionReturn(0);
111 }
112 
113 /* This could restrict auxiliary information to the coarse level.
114  */
115 static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
116 {
117   PetscFunctionBegin;
118   PetscFunctionReturn(0);
119 }
120 
121 static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
122 {
123   PetscFunctionBegin;
124   PetscCall(DMCopyDMTS(dm,subdm));
125   PetscFunctionReturn(0);
126 }
127 
128 /* This could restrict auxiliary information to the coarse level.
129  */
130 static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
131 {
132   PetscFunctionBegin;
133   PetscFunctionReturn(0);
134 }
135 
136 /*@C
137    DMTSCopy - copies the information in a DMTS to another DMTS
138 
139    Not Collective
140 
141    Input Parameters:
142 +  kdm - Original DMTS
143 -  nkdm - DMTS to receive the data, should have been created with DMTSCreate()
144 
145    Level: developer
146 
147 .seealso: `DMTSCreate()`, `DMTSDestroy()`
148 @*/
149 PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
150 {
151   PetscFunctionBegin;
152   PetscValidHeaderSpecific(kdm,DMTS_CLASSID,1);
153   PetscValidHeaderSpecific(nkdm,DMTS_CLASSID,2);
154   nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
155   nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
156   nkdm->ops->ifunction   = kdm->ops->ifunction;
157   nkdm->ops->ijacobian   = kdm->ops->ijacobian;
158   nkdm->ops->i2function  = kdm->ops->i2function;
159   nkdm->ops->i2jacobian  = kdm->ops->i2jacobian;
160   nkdm->ops->solution    = kdm->ops->solution;
161   nkdm->ops->destroy     = kdm->ops->destroy;
162   nkdm->ops->duplicate   = kdm->ops->duplicate;
163 
164   nkdm->rhsfunctionctx = kdm->rhsfunctionctx;
165   nkdm->rhsjacobianctx = kdm->rhsjacobianctx;
166   nkdm->ifunctionctx   = kdm->ifunctionctx;
167   nkdm->ijacobianctx   = kdm->ijacobianctx;
168   nkdm->i2functionctx  = kdm->i2functionctx;
169   nkdm->i2jacobianctx  = kdm->i2jacobianctx;
170   nkdm->solutionctx    = kdm->solutionctx;
171 
172   nkdm->data = kdm->data;
173 
174   /*
175   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
176   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
177   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
178   */
179 
180   /* implementation specific copy hooks */
181   if (kdm->ops->duplicate) PetscCall((*kdm->ops->duplicate)(kdm,nkdm));
182   PetscFunctionReturn(0);
183 }
184 
185 /*@C
186    DMGetDMTS - get read-only private DMTS context from a DM
187 
188    Not Collective
189 
190    Input Parameter:
191 .  dm - DM to be used with TS
192 
193    Output Parameter:
194 .  tsdm - private DMTS context
195 
196    Level: developer
197 
198    Notes:
199    Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.
200 
201 .seealso: `DMGetDMTSWrite()`
202 @*/
203 PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
204 {
205   PetscFunctionBegin;
206   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
207   *tsdm = (DMTS) dm->dmts;
208   if (!*tsdm) {
209     PetscCall(PetscInfo(dm,"Creating new DMTS\n"));
210     PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm));
211     dm->dmts = (PetscObject) *tsdm;
212     (*tsdm)->originaldm = dm;
213     PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL));
214     PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL));
215   }
216   PetscFunctionReturn(0);
217 }
218 
219 /*@C
220    DMGetDMTSWrite - get write access to private DMTS context from a DM
221 
222    Not Collective
223 
224    Input Parameter:
225 .  dm - DM to be used with TS
226 
227    Output Parameter:
228 .  tsdm - private DMTS context
229 
230    Level: developer
231 
232 .seealso: `DMGetDMTS()`
233 @*/
234 PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
235 {
236   DMTS           sdm;
237 
238   PetscFunctionBegin;
239   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
240   PetscCall(DMGetDMTS(dm,&sdm));
241   PetscCheck(sdm->originaldm,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DMTS has a NULL originaldm");
242   if (sdm->originaldm != dm) {  /* Copy on write */
243     DMTS oldsdm = sdm;
244     PetscCall(PetscInfo(dm,"Copying DMTS due to write\n"));
245     PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm));
246     PetscCall(DMTSCopy(oldsdm,sdm));
247     PetscCall(DMTSDestroy((DMTS*)&dm->dmts));
248     dm->dmts = (PetscObject) sdm;
249     sdm->originaldm = dm;
250   }
251   *tsdm = sdm;
252   PetscFunctionReturn(0);
253 }
254 
255 /*@C
256    DMCopyDMTS - copies a DM context to a new DM
257 
258    Logically Collective
259 
260    Input Parameters:
261 +  dmsrc - DM to obtain context from
262 -  dmdest - DM to add context to
263 
264    Level: developer
265 
266    Note:
267    The context is copied by reference. This function does not ensure that a context exists.
268 
269 .seealso: `DMGetDMTS()`, `TSSetDM()`
270 @*/
271 PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
272 {
273   PetscFunctionBegin;
274   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
275   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
276   PetscCall(DMTSDestroy((DMTS*)&dmdest->dmts));
277   dmdest->dmts = dmsrc->dmts;
278   PetscCall(PetscObjectReference(dmdest->dmts));
279   PetscCall(DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL));
280   PetscCall(DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL));
281   PetscFunctionReturn(0);
282 }
283 
284 /*@C
285    DMTSSetIFunction - set TS implicit function evaluation function
286 
287    Not Collective
288 
289    Input Parameters:
290 +  dm - DM to be used with TS
291 .  func - function evaluating f(t,u,u_t)
292 -  ctx - context for residual evaluation
293 
294    Calling sequence of func:
295 $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
296 
297 +  t   - time at step/stage being solved
298 .  u   - state vector
299 .  u_t - time derivative of state vector
300 .  F   - function vector
301 -  ctx - [optional] user-defined context for matrix evaluation routine
302 
303    Level: advanced
304 
305    Note:
306    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
307    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
308    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
309 
310 .seealso: `DMTSSetContext()`, `TSSetIFunction()`, `DMTSSetJacobian()`
311 @*/
312 PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
313 {
314   DMTS           tsdm;
315 
316   PetscFunctionBegin;
317   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
318   PetscCall(DMGetDMTSWrite(dm,&tsdm));
319   if (func) tsdm->ops->ifunction = func;
320   if (ctx)  tsdm->ifunctionctx = ctx;
321   PetscFunctionReturn(0);
322 }
323 
324 /*@C
325    DMTSGetIFunction - get TS implicit residual evaluation function
326 
327    Not Collective
328 
329    Input Parameter:
330 .  dm - DM to be used with TS
331 
332    Output Parameters:
333 +  func - function evaluation function, see TSSetIFunction() for calling sequence
334 -  ctx - context for residual evaluation
335 
336    Level: advanced
337 
338    Note:
339    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
340    associated with the DM.
341 
342 .seealso: `DMTSSetContext()`, `DMTSSetFunction()`, `TSSetFunction()`
343 @*/
344 PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
345 {
346   DMTS           tsdm;
347 
348   PetscFunctionBegin;
349   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
350   PetscCall(DMGetDMTS(dm,&tsdm));
351   if (func) *func = tsdm->ops->ifunction;
352   if (ctx)  *ctx = tsdm->ifunctionctx;
353   PetscFunctionReturn(0);
354 }
355 
356 /*@C
357    DMTSSetI2Function - set TS implicit function evaluation function for 2nd order systems
358 
359    Not Collective
360 
361    Input Parameters:
362 +  dm - DM to be used with TS
363 .  fun - function evaluation routine
364 -  ctx - context for residual evaluation
365 
366    Calling sequence of fun:
367 $     PetscErrorCode fun(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,Vec F,ctx);
368 
369 +  t    - time at step/stage being solved
370 .  U    - state vector
371 .  U_t  - time derivative of state vector
372 .  U_tt - second time derivative of state vector
373 .  F    - function vector
374 -  ctx  - [optional] user-defined context for matrix evaluation routine (may be NULL)
375 
376    Level: advanced
377 
378    Note:
379    TSSetI2Function() is normally used, but it calls this function internally because the user context is actually
380    associated with the DM.
381 
382 .seealso: `TSSetI2Function()`
383 @*/
384 PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx)
385 {
386   DMTS           tsdm;
387 
388   PetscFunctionBegin;
389   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
390   PetscCall(DMGetDMTSWrite(dm,&tsdm));
391   if (fun) tsdm->ops->i2function = fun;
392   if (ctx) tsdm->i2functionctx   = ctx;
393   PetscFunctionReturn(0);
394 }
395 
396 /*@C
397    DMTSGetI2Function - get TS implicit residual evaluation function for 2nd order systems
398 
399    Not Collective
400 
401    Input Parameter:
402 .  dm - DM to be used with TS
403 
404    Output Parameters:
405 +  fun - function evaluation function, see TSSetI2Function() for calling sequence
406 -  ctx - context for residual evaluation
407 
408    Level: advanced
409 
410    Note:
411    TSGetI2Function() is normally used, but it calls this function internally because the user context is actually
412    associated with the DM.
413 
414 .seealso: `DMTSSetI2Function(),TSGetI2Function()`
415 @*/
416 PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx)
417 {
418   DMTS           tsdm;
419 
420   PetscFunctionBegin;
421   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
422   PetscCall(DMGetDMTS(dm,&tsdm));
423   if (fun) *fun = tsdm->ops->i2function;
424   if (ctx) *ctx = tsdm->i2functionctx;
425   PetscFunctionReturn(0);
426 }
427 
428 /*@C
429    DMTSSetI2Jacobian - set TS implicit Jacobian evaluation function for 2nd order systems
430 
431    Not Collective
432 
433    Input Parameters:
434 +  dm - DM to be used with TS
435 .  fun - Jacobian evaluation routine
436 -  ctx - context for Jacobian evaluation
437 
438    Calling sequence of jac:
439 $    PetscErrorCode jac(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,PetscReal v,PetscReal a,Mat J,Mat P,void *ctx);
440 
441 +  t    - time at step/stage being solved
442 .  U    - state vector
443 .  U_t  - time derivative of state vector
444 .  U_tt - second time derivative of state vector
445 .  v    - shift for U_t
446 .  a    - shift for U_tt
447 .  J    - Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t  + a*dF/dU_tt
448 .  P    - preconditioning matrix for J, may be same as J
449 -  ctx  - [optional] user-defined context for matrix evaluation routine
450 
451    Level: advanced
452 
453    Note:
454    TSSetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
455    associated with the DM.
456 
457 .seealso: `TSSetI2Jacobian()`
458 @*/
459 PetscErrorCode DMTSSetI2Jacobian(DM dm,TSI2Jacobian jac,void *ctx)
460 {
461   DMTS           tsdm;
462 
463   PetscFunctionBegin;
464   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
465   PetscCall(DMGetDMTSWrite(dm,&tsdm));
466   if (jac) tsdm->ops->i2jacobian = jac;
467   if (ctx) tsdm->i2jacobianctx   = ctx;
468   PetscFunctionReturn(0);
469 }
470 
471 /*@C
472    DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems
473 
474    Not Collective
475 
476    Input Parameter:
477 .  dm - DM to be used with TS
478 
479    Output Parameters:
480 +  jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
481 -  ctx - context for Jacobian evaluation
482 
483    Level: advanced
484 
485    Note:
486    TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
487    associated with the DM.
488 
489 .seealso: `DMTSSetI2Jacobian(),TSGetI2Jacobian()`
490 @*/
491 PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx)
492 {
493   DMTS           tsdm;
494 
495   PetscFunctionBegin;
496   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
497   PetscCall(DMGetDMTS(dm,&tsdm));
498   if (jac) *jac = tsdm->ops->i2jacobian;
499   if (ctx) *ctx = tsdm->i2jacobianctx;
500   PetscFunctionReturn(0);
501 }
502 
503 /*@C
504    DMTSSetRHSFunction - set TS explicit residual evaluation function
505 
506    Not Collective
507 
508    Input Parameters:
509 +  dm - DM to be used with TS
510 .  func - RHS function evaluation routine
511 -  ctx - context for residual evaluation
512 
513     Calling sequence of func:
514 $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec F,void *ctx);
515 
516 +   ts - timestep context
517 .   t - current timestep
518 .   u - input vector
519 .   F - function vector
520 -   ctx - [optional] user-defined function context
521 
522    Level: advanced
523 
524    Note:
525    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
526    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
527    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
528 
529 .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetJacobian()`
530 @*/
531 PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
532 {
533   DMTS           tsdm;
534 
535   PetscFunctionBegin;
536   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
537   PetscCall(DMGetDMTSWrite(dm,&tsdm));
538   if (func) tsdm->ops->rhsfunction = func;
539   if (ctx)  tsdm->rhsfunctionctx = ctx;
540   PetscFunctionReturn(0);
541 }
542 
543 /*@C
544    DMTSSetTransientVariable - sets function to transform from state to transient variables
545 
546    Logically Collective
547 
548    Input Parameters:
549 +  dm - DM to be used with TS
550 .  tvar - a function that transforms to transient variables
551 -  ctx - a context for tvar
552 
553     Calling sequence of tvar:
554 $     PetscErrorCode tvar(TS ts,Vec p,Vec c,void *ctx);
555 
556 +   ts - timestep context
557 .   p - input vector (primative form)
558 .   c - output vector, transient variables (conservative form)
559 -   ctx - [optional] user-defined function context
560 
561    Level: advanced
562 
563    Notes:
564    This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., TSBDF)
565    can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
566    well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
567    C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
568    evaluated via the chain rule, as in
569 
570      dF/dP + shift * dF/dCdot dC/dP.
571 
572 .seealso: `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()`
573 @*/
574 PetscErrorCode DMTSSetTransientVariable(DM dm,TSTransientVariable tvar,void *ctx)
575 {
576   DMTS           dmts;
577 
578   PetscFunctionBegin;
579   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
580   PetscCall(DMGetDMTSWrite(dm,&dmts));
581   dmts->ops->transientvar = tvar;
582   dmts->transientvarctx = ctx;
583   PetscFunctionReturn(0);
584 }
585 
586 /*@C
587    DMTSGetTransientVariable - gets function to transform from state to transient variables
588 
589    Logically Collective
590 
591    Input Parameter:
592 .  dm - DM to be used with TS
593 
594    Output Parameters:
595 +  tvar - a function that transforms to transient variables
596 -  ctx - a context for tvar
597 
598    Level: advanced
599 
600 .seealso: `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()`
601 @*/
602 PetscErrorCode DMTSGetTransientVariable(DM dm,TSTransientVariable *tvar,void *ctx)
603 {
604   DMTS           dmts;
605 
606   PetscFunctionBegin;
607   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
608   PetscCall(DMGetDMTS(dm,&dmts));
609   if (tvar) *tvar = dmts->ops->transientvar;
610   if (ctx)  *(void**)ctx = dmts->transientvarctx;
611   PetscFunctionReturn(0);
612 }
613 
614 /*@C
615    DMTSGetSolutionFunction - gets the TS solution evaluation function
616 
617    Not Collective
618 
619    Input Parameter:
620 .  dm - DM to be used with TS
621 
622    Output Parameters:
623 +  func - solution function evaluation function, see TSSetSolution() for calling sequence
624 -  ctx - context for solution evaluation
625 
626    Level: advanced
627 
628 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSSetSolutionFunction()`
629 @*/
630 PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
631 {
632   DMTS           tsdm;
633 
634   PetscFunctionBegin;
635   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
636   PetscCall(DMGetDMTS(dm,&tsdm));
637   if (func) *func = tsdm->ops->solution;
638   if (ctx)  *ctx  = tsdm->solutionctx;
639   PetscFunctionReturn(0);
640 }
641 
642 /*@C
643    DMTSSetSolutionFunction - set TS solution evaluation function
644 
645    Not Collective
646 
647    Input Parameters:
648 +  dm - DM to be used with TS
649 .  func - solution function evaluation routine
650 -  ctx - context for solution evaluation
651 
652     Calling sequence of f:
653 $     PetscErrorCode f(TS ts,PetscReal t,Vec u,void *ctx);
654 
655 +   ts - timestep context
656 .   t - current timestep
657 .   u - output vector
658 -   ctx - [optional] user-defined function context
659 
660    Level: advanced
661 
662    Note:
663    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
664    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
665    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
666 
667 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSGetSolutionFunction()`
668 @*/
669 PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
670 {
671   DMTS           tsdm;
672 
673   PetscFunctionBegin;
674   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
675   PetscCall(DMGetDMTSWrite(dm,&tsdm));
676   if (func) tsdm->ops->solution = func;
677   if (ctx)  tsdm->solutionctx   = ctx;
678   PetscFunctionReturn(0);
679 }
680 
681 /*@C
682    DMTSSetForcingFunction - set TS forcing function evaluation function
683 
684    Not Collective
685 
686    Input Parameters:
687 +  dm - DM to be used with TS
688 .  f - forcing function evaluation routine
689 -  ctx - context for solution evaluation
690 
691     Calling sequence of func:
692 $     PetscErrorCode func (TS ts,PetscReal t,Vec f,void *ctx);
693 
694 +   ts - timestep context
695 .   t - current timestep
696 .   f - output vector
697 -   ctx - [optional] user-defined function context
698 
699    Level: advanced
700 
701    Note:
702    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
703    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
704    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
705 
706 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
707 @*/
708 PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx)
709 {
710   DMTS           tsdm;
711 
712   PetscFunctionBegin;
713   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
714   PetscCall(DMGetDMTSWrite(dm,&tsdm));
715   if (f)    tsdm->ops->forcing = f;
716   if (ctx)  tsdm->forcingctx   = ctx;
717   PetscFunctionReturn(0);
718 }
719 
720 /*@C
721    DMTSGetForcingFunction - get TS forcing function evaluation function
722 
723    Not Collective
724 
725    Input Parameter:
726 .   dm - DM to be used with TS
727 
728    Output Parameters:
729 +  f - forcing function evaluation function; see TSForcingFunction for details
730 -  ctx - context for solution evaluation
731 
732    Level: advanced
733 
734    Note:
735    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
736    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
737    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
738 
739 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
740 @*/
741 PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
742 {
743   DMTS           tsdm;
744 
745   PetscFunctionBegin;
746   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
747   PetscCall(DMGetDMTSWrite(dm,&tsdm));
748   if (f)   *f   = tsdm->ops->forcing;
749   if (ctx) *ctx = tsdm->forcingctx;
750   PetscFunctionReturn(0);
751 }
752 
753 /*@C
754    DMTSGetRHSFunction - get TS explicit residual evaluation function
755 
756    Not Collective
757 
758    Input Parameter:
759 .  dm - DM to be used with TS
760 
761    Output Parameters:
762 +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
763 -  ctx - context for residual evaluation
764 
765    Level: advanced
766 
767    Note:
768    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
769    associated with the DM.
770 
771 .seealso: `DMTSSetContext()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
772 @*/
773 PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
774 {
775   DMTS           tsdm;
776 
777   PetscFunctionBegin;
778   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
779   PetscCall(DMGetDMTS(dm,&tsdm));
780   if (func) *func = tsdm->ops->rhsfunction;
781   if (ctx)  *ctx = tsdm->rhsfunctionctx;
782   PetscFunctionReturn(0);
783 }
784 
785 /*@C
786    DMTSSetIJacobian - set TS Jacobian evaluation function
787 
788    Not Collective
789 
790    Input Parameters:
791 +  dm - DM to be used with TS
792 .  func - Jacobian evaluation routine
793 -  ctx - context for residual evaluation
794 
795    Calling sequence of f:
796 $    PetscErrorCode f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx);
797 
798 +  t    - time at step/stage being solved
799 .  U    - state vector
800 .  U_t  - time derivative of state vector
801 .  a    - shift
802 .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
803 .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
804 -  ctx  - [optional] user-defined context for matrix evaluation routine
805 
806    Level: advanced
807 
808    Note:
809    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
810    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
811    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
812 
813 .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSGetJacobian()`, `TSSetIJacobian()`, `TSSetIFunction()`
814 @*/
815 PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
816 {
817   DMTS           sdm;
818 
819   PetscFunctionBegin;
820   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
821   PetscCall(DMGetDMTSWrite(dm,&sdm));
822   if (func) sdm->ops->ijacobian = func;
823   if (ctx)  sdm->ijacobianctx   = ctx;
824   PetscFunctionReturn(0);
825 }
826 
827 /*@C
828    DMTSGetIJacobian - get TS Jacobian evaluation function
829 
830    Not Collective
831 
832    Input Parameter:
833 .  dm - DM to be used with TS
834 
835    Output Parameters:
836 +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
837 -  ctx - context for residual evaluation
838 
839    Level: advanced
840 
841    Note:
842    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
843    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
844    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
845 
846 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
847 @*/
848 PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
849 {
850   DMTS           tsdm;
851 
852   PetscFunctionBegin;
853   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
854   PetscCall(DMGetDMTS(dm,&tsdm));
855   if (func) *func = tsdm->ops->ijacobian;
856   if (ctx)  *ctx = tsdm->ijacobianctx;
857   PetscFunctionReturn(0);
858 }
859 
860 /*@C
861    DMTSSetRHSJacobian - set TS Jacobian evaluation function
862 
863    Not Collective
864 
865    Input Parameters:
866 +  dm - DM to be used with TS
867 .  func - Jacobian evaluation routine
868 -  ctx - context for residual evaluation
869 
870    Calling sequence of func:
871 $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);
872 
873 +  t - current timestep
874 .  u - input vector
875 .  Amat - (approximate) Jacobian matrix
876 .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
877 -  ctx - [optional] user-defined context for matrix evaluation routine
878 
879    Level: advanced
880 
881    Note:
882    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
883    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
884    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
885 
886 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSGetJacobian()`, `TSSetRHSJacobian()`
887 @*/
888 PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
889 {
890   DMTS           tsdm;
891 
892   PetscFunctionBegin;
893   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
894   PetscCall(DMGetDMTSWrite(dm,&tsdm));
895   if (func) tsdm->ops->rhsjacobian = func;
896   if (ctx)  tsdm->rhsjacobianctx = ctx;
897   PetscFunctionReturn(0);
898 }
899 
900 /*@C
901    DMTSGetRHSJacobian - get TS Jacobian evaluation function
902 
903    Not Collective
904 
905    Input Parameter:
906 .  dm - DM to be used with TS
907 
908    Output Parameters:
909 +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
910 -  ctx - context for residual evaluation
911 
912    Level: advanced
913 
914    Note:
915    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
916    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
917    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
918 
919 .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
920 @*/
921 PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
922 {
923   DMTS           tsdm;
924 
925   PetscFunctionBegin;
926   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
927   PetscCall(DMGetDMTS(dm,&tsdm));
928   if (func) *func = tsdm->ops->rhsjacobian;
929   if (ctx)  *ctx = tsdm->rhsjacobianctx;
930   PetscFunctionReturn(0);
931 }
932 
933 /*@C
934    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
935 
936    Not Collective
937 
938    Input Parameters:
939 +  dm - DM to be used with TS
940 .  view - viewer function
941 -  load - loading function
942 
943    Level: advanced
944 
945 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
946 @*/
947 PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
948 {
949   DMTS           tsdm;
950 
951   PetscFunctionBegin;
952   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
953   PetscCall(DMGetDMTSWrite(dm,&tsdm));
954   tsdm->ops->ifunctionview = view;
955   tsdm->ops->ifunctionload = load;
956   PetscFunctionReturn(0);
957 }
958 
959 /*@C
960    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
961 
962    Not Collective
963 
964    Input Parameters:
965 +  dm - DM to be used with TS
966 .  view - viewer function
967 -  load - loading function
968 
969    Level: advanced
970 
971 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
972 @*/
973 PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
974 {
975   DMTS           tsdm;
976 
977   PetscFunctionBegin;
978   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
979   PetscCall(DMGetDMTSWrite(dm,&tsdm));
980   tsdm->ops->ijacobianview = view;
981   tsdm->ops->ijacobianload = load;
982   PetscFunctionReturn(0);
983 }
984