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