xref: /petsc/src/ts/utils/dmts.c (revision dfbbaf821b4c49d07b4ce746493b0d955783fdf9)
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    Level: advanced
654 
655 .seealso: `TSSetI2JacobianContextDestroy()`, `DMTSSetI2Jacobian()`, `TSSetI2Jacobian()`
656 @*/
657 PetscErrorCode DMTSSetI2JacobianContextDestroy(DM dm,PetscErrorCode (*f)(void*))
658 {
659   DMTS           tsdm;
660 
661   PetscFunctionBegin;
662   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
663   PetscCall(DMGetDMTSWrite(dm,&tsdm));
664   if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2jacobianctxcontainer,f));
665   PetscFunctionReturn(0);
666 }
667 
668 PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM dm)
669 {
670   DMTS           tsdm;
671 
672   PetscFunctionBegin;
673   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
674   PetscCall(DMGetDMTSWrite(dm,&tsdm));
675   PetscCall(DMTSUnsetI2JacobianContext_DMTS(tsdm));
676   PetscFunctionReturn(0);
677 }
678 
679 /*@C
680    DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems
681 
682    Not Collective
683 
684    Input Parameter:
685 .  dm - DM to be used with TS
686 
687    Output Parameters:
688 +  jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
689 -  ctx - context for Jacobian evaluation
690 
691    Level: advanced
692 
693    Note:
694    TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
695    associated with the DM.
696 
697 .seealso: `DMTSSetI2Jacobian()`, `TSGetI2Jacobian()`
698 @*/
699 PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx)
700 {
701   DMTS           tsdm;
702 
703   PetscFunctionBegin;
704   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
705   PetscCall(DMGetDMTS(dm,&tsdm));
706   if (jac) *jac = tsdm->ops->i2jacobian;
707   if (ctx) {
708     if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2jacobianctxcontainer,ctx));
709     else *ctx = NULL;
710   }
711   PetscFunctionReturn(0);
712 }
713 
714 /*@C
715    DMTSSetRHSFunction - set TS explicit residual evaluation function
716 
717    Not Collective
718 
719    Input Parameters:
720 +  dm - DM to be used with TS
721 .  func - RHS function evaluation routine
722 -  ctx - context for residual evaluation
723 
724     Calling sequence of func:
725 $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec F,void *ctx);
726 
727 +   ts - timestep context
728 .   t - current timestep
729 .   u - input vector
730 .   F - function vector
731 -   ctx - [optional] user-defined function context
732 
733    Level: advanced
734 
735    Note:
736    TSSetRHSFunction() is normally used, but it calls this function internally because the user context is actually
737    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
738    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
739 
740 .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetJacobian()`
741 @*/
742 PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
743 {
744   DMTS           tsdm;
745 
746   PetscFunctionBegin;
747   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
748   PetscCall(DMGetDMTSWrite(dm,&tsdm));
749   if (func) tsdm->ops->rhsfunction = func;
750   if (ctx) {
751     PetscContainer ctxcontainer;
752     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm),&ctxcontainer));
753     PetscCall(PetscContainerSetPointer(ctxcontainer,ctx));
754     PetscCall(PetscObjectCompose((PetscObject)tsdm,"rhs function ctx",(PetscObject)ctxcontainer));
755     tsdm->rhsfunctionctxcontainer = ctxcontainer;
756     PetscCall(PetscContainerDestroy(&ctxcontainer));
757   }
758   PetscFunctionReturn(0);
759 }
760 
761 /*@C
762    DMTSSetRHSFunctionContextDestroy - set `TS` explicit residual evaluation context destroy function
763 
764    Not Collective
765 
766    Input Parameters:
767 +  dm - `DM` to be used with `TS`
768 -  f - explicit evaluation context destroy function
769 
770    Level: advanced
771 
772    Note:
773    `TSSetRHSFunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
774    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
775    not.
776 
777    Developer Note:
778    If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
779 
780 .seealso: `TSSetRHSFunctionContextDestroy()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
781 @*/
782 PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM dm,PetscErrorCode (*f)(void*))
783 {
784   DMTS           tsdm;
785 
786   PetscFunctionBegin;
787   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
788   PetscCall(DMGetDMTSWrite(dm,&tsdm));
789   if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsfunctionctxcontainer,f));
790   PetscFunctionReturn(0);
791 }
792 
793 PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM dm)
794 {
795   DMTS           tsdm;
796 
797   PetscFunctionBegin;
798   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
799   PetscCall(DMGetDMTSWrite(dm,&tsdm));
800   PetscCall(DMTSUnsetRHSFunctionContext_DMTS(tsdm));
801   tsdm->rhsfunctionctxcontainer = NULL;
802   PetscFunctionReturn(0);
803 }
804 
805 /*@C
806    DMTSSetTransientVariable - sets function to transform from state to transient variables
807 
808    Logically Collective
809 
810    Input Parameters:
811 +  dm - DM to be used with TS
812 .  tvar - a function that transforms to transient variables
813 -  ctx - a context for tvar
814 
815     Calling sequence of tvar:
816 $     PetscErrorCode tvar(TS ts,Vec p,Vec c,void *ctx);
817 
818 +   ts - timestep context
819 .   p - input vector (primitive form)
820 .   c - output vector, transient variables (conservative form)
821 -   ctx - [optional] user-defined function context
822 
823    Level: advanced
824 
825    Notes:
826    This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., TSBDF)
827    can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
828    well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
829    C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
830    evaluated via the chain rule, as in
831 
832      dF/dP + shift * dF/dCdot dC/dP.
833 
834 .seealso: `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()`
835 @*/
836 PetscErrorCode DMTSSetTransientVariable(DM dm,TSTransientVariable tvar,void *ctx)
837 {
838   DMTS           dmts;
839 
840   PetscFunctionBegin;
841   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
842   PetscCall(DMGetDMTSWrite(dm,&dmts));
843   dmts->ops->transientvar = tvar;
844   dmts->transientvarctx = ctx;
845   PetscFunctionReturn(0);
846 }
847 
848 /*@C
849    DMTSGetTransientVariable - gets function to transform from state to transient variables
850 
851    Logically Collective
852 
853    Input Parameter:
854 .  dm - DM to be used with TS
855 
856    Output Parameters:
857 +  tvar - a function that transforms to transient variables
858 -  ctx - a context for tvar
859 
860    Level: advanced
861 
862 .seealso: `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()`
863 @*/
864 PetscErrorCode DMTSGetTransientVariable(DM dm,TSTransientVariable *tvar,void *ctx)
865 {
866   DMTS           dmts;
867 
868   PetscFunctionBegin;
869   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
870   PetscCall(DMGetDMTS(dm,&dmts));
871   if (tvar) *tvar = dmts->ops->transientvar;
872   if (ctx)  *(void**)ctx = dmts->transientvarctx;
873   PetscFunctionReturn(0);
874 }
875 
876 /*@C
877    DMTSGetSolutionFunction - gets the TS solution evaluation function
878 
879    Not Collective
880 
881    Input Parameter:
882 .  dm - DM to be used with TS
883 
884    Output Parameters:
885 +  func - solution function evaluation function, see TSSetSolution() for calling sequence
886 -  ctx - context for solution evaluation
887 
888    Level: advanced
889 
890 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSSetSolutionFunction()`
891 @*/
892 PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
893 {
894   DMTS           tsdm;
895 
896   PetscFunctionBegin;
897   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
898   PetscCall(DMGetDMTS(dm,&tsdm));
899   if (func) *func = tsdm->ops->solution;
900   if (ctx)  *ctx  = tsdm->solutionctx;
901   PetscFunctionReturn(0);
902 }
903 
904 /*@C
905    DMTSSetSolutionFunction - set TS solution evaluation function
906 
907    Not Collective
908 
909    Input Parameters:
910 +  dm - DM to be used with TS
911 .  func - solution function evaluation routine
912 -  ctx - context for solution evaluation
913 
914     Calling sequence of f:
915 $     PetscErrorCode f(TS ts,PetscReal t,Vec u,void *ctx);
916 
917 +   ts - timestep context
918 .   t - current timestep
919 .   u - output vector
920 -   ctx - [optional] user-defined function context
921 
922    Level: advanced
923 
924    Note:
925    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
926    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
927    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
928 
929 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSGetSolutionFunction()`
930 @*/
931 PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
932 {
933   DMTS           tsdm;
934 
935   PetscFunctionBegin;
936   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
937   PetscCall(DMGetDMTSWrite(dm,&tsdm));
938   if (func) tsdm->ops->solution = func;
939   if (ctx)  tsdm->solutionctx   = ctx;
940   PetscFunctionReturn(0);
941 }
942 
943 /*@C
944    DMTSSetForcingFunction - set TS forcing function evaluation function
945 
946    Not Collective
947 
948    Input Parameters:
949 +  dm - DM to be used with TS
950 .  f - forcing function evaluation routine
951 -  ctx - context for solution evaluation
952 
953     Calling sequence of func:
954 $     PetscErrorCode func (TS ts,PetscReal t,Vec f,void *ctx);
955 
956 +   ts - timestep context
957 .   t - current timestep
958 .   f - output vector
959 -   ctx - [optional] user-defined function context
960 
961    Level: advanced
962 
963    Note:
964    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
965    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
966    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
967 
968 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
969 @*/
970 PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx)
971 {
972   DMTS           tsdm;
973 
974   PetscFunctionBegin;
975   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
976   PetscCall(DMGetDMTSWrite(dm,&tsdm));
977   if (f)    tsdm->ops->forcing = f;
978   if (ctx)  tsdm->forcingctx   = ctx;
979   PetscFunctionReturn(0);
980 }
981 
982 /*@C
983    DMTSGetForcingFunction - get TS forcing function evaluation function
984 
985    Not Collective
986 
987    Input Parameter:
988 .   dm - DM to be used with TS
989 
990    Output Parameters:
991 +  f - forcing function evaluation function; see TSForcingFunction for details
992 -  ctx - context for solution evaluation
993 
994    Level: advanced
995 
996    Note:
997    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
998    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
999    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
1000 
1001 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
1002 @*/
1003 PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
1004 {
1005   DMTS           tsdm;
1006 
1007   PetscFunctionBegin;
1008   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1009   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1010   if (f)   *f   = tsdm->ops->forcing;
1011   if (ctx) *ctx = tsdm->forcingctx;
1012   PetscFunctionReturn(0);
1013 }
1014 
1015 /*@C
1016    DMTSGetRHSFunction - get TS explicit residual evaluation function
1017 
1018    Not Collective
1019 
1020    Input Parameter:
1021 .  dm - DM to be used with TS
1022 
1023    Output Parameters:
1024 +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
1025 -  ctx - context for residual evaluation
1026 
1027    Level: advanced
1028 
1029    Note:
1030    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
1031    associated with the DM.
1032 
1033 .seealso: `DMTSSetContext()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
1034 @*/
1035 PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
1036 {
1037   DMTS           tsdm;
1038 
1039   PetscFunctionBegin;
1040   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1041   PetscCall(DMGetDMTS(dm,&tsdm));
1042   if (func) *func = tsdm->ops->rhsfunction;
1043   if (ctx) {
1044     if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsfunctionctxcontainer,ctx));
1045     else *ctx = NULL;
1046   }
1047   PetscFunctionReturn(0);
1048 }
1049 
1050 /*@C
1051    DMTSSetIJacobian - set TS Jacobian evaluation function
1052 
1053    Not Collective
1054 
1055    Input Parameters:
1056 +  dm - DM to be used with TS
1057 .  func - Jacobian evaluation routine
1058 -  ctx - context for residual evaluation
1059 
1060    Calling sequence of f:
1061 $    PetscErrorCode f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx);
1062 
1063 +  t    - time at step/stage being solved
1064 .  U    - state vector
1065 .  U_t  - time derivative of state vector
1066 .  a    - shift
1067 .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
1068 .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
1069 -  ctx  - [optional] user-defined context for matrix evaluation routine
1070 
1071    Level: advanced
1072 
1073    Note:
1074    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
1075    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
1076    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
1077 
1078 .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSGetJacobian()`, `TSSetIJacobian()`, `TSSetIFunction()`
1079 @*/
1080 PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
1081 {
1082   DMTS           tsdm;
1083 
1084   PetscFunctionBegin;
1085   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1086   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1087   if (func) tsdm->ops->ijacobian = func;
1088   if (ctx) {
1089     PetscContainer ctxcontainer;
1090     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm),&ctxcontainer));
1091     PetscCall(PetscContainerSetPointer(ctxcontainer,ctx));
1092     PetscCall(PetscObjectCompose((PetscObject)tsdm,"ijacobian ctx",(PetscObject)ctxcontainer));
1093     tsdm->ijacobianctxcontainer = ctxcontainer;
1094     PetscCall(PetscContainerDestroy(&ctxcontainer));
1095   }
1096   PetscFunctionReturn(0);
1097 }
1098 
1099 /*@C
1100    DMTSSetIJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function
1101 
1102    Not Collective
1103 
1104    Input Parameters:
1105 +  dm - `DM` to be used with `TS`
1106 -  f - Jacobian evaluation context destroy function
1107 
1108    Level: advanced
1109 
1110    Note:
1111    `TSSetIJacobianContextDestroy()` is normally used, but it calls this function internally because the user context is actually
1112    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1113    not.
1114 
1115    Developer Note:
1116    If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1117 
1118 .seealso: `TSSetIJacobianContextDestroy()`, `TSSetI2JacobianContextDestroy()`, `DMTSSetIJacobian()`, `TSSetIJacobian()`
1119 @*/
1120 PetscErrorCode DMTSSetIJacobianContextDestroy(DM dm,PetscErrorCode (*f)(void*))
1121 {
1122   DMTS           tsdm;
1123 
1124   PetscFunctionBegin;
1125   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1126   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1127   if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ijacobianctxcontainer,f));
1128   PetscFunctionReturn(0);
1129 }
1130 
1131 PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM dm)
1132 {
1133   DMTS           tsdm;
1134 
1135   PetscFunctionBegin;
1136   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1137   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1138   PetscCall(DMTSUnsetIJacobianContext_DMTS(tsdm));
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 /*@C
1143    DMTSGetIJacobian - get TS Jacobian evaluation function
1144 
1145    Not Collective
1146 
1147    Input Parameter:
1148 .  dm - DM to be used with TS
1149 
1150    Output Parameters:
1151 +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
1152 -  ctx - context for residual evaluation
1153 
1154    Level: advanced
1155 
1156    Note:
1157    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
1158    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
1159    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
1160 
1161 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1162 @*/
1163 PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
1164 {
1165   DMTS           tsdm;
1166 
1167   PetscFunctionBegin;
1168   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1169   PetscCall(DMGetDMTS(dm,&tsdm));
1170   if (func) *func = tsdm->ops->ijacobian;
1171   if (ctx) {
1172     if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ijacobianctxcontainer,ctx));
1173     else *ctx = NULL;
1174   }
1175   PetscFunctionReturn(0);
1176 }
1177 
1178 /*@C
1179    DMTSSetRHSJacobian - set TS Jacobian evaluation function
1180 
1181    Not Collective
1182 
1183    Input Parameters:
1184 +  dm - DM to be used with TS
1185 .  func - Jacobian evaluation routine
1186 -  ctx - context for residual evaluation
1187 
1188    Calling sequence of func:
1189 $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);
1190 
1191 +  t - current timestep
1192 .  u - input vector
1193 .  Amat - (approximate) Jacobian matrix
1194 .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1195 -  ctx - [optional] user-defined context for matrix evaluation routine
1196 
1197    Level: advanced
1198 
1199    Note:
1200    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
1201    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
1202    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
1203 
1204 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSGetJacobian()`, `TSSetRHSJacobian()`
1205 @*/
1206 PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
1207 {
1208   DMTS           tsdm;
1209 
1210   PetscFunctionBegin;
1211   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1212   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1213   if (func) tsdm->ops->rhsjacobian = func;
1214   if (ctx) {
1215     PetscContainer ctxcontainer;
1216     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm),&ctxcontainer));
1217     PetscCall(PetscContainerSetPointer(ctxcontainer,ctx));
1218     PetscCall(PetscObjectCompose((PetscObject)tsdm,"rhs jacobian ctx",(PetscObject)ctxcontainer));
1219     tsdm->rhsjacobianctxcontainer = ctxcontainer;
1220     PetscCall(PetscContainerDestroy(&ctxcontainer));
1221   }
1222   PetscFunctionReturn(0);
1223 }
1224 
1225 /*@C
1226    DMTSSetRHSJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function
1227 
1228    Not Collective
1229 
1230    Input Parameters:
1231 +  dm - `DM` to be used with `TS`
1232 -  f - Jacobian evaluation context destroy function
1233 
1234    Level: advanced
1235 
1236    Note:
1237    The user usually calls `TSSetRHSJacobianContextDestroy()` which calls this routine
1238 
1239    Level: advanced
1240 
1241 .seealso: `TSSetRHSJacobianContextDestroy()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
1242 @*/
1243 PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM dm,PetscErrorCode (*f)(void*))
1244 {
1245   DMTS           tsdm;
1246 
1247   PetscFunctionBegin;
1248   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1249   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1250   if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsjacobianctxcontainer,f));
1251   PetscFunctionReturn(0);
1252 }
1253 
1254 PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM dm)
1255 {
1256   DMTS           tsdm;
1257 
1258   PetscFunctionBegin;
1259   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1260   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1261   PetscCall(DMTSUnsetRHSJacobianContext_DMTS(tsdm));
1262   PetscFunctionReturn(0);
1263 }
1264 
1265 /*@C
1266    DMTSGetRHSJacobian - get TS Jacobian evaluation function
1267 
1268    Not Collective
1269 
1270    Input Parameter:
1271 .  dm - DM to be used with TS
1272 
1273    Output Parameters:
1274 +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
1275 -  ctx - context for residual evaluation
1276 
1277    Level: advanced
1278 
1279    Note:
1280    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
1281    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
1282    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
1283 
1284 .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
1285 @*/
1286 PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
1287 {
1288   DMTS           tsdm;
1289 
1290   PetscFunctionBegin;
1291   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1292   PetscCall(DMGetDMTS(dm,&tsdm));
1293   if (func) *func = tsdm->ops->rhsjacobian;
1294   if (ctx) {
1295     if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsjacobianctxcontainer,ctx));
1296     else *ctx = NULL;
1297   }
1298   PetscFunctionReturn(0);
1299 }
1300 
1301 /*@C
1302    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
1303 
1304    Not Collective
1305 
1306    Input Parameters:
1307 +  dm - DM to be used with TS
1308 .  view - viewer function
1309 -  load - loading function
1310 
1311    Level: advanced
1312 
1313 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1314 @*/
1315 PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
1316 {
1317   DMTS           tsdm;
1318 
1319   PetscFunctionBegin;
1320   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1321   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1322   tsdm->ops->ifunctionview = view;
1323   tsdm->ops->ifunctionload = load;
1324   PetscFunctionReturn(0);
1325 }
1326 
1327 /*@C
1328    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
1329 
1330    Not Collective
1331 
1332    Input Parameters:
1333 +  dm - DM to be used with TS
1334 .  view - viewer function
1335 -  load - loading function
1336 
1337    Level: advanced
1338 
1339 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1340 @*/
1341 PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
1342 {
1343   DMTS           tsdm;
1344 
1345   PetscFunctionBegin;
1346   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1347   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1348   tsdm->ops->ijacobianview = view;
1349   tsdm->ops->ijacobianload = load;
1350   PetscFunctionReturn(0);
1351 }
1352