xref: /petsc/src/ts/utils/dmts.c (revision fbf9dbe564678ed6eff1806adbc4c4f01b9743f4)
1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
2 #include <petsc/private/dmimpl.h>
3 
4 static PetscErrorCode DMTSUnsetRHSFunctionContext_DMTS(DMTS tsdm)
5 {
6   PetscFunctionBegin;
7   PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs function ctx", NULL));
8   tsdm->rhsfunctionctxcontainer = NULL;
9   PetscFunctionReturn(PETSC_SUCCESS);
10 }
11 
12 static PetscErrorCode DMTSUnsetRHSJacobianContext_DMTS(DMTS tsdm)
13 {
14   PetscFunctionBegin;
15   PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs jacobian ctx", NULL));
16   tsdm->rhsjacobianctxcontainer = NULL;
17   PetscFunctionReturn(PETSC_SUCCESS);
18 }
19 
20 static PetscErrorCode DMTSUnsetIFunctionContext_DMTS(DMTS tsdm)
21 {
22   PetscFunctionBegin;
23   PetscCall(PetscObjectCompose((PetscObject)tsdm, "ifunction ctx", NULL));
24   tsdm->ifunctionctxcontainer = NULL;
25   PetscFunctionReturn(PETSC_SUCCESS);
26 }
27 
28 static PetscErrorCode DMTSUnsetIJacobianContext_DMTS(DMTS tsdm)
29 {
30   PetscFunctionBegin;
31   PetscCall(PetscObjectCompose((PetscObject)tsdm, "ijacobian ctx", NULL));
32   tsdm->ijacobianctxcontainer = NULL;
33   PetscFunctionReturn(PETSC_SUCCESS);
34 }
35 
36 static PetscErrorCode DMTSUnsetI2FunctionContext_DMTS(DMTS tsdm)
37 {
38   PetscFunctionBegin;
39   PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2function ctx", NULL));
40   tsdm->i2functionctxcontainer = NULL;
41   PetscFunctionReturn(PETSC_SUCCESS);
42 }
43 
44 static PetscErrorCode DMTSUnsetI2JacobianContext_DMTS(DMTS tsdm)
45 {
46   PetscFunctionBegin;
47   PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2jacobian ctx", NULL));
48   tsdm->i2jacobianctxcontainer = NULL;
49   PetscFunctionReturn(PETSC_SUCCESS);
50 }
51 
52 static PetscErrorCode DMTSDestroy(DMTS *kdm)
53 {
54   PetscFunctionBegin;
55   if (!*kdm) PetscFunctionReturn(PETSC_SUCCESS);
56   PetscValidHeaderSpecific((*kdm), DMTS_CLASSID, 1);
57   if (--((PetscObject)(*kdm))->refct > 0) {
58     *kdm = NULL;
59     PetscFunctionReturn(PETSC_SUCCESS);
60   }
61   PetscCall(DMTSUnsetRHSFunctionContext_DMTS(*kdm));
62   PetscCall(DMTSUnsetRHSJacobianContext_DMTS(*kdm));
63   PetscCall(DMTSUnsetIFunctionContext_DMTS(*kdm));
64   PetscCall(DMTSUnsetIJacobianContext_DMTS(*kdm));
65   PetscCall(DMTSUnsetI2FunctionContext_DMTS(*kdm));
66   PetscCall(DMTSUnsetI2JacobianContext_DMTS(*kdm));
67   PetscTryTypeMethod(*kdm, destroy);
68   PetscCall(PetscHeaderDestroy(kdm));
69   PetscFunctionReturn(PETSC_SUCCESS);
70 }
71 
72 PetscErrorCode DMTSLoad(DMTS kdm, PetscViewer viewer)
73 {
74   PetscFunctionBegin;
75   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunction, 1, NULL, PETSC_FUNCTION));
76   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunctionview, 1, NULL, PETSC_FUNCTION));
77   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunctionload, 1, NULL, PETSC_FUNCTION));
78   if (kdm->ops->ifunctionload) {
79     void *ctx;
80 
81     PetscCall(PetscContainerGetPointer(kdm->ifunctionctxcontainer, &ctx));
82     PetscCall((*kdm->ops->ifunctionload)(&ctx, viewer));
83   }
84   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobian, 1, NULL, PETSC_FUNCTION));
85   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobianview, 1, NULL, PETSC_FUNCTION));
86   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobianload, 1, NULL, PETSC_FUNCTION));
87   if (kdm->ops->ijacobianload) {
88     void *ctx;
89 
90     PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer, &ctx));
91     PetscCall((*kdm->ops->ijacobianload)(&ctx, viewer));
92   }
93   PetscFunctionReturn(PETSC_SUCCESS);
94 }
95 
96 PetscErrorCode DMTSView(DMTS kdm, PetscViewer viewer)
97 {
98   PetscBool isascii, isbinary;
99 
100   PetscFunctionBegin;
101   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
102   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
103   if (isascii) {
104 #if defined(PETSC_SERIALIZE_FUNCTIONS)
105     const char *fname;
106 
107     PetscCall(PetscFPTFind(kdm->ops->ifunction, &fname));
108     if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, "  IFunction used by TS: %s\n", fname));
109     PetscCall(PetscFPTFind(kdm->ops->ijacobian, &fname));
110     if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, "  IJacobian function used by TS: %s\n", fname));
111 #endif
112   } else if (isbinary) {
113     struct {
114       TSIFunction ifunction;
115     } funcstruct;
116     struct {
117       PetscErrorCode (*ifunctionview)(void *, PetscViewer);
118     } funcviewstruct;
119     struct {
120       PetscErrorCode (*ifunctionload)(void **, PetscViewer);
121     } funcloadstruct;
122     struct {
123       TSIJacobian ijacobian;
124     } jacstruct;
125     struct {
126       PetscErrorCode (*ijacobianview)(void *, PetscViewer);
127     } jacviewstruct;
128     struct {
129       PetscErrorCode (*ijacobianload)(void **, PetscViewer);
130     } jacloadstruct;
131 
132     funcstruct.ifunction         = kdm->ops->ifunction;
133     funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
134     funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
135     PetscCall(PetscViewerBinaryWrite(viewer, &funcstruct, 1, PETSC_FUNCTION));
136     PetscCall(PetscViewerBinaryWrite(viewer, &funcviewstruct, 1, PETSC_FUNCTION));
137     PetscCall(PetscViewerBinaryWrite(viewer, &funcloadstruct, 1, PETSC_FUNCTION));
138     if (kdm->ops->ifunctionview) {
139       void *ctx;
140 
141       PetscCall(PetscContainerGetPointer(kdm->ifunctionctxcontainer, &ctx));
142       PetscCall((*kdm->ops->ifunctionview)(ctx, viewer));
143     }
144     jacstruct.ijacobian         = kdm->ops->ijacobian;
145     jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
146     jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
147     PetscCall(PetscViewerBinaryWrite(viewer, &jacstruct, 1, PETSC_FUNCTION));
148     PetscCall(PetscViewerBinaryWrite(viewer, &jacviewstruct, 1, PETSC_FUNCTION));
149     PetscCall(PetscViewerBinaryWrite(viewer, &jacloadstruct, 1, PETSC_FUNCTION));
150     if (kdm->ops->ijacobianview) {
151       void *ctx;
152 
153       PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer, &ctx));
154       PetscCall((*kdm->ops->ijacobianview)(ctx, viewer));
155     }
156   }
157   PetscFunctionReturn(PETSC_SUCCESS);
158 }
159 
160 static PetscErrorCode DMTSCreate(MPI_Comm comm, DMTS *kdm)
161 {
162   PetscFunctionBegin;
163   PetscCall(TSInitializePackage());
164   PetscCall(PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView));
165   PetscFunctionReturn(PETSC_SUCCESS);
166 }
167 
168 /* Attaches the DMTS to the coarse level.
169  * Under what conditions should we copy versus duplicate?
170  */
171 static PetscErrorCode DMCoarsenHook_DMTS(DM dm, DM dmc, void *ctx)
172 {
173   PetscFunctionBegin;
174   PetscCall(DMCopyDMTS(dm, dmc));
175   PetscFunctionReturn(PETSC_SUCCESS);
176 }
177 
178 /* This could restrict auxiliary information to the coarse level.
179  */
180 static PetscErrorCode DMRestrictHook_DMTS(DM dm, Mat Restrict, Vec rscale, Mat Inject, DM dmc, void *ctx)
181 {
182   PetscFunctionBegin;
183   PetscFunctionReturn(PETSC_SUCCESS);
184 }
185 
186 static PetscErrorCode DMSubDomainHook_DMTS(DM dm, DM subdm, void *ctx)
187 {
188   PetscFunctionBegin;
189   PetscCall(DMCopyDMTS(dm, subdm));
190   PetscFunctionReturn(PETSC_SUCCESS);
191 }
192 
193 /* This could restrict auxiliary information to the coarse level.
194  */
195 static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
196 {
197   PetscFunctionBegin;
198   PetscFunctionReturn(PETSC_SUCCESS);
199 }
200 
201 /*@C
202    DMTSCopy - copies the information in a `DMTS` to another `DMTS`
203 
204    Not Collective
205 
206    Input Parameters:
207 +  kdm - Original `DMTS`
208 -  nkdm - `DMTS` to receive the data, should have been created with `DMTSCreate()`
209 
210    Level: developer
211 
212 .seealso: [](ch_ts), `DMTSCreate()`, `DMTSDestroy()`
213 @*/
214 PetscErrorCode DMTSCopy(DMTS kdm, DMTS nkdm)
215 {
216   PetscFunctionBegin;
217   PetscValidHeaderSpecific(kdm, DMTS_CLASSID, 1);
218   PetscValidHeaderSpecific(nkdm, DMTS_CLASSID, 2);
219   nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
220   nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
221   nkdm->ops->ifunction   = kdm->ops->ifunction;
222   nkdm->ops->ijacobian   = kdm->ops->ijacobian;
223   nkdm->ops->i2function  = kdm->ops->i2function;
224   nkdm->ops->i2jacobian  = kdm->ops->i2jacobian;
225   nkdm->ops->solution    = kdm->ops->solution;
226   nkdm->ops->destroy     = kdm->ops->destroy;
227   nkdm->ops->duplicate   = kdm->ops->duplicate;
228 
229   nkdm->solutionctx             = kdm->solutionctx;
230   nkdm->rhsfunctionctxcontainer = kdm->rhsfunctionctxcontainer;
231   nkdm->rhsjacobianctxcontainer = kdm->rhsjacobianctxcontainer;
232   nkdm->ifunctionctxcontainer   = kdm->ifunctionctxcontainer;
233   nkdm->ijacobianctxcontainer   = kdm->ijacobianctxcontainer;
234   nkdm->i2functionctxcontainer  = kdm->i2functionctxcontainer;
235   nkdm->i2jacobianctxcontainer  = kdm->i2jacobianctxcontainer;
236   if (nkdm->rhsfunctionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "rhs function ctx", (PetscObject)nkdm->rhsfunctionctxcontainer));
237   if (nkdm->rhsjacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "rhs jacobian ctx", (PetscObject)nkdm->rhsjacobianctxcontainer));
238   if (nkdm->ifunctionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "ifunction ctx", (PetscObject)nkdm->ifunctionctxcontainer));
239   if (nkdm->ijacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "ijacobian ctx", (PetscObject)nkdm->ijacobianctxcontainer));
240   if (nkdm->i2functionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "i2function ctx", (PetscObject)nkdm->i2functionctxcontainer));
241   if (nkdm->i2jacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "i2jacobian ctx", (PetscObject)nkdm->i2jacobianctxcontainer));
242 
243   nkdm->data = kdm->data;
244 
245   /*
246   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
247   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
248   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
249   */
250 
251   /* implementation specific copy hooks */
252   PetscTryTypeMethod(kdm, duplicate, nkdm);
253   PetscFunctionReturn(PETSC_SUCCESS);
254 }
255 
256 /*@C
257    DMGetDMTS - get read-only private `DMTS` context from a `DM`
258 
259    Not Collective
260 
261    Input Parameter:
262 .  dm - `DM` to be used with `TS`
263 
264    Output Parameter:
265 .  tsdm - private `DMTS` context
266 
267    Level: developer
268 
269    Notes:
270    Use `DMGetDMTSWrite()` if write access is needed. The `DMTSSetXXX()` API should be used wherever possible.
271 
272 .seealso: [](ch_ts), `DMTS`, `DMGetDMTSWrite()`
273 @*/
274 PetscErrorCode DMGetDMTS(DM dm, DMTS *tsdm)
275 {
276   PetscFunctionBegin;
277   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
278   *tsdm = (DMTS)dm->dmts;
279   if (!*tsdm) {
280     PetscCall(PetscInfo(dm, "Creating new DMTS\n"));
281     PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm), tsdm));
282     dm->dmts            = (PetscObject)*tsdm;
283     (*tsdm)->originaldm = dm;
284     PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_DMTS, DMRestrictHook_DMTS, NULL));
285     PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_DMTS, DMSubDomainRestrictHook_DMTS, NULL));
286   }
287   PetscFunctionReturn(PETSC_SUCCESS);
288 }
289 
290 /*@C
291    DMGetDMTSWrite - get write access to private `DMTS` context from a `DM`
292 
293    Not Collective
294 
295    Input Parameter:
296 .  dm - `DM` to be used with `TS`
297 
298    Output Parameter:
299 .  tsdm - private `DMTS` context
300 
301    Level: developer
302 
303 .seealso: [](ch_ts), `DMTS`, `DMGetDMTS()`
304 @*/
305 PetscErrorCode DMGetDMTSWrite(DM dm, DMTS *tsdm)
306 {
307   DMTS sdm;
308 
309   PetscFunctionBegin;
310   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
311   PetscCall(DMGetDMTS(dm, &sdm));
312   PetscCheck(sdm->originaldm, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DMTS has a NULL originaldm");
313   if (sdm->originaldm != dm) { /* Copy on write */
314     DMTS oldsdm = sdm;
315     PetscCall(PetscInfo(dm, "Copying DMTS due to write\n"));
316     PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm), &sdm));
317     PetscCall(DMTSCopy(oldsdm, sdm));
318     PetscCall(DMTSDestroy((DMTS *)&dm->dmts));
319     dm->dmts        = (PetscObject)sdm;
320     sdm->originaldm = dm;
321   }
322   *tsdm = sdm;
323   PetscFunctionReturn(PETSC_SUCCESS);
324 }
325 
326 /*@C
327    DMCopyDMTS - copies a `DM` context to a new `DM`
328 
329    Logically Collective
330 
331    Input Parameters:
332 +  dmsrc - `DM` to obtain context from
333 -  dmdest - `DM` to add context to
334 
335    Level: developer
336 
337    Note:
338    The context is copied by reference. This function does not ensure that a context exists.
339 
340 .seealso: [](ch_ts), `DMTS`, `DMGetDMTS()`, `TSSetDM()`
341 @*/
342 PetscErrorCode DMCopyDMTS(DM dmsrc, DM dmdest)
343 {
344   PetscFunctionBegin;
345   PetscValidHeaderSpecific(dmsrc, DM_CLASSID, 1);
346   PetscValidHeaderSpecific(dmdest, DM_CLASSID, 2);
347   PetscCall(DMTSDestroy((DMTS *)&dmdest->dmts));
348   dmdest->dmts = dmsrc->dmts;
349   PetscCall(PetscObjectReference(dmdest->dmts));
350   PetscCall(DMCoarsenHookAdd(dmdest, DMCoarsenHook_DMTS, DMRestrictHook_DMTS, NULL));
351   PetscCall(DMSubDomainHookAdd(dmdest, DMSubDomainHook_DMTS, DMSubDomainRestrictHook_DMTS, NULL));
352   PetscFunctionReturn(PETSC_SUCCESS);
353 }
354 
355 /*@C
356    DMTSSetIFunction - set `TS` implicit function evaluation function
357 
358    Not Collective
359 
360    Input Parameters:
361 +  dm - `DM` to be used with `TS`
362 .  func - function evaluating f(t,u,u_t)
363 -  ctx - context for residual evaluation
364 
365    Calling sequence of `func`:
366 $     PetscErrorCode func(TS ts, PetscReal t, Vec u, Vec u_t, Vec F,void *ctx);
367 +  ts  - the `TS` context obtained from `TSCreate()`
368 .  t   - time at step/stage being solved
369 .  u   - state vector
370 .  u_t - time derivative of state vector
371 .  F   - function vector
372 -  ctx - [optional] user-defined context for matrix evaluation routine
373 
374    Level: advanced
375 
376    Note:
377    `TSSetFunction()` is normally used, but it calls this function internally because the user context is actually
378    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
379    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
380 
381 .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `TSSetIFunction()`, `DMTSSetJacobian()`
382 @*/
383 PetscErrorCode DMTSSetIFunction(DM dm, TSIFunction func, void *ctx)
384 {
385   DMTS tsdm;
386 
387   PetscFunctionBegin;
388   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
389   PetscCall(DMGetDMTSWrite(dm, &tsdm));
390   if (func) tsdm->ops->ifunction = func;
391   if (ctx) {
392     PetscContainer ctxcontainer;
393     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
394     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
395     PetscCall(PetscObjectCompose((PetscObject)tsdm, "ifunction ctx", (PetscObject)ctxcontainer));
396     tsdm->ifunctionctxcontainer = ctxcontainer;
397     PetscCall(PetscContainerDestroy(&ctxcontainer));
398   }
399   PetscFunctionReturn(PETSC_SUCCESS);
400 }
401 
402 /*@C
403    DMTSSetIFunctionContextDestroy - set `TS` implicit evaluation context destroy function
404 
405    Not Collective
406 
407    Input Parameters:
408 +  dm - `DM` to be used with `TS`
409 -  f - implicit evaluation context destroy function
410 
411    Level: advanced
412 
413 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetIFunction()`, `TSSetIFunction()`
414 @*/
415 PetscErrorCode DMTSSetIFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
416 {
417   DMTS tsdm;
418 
419   PetscFunctionBegin;
420   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
421   PetscCall(DMGetDMTSWrite(dm, &tsdm));
422   if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ifunctionctxcontainer, f));
423   PetscFunctionReturn(PETSC_SUCCESS);
424 }
425 
426 PetscErrorCode DMTSUnsetIFunctionContext_Internal(DM dm)
427 {
428   DMTS tsdm;
429 
430   PetscFunctionBegin;
431   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
432   PetscCall(DMGetDMTSWrite(dm, &tsdm));
433   PetscCall(DMTSUnsetIFunctionContext_DMTS(tsdm));
434   PetscFunctionReturn(PETSC_SUCCESS);
435 }
436 
437 /*@C
438    DMTSGetIFunction - get `TS` implicit residual evaluation function
439 
440    Not Collective
441 
442    Input Parameter:
443 .  dm - `DM` to be used with `TS`
444 
445    Output Parameters:
446 +  func - function evaluation function, for calling sequence see `TSSetIFunction()`
447 -  ctx - context for residual evaluation
448 
449    Level: advanced
450 
451    Note:
452    `TSGetFunction()` is normally used, but it calls this function internally because the user context is actually
453    associated with the `DM`.
454 
455 .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `DMTSSetFunction()`, `TSSetFunction()`
456 @*/
457 PetscErrorCode DMTSGetIFunction(DM dm, TSIFunction *func, void **ctx)
458 {
459   DMTS tsdm;
460 
461   PetscFunctionBegin;
462   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
463   PetscCall(DMGetDMTS(dm, &tsdm));
464   if (func) *func = tsdm->ops->ifunction;
465   if (ctx) {
466     if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ifunctionctxcontainer, ctx));
467     else *ctx = NULL;
468   }
469   PetscFunctionReturn(PETSC_SUCCESS);
470 }
471 
472 /*@C
473    DMTSSetI2Function - set `TS` implicit function evaluation function for 2nd order systems
474 
475    Not Collective
476 
477    Input Parameters:
478 +  dm - `DM` to be used with `TS`
479 .  fun - function evaluation routine
480 -  ctx - context for residual evaluation
481 
482    Calling sequence of `fun`:
483 $  PetscErrorCode fun(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F,void *ctx);
484 +  ts  - the `TS` context obtained from `TSCreate()`
485 .  t    - time at step/stage being solved
486 .  U    - state vector
487 .  U_t  - time derivative of state vector
488 .  U_tt - second time derivative of state vector
489 .  F    - function vector
490 -  ctx  - [optional] user-defined context for matrix evaluation routine (may be `NULL`)
491 
492    Level: advanced
493 
494    Note:
495    `TSSetI2Function()` is normally used, but it calls this function internally because the user context is actually
496    associated with the `DM`.
497 
498 .seealso: [](ch_ts), `DM`, `TS`, `TSSetI2Function()`
499 @*/
500 PetscErrorCode DMTSSetI2Function(DM dm, TSI2Function fun, void *ctx)
501 {
502   DMTS tsdm;
503 
504   PetscFunctionBegin;
505   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
506   PetscCall(DMGetDMTSWrite(dm, &tsdm));
507   if (fun) tsdm->ops->i2function = fun;
508   if (ctx) {
509     PetscContainer ctxcontainer;
510     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
511     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
512     PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2function ctx", (PetscObject)ctxcontainer));
513     tsdm->i2functionctxcontainer = ctxcontainer;
514     PetscCall(PetscContainerDestroy(&ctxcontainer));
515   }
516   PetscFunctionReturn(PETSC_SUCCESS);
517 }
518 
519 /*@C
520    DMTSSetI2FunctionContextDestroy - set `TS` implicit evaluation for 2nd order systems context destroy
521 
522    Not Collective
523 
524    Input Parameters:
525 +  dm - `DM` to be used with `TS`
526 -  f - implicit evaluation context destroy function
527 
528    Level: advanced
529 
530    Note:
531    `TSSetI2FunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
532    associated with the `DM`.
533 
534 .seealso: [](ch_ts), `TSSetI2FunctionContextDestroy()`, `DMTSSetI2Function()`, `TSSetI2Function()`
535 @*/
536 PetscErrorCode DMTSSetI2FunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
537 {
538   DMTS tsdm;
539 
540   PetscFunctionBegin;
541   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
542   PetscCall(DMGetDMTSWrite(dm, &tsdm));
543   if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2functionctxcontainer, f));
544   PetscFunctionReturn(PETSC_SUCCESS);
545 }
546 
547 PetscErrorCode DMTSUnsetI2FunctionContext_Internal(DM dm)
548 {
549   DMTS tsdm;
550 
551   PetscFunctionBegin;
552   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
553   PetscCall(DMGetDMTSWrite(dm, &tsdm));
554   PetscCall(DMTSUnsetI2FunctionContext_DMTS(tsdm));
555   PetscFunctionReturn(PETSC_SUCCESS);
556 }
557 
558 /*@C
559    DMTSGetI2Function - get `TS` implicit residual evaluation function for 2nd order systems
560 
561    Not Collective
562 
563    Input Parameter:
564 .  dm - `DM` to be used with `TS`
565 
566    Output Parameters:
567 +  fun - function evaluation function, for calling sequence see `TSSetI2Function()`
568 -  ctx - context for residual evaluation
569 
570    Level: advanced
571 
572    Note:
573    `TSGetI2Function()` is normally used, but it calls this function internally because the user context is actually
574    associated with the `DM`.
575 
576 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetI2Function()`, `TSGetI2Function()`
577 @*/
578 PetscErrorCode DMTSGetI2Function(DM dm, TSI2Function *fun, void **ctx)
579 {
580   DMTS tsdm;
581 
582   PetscFunctionBegin;
583   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
584   PetscCall(DMGetDMTS(dm, &tsdm));
585   if (fun) *fun = tsdm->ops->i2function;
586   if (ctx) {
587     if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2functionctxcontainer, ctx));
588     else *ctx = NULL;
589   }
590   PetscFunctionReturn(PETSC_SUCCESS);
591 }
592 
593 /*@C
594    DMTSSetI2Jacobian - set `TS` implicit Jacobian evaluation function for 2nd order systems
595 
596    Not Collective
597 
598    Input Parameters:
599 +  dm - `DM` to be used with `TS`
600 .  fun - Jacobian evaluation routine
601 -  ctx - context for Jacobian evaluation
602 
603    Calling sequence of `jac`:
604 $    PetscErrorCode jac(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, PetscReal v, PetscReal a, Mat J, Mat P, void *ctx)
605 +  ts  - the `TS` context obtained from `TSCreate()`
606 .  t    - time at step/stage being solved
607 .  U    - state vector
608 .  U_t  - time derivative of state vector
609 .  U_tt - second time derivative of state vector
610 .  v    - shift for U_t
611 .  a    - shift for U_tt
612 .  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
613 .  P    - preconditioning matrix for J, may be same as `J`
614 -  ctx  - [optional] user-defined context for matrix evaluation routine
615 
616    Level: advanced
617 
618    Note:
619    `TSSetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
620    associated with the `DM`.
621 
622 .seealso: [](ch_ts), `DM`, `TS`, `TSSetI2Jacobian()`
623 @*/
624 PetscErrorCode DMTSSetI2Jacobian(DM dm, TSI2Jacobian jac, void *ctx)
625 {
626   DMTS tsdm;
627 
628   PetscFunctionBegin;
629   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
630   PetscCall(DMGetDMTSWrite(dm, &tsdm));
631   if (jac) tsdm->ops->i2jacobian = jac;
632   if (ctx) {
633     PetscContainer ctxcontainer;
634     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
635     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
636     PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2jacobian ctx", (PetscObject)ctxcontainer));
637     tsdm->i2jacobianctxcontainer = ctxcontainer;
638     PetscCall(PetscContainerDestroy(&ctxcontainer));
639   }
640   PetscFunctionReturn(PETSC_SUCCESS);
641 }
642 
643 /*@C
644    DMTSSetI2JacobianContextDestroy - set `TS` implicit Jacobian evaluation for 2nd order systems context destroy function
645 
646    Not Collective
647 
648    Input Parameters:
649 +  dm - `DM` to be used with `TS`
650 -  f - implicit Jacobian evaluation context destroy function
651 
652    Level: advanced
653 
654 .seealso: [](ch_ts), `DM`, `TS`, `TSSetI2JacobianContextDestroy()`, `DMTSSetI2Jacobian()`, `TSSetI2Jacobian()`
655 @*/
656 PetscErrorCode DMTSSetI2JacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
657 {
658   DMTS tsdm;
659 
660   PetscFunctionBegin;
661   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
662   PetscCall(DMGetDMTSWrite(dm, &tsdm));
663   if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2jacobianctxcontainer, f));
664   PetscFunctionReturn(PETSC_SUCCESS);
665 }
666 
667 PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM dm)
668 {
669   DMTS tsdm;
670 
671   PetscFunctionBegin;
672   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
673   PetscCall(DMGetDMTSWrite(dm, &tsdm));
674   PetscCall(DMTSUnsetI2JacobianContext_DMTS(tsdm));
675   PetscFunctionReturn(PETSC_SUCCESS);
676 }
677 
678 /*@C
679    DMTSGetI2Jacobian - get `TS` implicit Jacobian evaluation function for 2nd order systems
680 
681    Not Collective
682 
683    Input Parameter:
684 .  dm - `DM` to be used with `TS`
685 
686    Output Parameters:
687 +  jac - Jacobian evaluation function,  for calling sequence see `TSSetI2Jacobian()`
688 -  ctx - context for Jacobian evaluation
689 
690    Level: advanced
691 
692    Note:
693    `TSGetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
694    associated with the `DM`.
695 
696 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetI2Jacobian()`, `TSGetI2Jacobian()`
697 @*/
698 PetscErrorCode DMTSGetI2Jacobian(DM dm, TSI2Jacobian *jac, void **ctx)
699 {
700   DMTS tsdm;
701 
702   PetscFunctionBegin;
703   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
704   PetscCall(DMGetDMTS(dm, &tsdm));
705   if (jac) *jac = tsdm->ops->i2jacobian;
706   if (ctx) {
707     if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2jacobianctxcontainer, ctx));
708     else *ctx = NULL;
709   }
710   PetscFunctionReturn(PETSC_SUCCESS);
711 }
712 
713 /*@C
714    DMTSSetRHSFunction - set `TS` explicit residual evaluation function
715 
716    Not Collective
717 
718    Input Parameters:
719 +  dm - `DM` to be used with `TS`
720 .  func - RHS function evaluation routine
721 -  ctx - context for residual evaluation
722 
723     Calling sequence of `func`:
724 $   PetscErrorCode func(TS ts, PetscReal t, Vec u, Vec F, void *ctx);
725 +   ts - timestep context
726 .   t - current timestep
727 .   u - input vector
728 .   F - function vector
729 -   ctx - [optional] user-defined function context
730 
731    Level: advanced
732 
733    Note:
734    `TSSetRHSFunction()` is normally used, but it calls this function internally because the user context is actually
735    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
736    not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
737 
738 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetJacobian()`
739 @*/
740 PetscErrorCode DMTSSetRHSFunction(DM dm, TSRHSFunction func, void *ctx)
741 {
742   DMTS tsdm;
743 
744   PetscFunctionBegin;
745   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
746   PetscCall(DMGetDMTSWrite(dm, &tsdm));
747   if (func) tsdm->ops->rhsfunction = func;
748   if (ctx) {
749     PetscContainer ctxcontainer;
750     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
751     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
752     PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs function ctx", (PetscObject)ctxcontainer));
753     tsdm->rhsfunctionctxcontainer = ctxcontainer;
754     PetscCall(PetscContainerDestroy(&ctxcontainer));
755   }
756   PetscFunctionReturn(PETSC_SUCCESS);
757 }
758 
759 /*@C
760    DMTSSetRHSFunctionContextDestroy - set `TS` explicit residual evaluation context destroy function
761 
762    Not Collective
763 
764    Input Parameters:
765 +  dm - `DM` to be used with `TS`
766 -  f - explicit evaluation context destroy function
767 
768    Level: advanced
769 
770    Note:
771    `TSSetRHSFunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
772    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
773    not.
774 
775    Developer Note:
776    If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
777 
778 .seealso: [](ch_ts), `TSSetRHSFunctionContextDestroy()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
779 @*/
780 PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
781 {
782   DMTS tsdm;
783 
784   PetscFunctionBegin;
785   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
786   PetscCall(DMGetDMTSWrite(dm, &tsdm));
787   if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsfunctionctxcontainer, f));
788   PetscFunctionReturn(PETSC_SUCCESS);
789 }
790 
791 PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM dm)
792 {
793   DMTS tsdm;
794 
795   PetscFunctionBegin;
796   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
797   PetscCall(DMGetDMTSWrite(dm, &tsdm));
798   PetscCall(DMTSUnsetRHSFunctionContext_DMTS(tsdm));
799   tsdm->rhsfunctionctxcontainer = NULL;
800   PetscFunctionReturn(PETSC_SUCCESS);
801 }
802 
803 /*@C
804    DMTSSetTransientVariable - sets function to transform from state to transient variables
805 
806    Logically Collective
807 
808    Input Parameters:
809 +  dm - `DM` to be used with `TS`
810 .  tvar - a function that transforms to transient variables
811 -  ctx - a context for tvar
812 
813     Calling sequence of `tvar`:
814 $   PetscErrorCode tvar(TS ts, Vec p, Vec c, void *ctx);
815 +   ts - timestep context
816 .   p - input vector (primitive form)
817 .   c - output vector, transient variables (conservative form)
818 -   ctx - [optional] user-defined function context
819 
820    Level: advanced
821 
822    Notes:
823    This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., `TSBDF`)
824    can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
825    well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
826    C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
827    evaluated via the chain rule, as in
828 
829      dF/dP + shift * dF/dCdot dC/dP.
830 
831 .seealso: [](ch_ts), `TS`, `TSBDF`, `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()`
832 @*/
833 PetscErrorCode DMTSSetTransientVariable(DM dm, TSTransientVariable tvar, void *ctx)
834 {
835   DMTS dmts;
836 
837   PetscFunctionBegin;
838   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
839   PetscCall(DMGetDMTSWrite(dm, &dmts));
840   dmts->ops->transientvar = tvar;
841   dmts->transientvarctx   = ctx;
842   PetscFunctionReturn(PETSC_SUCCESS);
843 }
844 
845 /*@C
846    DMTSGetTransientVariable - gets function to transform from state to transient variables set with `DMTSSetTransientVariable()`
847 
848    Logically Collective
849 
850    Input Parameter:
851 .  dm - `DM` to be used with `TS`
852 
853    Output Parameters:
854 +  tvar - a function that transforms to transient variables
855 -  ctx - a context for tvar
856 
857    Level: advanced
858 
859 .seealso: [](ch_ts), `DM`, `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()`
860 @*/
861 PetscErrorCode DMTSGetTransientVariable(DM dm, TSTransientVariable *tvar, void *ctx)
862 {
863   DMTS dmts;
864 
865   PetscFunctionBegin;
866   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
867   PetscCall(DMGetDMTS(dm, &dmts));
868   if (tvar) *tvar = dmts->ops->transientvar;
869   if (ctx) *(void **)ctx = dmts->transientvarctx;
870   PetscFunctionReturn(PETSC_SUCCESS);
871 }
872 
873 /*@C
874    DMTSGetSolutionFunction - gets the `TS` solution evaluation function
875 
876    Not Collective
877 
878    Input Parameter:
879 .  dm - `DM` to be used with `TS`
880 
881    Output Parameters:
882 +  func - solution function evaluation function, for calling sequence see `TSSetSolution()`
883 -  ctx - context for solution evaluation
884 
885    Level: advanced
886 
887 .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSSetSolutionFunction()`
888 @*/
889 PetscErrorCode DMTSGetSolutionFunction(DM dm, TSSolutionFunction *func, void **ctx)
890 {
891   DMTS tsdm;
892 
893   PetscFunctionBegin;
894   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
895   PetscCall(DMGetDMTS(dm, &tsdm));
896   if (func) *func = tsdm->ops->solution;
897   if (ctx) *ctx = tsdm->solutionctx;
898   PetscFunctionReturn(PETSC_SUCCESS);
899 }
900 
901 /*@C
902    DMTSSetSolutionFunction - set `TS` solution evaluation function
903 
904    Not Collective
905 
906    Input Parameters:
907 +  dm - `DM` to be used with `TS`
908 .  func - solution function evaluation routine
909 -  ctx - context for solution evaluation
910 
911     Calling sequence of `f`:
912 $   PetscErrorCode f(TS ts, PetscReal t, Vec u, void *ctx);
913 +   ts - timestep context
914 .   t - current timestep
915 .   u - output vector
916 -   ctx - [optional] user-defined function context
917 
918    Level: advanced
919 
920    Note:
921    `TSSetSolutionFunction()` is normally used, but it calls this function internally because the user context is actually
922    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
923    not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
924 
925 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSGetSolutionFunction()`
926 @*/
927 PetscErrorCode DMTSSetSolutionFunction(DM dm, TSSolutionFunction func, void *ctx)
928 {
929   DMTS tsdm;
930 
931   PetscFunctionBegin;
932   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
933   PetscCall(DMGetDMTSWrite(dm, &tsdm));
934   if (func) tsdm->ops->solution = func;
935   if (ctx) tsdm->solutionctx = ctx;
936   PetscFunctionReturn(PETSC_SUCCESS);
937 }
938 
939 /*@C
940    DMTSSetForcingFunction - set `TS` forcing function evaluation function
941 
942    Not Collective
943 
944    Input Parameters:
945 +  dm - `DM` to be used with `TS`
946 .  func - forcing function evaluation routine
947 -  ctx - context for solution evaluation
948 
949     Calling sequence of `func`:
950 $     PetscErrorCode func (TS ts, PetscReal t, Vec f,void *ctx)
951 +   ts - timestep context
952 .   t - current timestep
953 .   f - output vector
954 -   ctx - [optional] user-defined function context
955 
956    Level: advanced
957 
958    Note:
959    `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually
960    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
961    not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
962 
963 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
964 @*/
965 PetscErrorCode DMTSSetForcingFunction(DM dm, TSForcingFunction func, void *ctx)
966 {
967   DMTS tsdm;
968 
969   PetscFunctionBegin;
970   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
971   PetscCall(DMGetDMTSWrite(dm, &tsdm));
972   if (func) tsdm->ops->forcing = func;
973   if (ctx) tsdm->forcingctx = ctx;
974   PetscFunctionReturn(PETSC_SUCCESS);
975 }
976 
977 /*@C
978    DMTSGetForcingFunction - get `TS` forcing function evaluation function
979 
980    Not Collective
981 
982    Input Parameter:
983 .   dm - `DM` to be used with `TS`
984 
985    Output Parameters:
986 +  f - forcing function evaluation function; see `TSForcingFunction` for details
987 -  ctx - context for solution evaluation
988 
989    Level: advanced
990 
991    Note:
992    `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually
993    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
994    not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
995 
996 .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
997 @*/
998 PetscErrorCode DMTSGetForcingFunction(DM dm, TSForcingFunction *f, void **ctx)
999 {
1000   DMTS tsdm;
1001 
1002   PetscFunctionBegin;
1003   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1004   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1005   if (f) *f = tsdm->ops->forcing;
1006   if (ctx) *ctx = tsdm->forcingctx;
1007   PetscFunctionReturn(PETSC_SUCCESS);
1008 }
1009 
1010 /*@C
1011    DMTSGetRHSFunction - get `TS` explicit residual evaluation function
1012 
1013    Not Collective
1014 
1015    Input Parameter:
1016 .  dm - `DM` to be used with `TS`
1017 
1018    Output Parameters:
1019 +  func - residual evaluation function, for calling sequence see `TSSetRHSFunction()`
1020 -  ctx - context for residual evaluation
1021 
1022    Level: advanced
1023 
1024    Note:
1025    `TSGetFunction()` is normally used, but it calls this function internally because the user context is actually
1026    associated with the DM.
1027 
1028 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
1029 @*/
1030 PetscErrorCode DMTSGetRHSFunction(DM dm, TSRHSFunction *func, void **ctx)
1031 {
1032   DMTS tsdm;
1033 
1034   PetscFunctionBegin;
1035   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1036   PetscCall(DMGetDMTS(dm, &tsdm));
1037   if (func) *func = tsdm->ops->rhsfunction;
1038   if (ctx) {
1039     if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsfunctionctxcontainer, ctx));
1040     else *ctx = NULL;
1041   }
1042   PetscFunctionReturn(PETSC_SUCCESS);
1043 }
1044 
1045 /*@C
1046    DMTSSetIJacobian - set `TS` Jacobian evaluation function
1047 
1048    Not Collective
1049 
1050    Input Parameters:
1051 +  dm - `DM` to be used with `TS`
1052 .  func - Jacobian evaluation routine
1053 -  ctx - context for residual evaluation
1054 
1055    Calling sequence of `f`:
1056 $    PetscErrorCode f(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx);
1057 +  ts  - the `TS` context obtained from `TSCreate()`
1058 .  t    - time at step/stage being solved
1059 .  U    - state vector
1060 .  U_t  - time derivative of state vector
1061 .  a    - shift
1062 .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
1063 .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
1064 -  ctx  - [optional] user-defined context for matrix evaluation routine
1065 
1066    Level: advanced
1067 
1068    Note:
1069    `TSSetJacobian()` is normally used, but it calls this function internally because the user context is actually
1070    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1071    not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1072 
1073 .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSGetJacobian()`, `TSSetIJacobian()`, `TSSetIFunction()`
1074 @*/
1075 PetscErrorCode DMTSSetIJacobian(DM dm, TSIJacobian func, void *ctx)
1076 {
1077   DMTS tsdm;
1078 
1079   PetscFunctionBegin;
1080   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1081   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1082   if (func) tsdm->ops->ijacobian = func;
1083   if (ctx) {
1084     PetscContainer ctxcontainer;
1085     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
1086     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
1087     PetscCall(PetscObjectCompose((PetscObject)tsdm, "ijacobian ctx", (PetscObject)ctxcontainer));
1088     tsdm->ijacobianctxcontainer = ctxcontainer;
1089     PetscCall(PetscContainerDestroy(&ctxcontainer));
1090   }
1091   PetscFunctionReturn(PETSC_SUCCESS);
1092 }
1093 
1094 /*@C
1095    DMTSSetIJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function
1096 
1097    Not Collective
1098 
1099    Input Parameters:
1100 +  dm - `DM` to be used with `TS`
1101 -  f - Jacobian evaluation context destroy function
1102 
1103    Level: advanced
1104 
1105    Note:
1106    `TSSetIJacobianContextDestroy()` is normally used, but it calls this function internally because the user context is actually
1107    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1108    not.
1109 
1110    Developer Note:
1111    If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1112 
1113 .seealso: [](ch_ts), `TSSetIJacobianContextDestroy()`, `TSSetI2JacobianContextDestroy()`, `DMTSSetIJacobian()`, `TSSetIJacobian()`
1114 @*/
1115 PetscErrorCode DMTSSetIJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
1116 {
1117   DMTS tsdm;
1118 
1119   PetscFunctionBegin;
1120   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1121   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1122   if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ijacobianctxcontainer, f));
1123   PetscFunctionReturn(PETSC_SUCCESS);
1124 }
1125 
1126 PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM dm)
1127 {
1128   DMTS tsdm;
1129 
1130   PetscFunctionBegin;
1131   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1132   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1133   PetscCall(DMTSUnsetIJacobianContext_DMTS(tsdm));
1134   PetscFunctionReturn(PETSC_SUCCESS);
1135 }
1136 
1137 /*@C
1138    DMTSGetIJacobian - get `TS` Jacobian evaluation function
1139 
1140    Not Collective
1141 
1142    Input Parameter:
1143 .  dm - `DM` to be used with `TS`
1144 
1145    Output Parameters:
1146 +  func - Jacobian evaluation function, for calling sequence see `TSSetIJacobian()`
1147 -  ctx - context for residual evaluation
1148 
1149    Level: advanced
1150 
1151    Note:
1152    `TSGetJacobian()` is normally used, but it calls this function internally because the user context is actually
1153    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1154    not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1155 
1156 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1157 @*/
1158 PetscErrorCode DMTSGetIJacobian(DM dm, TSIJacobian *func, void **ctx)
1159 {
1160   DMTS tsdm;
1161 
1162   PetscFunctionBegin;
1163   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1164   PetscCall(DMGetDMTS(dm, &tsdm));
1165   if (func) *func = tsdm->ops->ijacobian;
1166   if (ctx) {
1167     if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ijacobianctxcontainer, ctx));
1168     else *ctx = NULL;
1169   }
1170   PetscFunctionReturn(PETSC_SUCCESS);
1171 }
1172 
1173 /*@C
1174    DMTSSetRHSJacobian - set `TS` Jacobian evaluation function
1175 
1176    Not Collective
1177 
1178    Input Parameters:
1179 +  dm - `DM` to be used with `TS`
1180 .  func - Jacobian evaluation routine
1181 -  ctx - context for residual evaluation
1182 
1183    Calling sequence of `func`:
1184 $     PetscErrorCode func(TS ts, PetscReal t, Vec u, Mat A, Mat B, void *ctx);
1185 +  ts  - the `TS` context obtained from `TSCreate()`
1186 .  t - current timestep
1187 .  u - input vector
1188 .  Amat - (approximate) Jacobian matrix
1189 .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1190 -  ctx - [optional] user-defined context for matrix evaluation routine
1191 
1192    Level: advanced
1193 
1194    Note:
1195    `TSSetJacobian()` is normally used, but it calls this function internally because the user context is actually
1196    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1197    not.
1198 
1199    Developer Note:
1200    If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1201 
1202 .seealso: [](ch_ts), `DMTSSetContext()`, `TSSetFunction()`, `DMTSGetJacobian()`, `TSSetRHSJacobian()`
1203 @*/
1204 PetscErrorCode DMTSSetRHSJacobian(DM dm, TSRHSJacobian func, void *ctx)
1205 {
1206   DMTS tsdm;
1207 
1208   PetscFunctionBegin;
1209   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1210   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1211   if (func) tsdm->ops->rhsjacobian = func;
1212   if (ctx) {
1213     PetscContainer ctxcontainer;
1214     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
1215     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
1216     PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs jacobian ctx", (PetscObject)ctxcontainer));
1217     tsdm->rhsjacobianctxcontainer = ctxcontainer;
1218     PetscCall(PetscContainerDestroy(&ctxcontainer));
1219   }
1220   PetscFunctionReturn(PETSC_SUCCESS);
1221 }
1222 
1223 /*@C
1224    DMTSSetRHSJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function
1225 
1226    Not Collective
1227 
1228    Input Parameters:
1229 +  dm - `DM` to be used with `TS`
1230 -  f - Jacobian evaluation context destroy function
1231 
1232    Level: advanced
1233 
1234    Note:
1235    The user usually calls `TSSetRHSJacobianContextDestroy()` which calls this routine
1236 
1237 .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianContextDestroy()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
1238 @*/
1239 PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
1240 {
1241   DMTS tsdm;
1242 
1243   PetscFunctionBegin;
1244   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1245   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1246   if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsjacobianctxcontainer, f));
1247   PetscFunctionReturn(PETSC_SUCCESS);
1248 }
1249 
1250 PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM dm)
1251 {
1252   DMTS tsdm;
1253 
1254   PetscFunctionBegin;
1255   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1256   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1257   PetscCall(DMTSUnsetRHSJacobianContext_DMTS(tsdm));
1258   PetscFunctionReturn(PETSC_SUCCESS);
1259 }
1260 
1261 /*@C
1262    DMTSGetRHSJacobian - get `TS` Jacobian evaluation function
1263 
1264    Not Collective
1265 
1266    Input Parameter:
1267 .  dm - `DM` to be used with `TS`
1268 
1269    Output Parameters:
1270 +  func - Jacobian evaluation function, for calling sequence see `TSSetRHSJacobian()`
1271 -  ctx - context for residual evaluation
1272 
1273    Level: advanced
1274 
1275    Note:
1276    `TSGetJacobian()` is normally used, but it calls this function internally because the user context is actually
1277    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1278    not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1279 
1280 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
1281 @*/
1282 PetscErrorCode DMTSGetRHSJacobian(DM dm, TSRHSJacobian *func, void **ctx)
1283 {
1284   DMTS tsdm;
1285 
1286   PetscFunctionBegin;
1287   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1288   PetscCall(DMGetDMTS(dm, &tsdm));
1289   if (func) *func = tsdm->ops->rhsjacobian;
1290   if (ctx) {
1291     if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsjacobianctxcontainer, ctx));
1292     else *ctx = NULL;
1293   }
1294   PetscFunctionReturn(PETSC_SUCCESS);
1295 }
1296 
1297 /*@C
1298    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
1299 
1300    Not Collective
1301 
1302    Input Parameters:
1303 +  dm - `DM` to be used with `TS`
1304 .  view - viewer function
1305 -  load - loading function
1306 
1307    Level: advanced
1308 
1309 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1310 @*/
1311 PetscErrorCode DMTSSetIFunctionSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1312 {
1313   DMTS tsdm;
1314 
1315   PetscFunctionBegin;
1316   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1317   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1318   tsdm->ops->ifunctionview = view;
1319   tsdm->ops->ifunctionload = load;
1320   PetscFunctionReturn(PETSC_SUCCESS);
1321 }
1322 
1323 /*@C
1324    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
1325 
1326    Not Collective
1327 
1328    Input Parameters:
1329 +  dm - `DM` to be used with `TS`
1330 .  view - viewer function
1331 -  load - loading function
1332 
1333    Level: advanced
1334 
1335 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1336 @*/
1337 PetscErrorCode DMTSSetIJacobianSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1338 {
1339   DMTS tsdm;
1340 
1341   PetscFunctionBegin;
1342   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1343   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1344   tsdm->ops->ijacobianview = view;
1345   tsdm->ops->ijacobianload = load;
1346   PetscFunctionReturn(PETSC_SUCCESS);
1347 }
1348