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