xref: /petsc/src/ts/utils/dmts.c (revision 14d0ab18b70ee075d422f67a0c1395817de67fab)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
2af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
324989b8cSPeter Brune 
4d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSUnsetRHSFunctionContext_DMTS(DMTS tsdm)
5d71ae5a4SJacob Faibussowitsch {
6800f99ffSJeremy L Thompson   PetscFunctionBegin;
7800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs function ctx", NULL));
8800f99ffSJeremy L Thompson   tsdm->rhsfunctionctxcontainer = NULL;
93ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10800f99ffSJeremy L Thompson }
11800f99ffSJeremy L Thompson 
12d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSUnsetRHSJacobianContext_DMTS(DMTS tsdm)
13d71ae5a4SJacob Faibussowitsch {
14800f99ffSJeremy L Thompson   PetscFunctionBegin;
15800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs jacobian ctx", NULL));
16800f99ffSJeremy L Thompson   tsdm->rhsjacobianctxcontainer = NULL;
173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18800f99ffSJeremy L Thompson }
19800f99ffSJeremy L Thompson 
20d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSUnsetIFunctionContext_DMTS(DMTS tsdm)
21d71ae5a4SJacob Faibussowitsch {
22800f99ffSJeremy L Thompson   PetscFunctionBegin;
23800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm, "ifunction ctx", NULL));
24800f99ffSJeremy L Thompson   tsdm->ifunctionctxcontainer = NULL;
253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26800f99ffSJeremy L Thompson }
27800f99ffSJeremy L Thompson 
28d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSUnsetIJacobianContext_DMTS(DMTS tsdm)
29d71ae5a4SJacob Faibussowitsch {
30800f99ffSJeremy L Thompson   PetscFunctionBegin;
31800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm, "ijacobian ctx", NULL));
32800f99ffSJeremy L Thompson   tsdm->ijacobianctxcontainer = NULL;
333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34800f99ffSJeremy L Thompson }
35800f99ffSJeremy L Thompson 
36d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSUnsetI2FunctionContext_DMTS(DMTS tsdm)
37d71ae5a4SJacob Faibussowitsch {
38800f99ffSJeremy L Thompson   PetscFunctionBegin;
39800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2function ctx", NULL));
40800f99ffSJeremy L Thompson   tsdm->i2functionctxcontainer = NULL;
413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42800f99ffSJeremy L Thompson }
43800f99ffSJeremy L Thompson 
44d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSUnsetI2JacobianContext_DMTS(DMTS tsdm)
45d71ae5a4SJacob Faibussowitsch {
46800f99ffSJeremy L Thompson   PetscFunctionBegin;
47800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2jacobian ctx", NULL));
48800f99ffSJeremy L Thompson   tsdm->i2jacobianctxcontainer = NULL;
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50800f99ffSJeremy L Thompson }
51800f99ffSJeremy L Thompson 
52d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSDestroy(DMTS *kdm)
53d71ae5a4SJacob Faibussowitsch {
54d74926cbSBarry Smith   PetscFunctionBegin;
553ba16761SJacob Faibussowitsch   if (!*kdm) PetscFunctionReturn(PETSC_SUCCESS);
56d74926cbSBarry Smith   PetscValidHeaderSpecific((*kdm), DMTS_CLASSID, 1);
579371c9d4SSatish Balay   if (--((PetscObject)(*kdm))->refct > 0) {
589371c9d4SSatish Balay     *kdm = NULL;
593ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
609371c9d4SSatish Balay   }
61800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetRHSFunctionContext_DMTS(*kdm));
62800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetRHSJacobianContext_DMTS(*kdm));
63800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIFunctionContext_DMTS(*kdm));
64800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIJacobianContext_DMTS(*kdm));
65800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2FunctionContext_DMTS(*kdm));
66800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2JacobianContext_DMTS(*kdm));
67dbbe0bcdSBarry Smith   PetscTryTypeMethod(*kdm, destroy);
689566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(kdm));
693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70d74926cbSBarry Smith }
71d74926cbSBarry Smith 
72d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSLoad(DMTS kdm, PetscViewer viewer)
73d71ae5a4SJacob Faibussowitsch {
742d53ad75SBarry Smith   PetscFunctionBegin;
759566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunction, 1, NULL, PETSC_FUNCTION));
769566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunctionview, 1, NULL, PETSC_FUNCTION));
779566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunctionload, 1, NULL, PETSC_FUNCTION));
78ad6bc421SBarry Smith   if (kdm->ops->ifunctionload) {
79800f99ffSJeremy L Thompson     void *ctx;
80800f99ffSJeremy L Thompson 
81800f99ffSJeremy L Thompson     PetscCall(PetscContainerGetPointer(kdm->ifunctionctxcontainer, &ctx));
82800f99ffSJeremy L Thompson     PetscCall((*kdm->ops->ifunctionload)(&ctx, viewer));
83ad6bc421SBarry Smith   }
849566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobian, 1, NULL, PETSC_FUNCTION));
859566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobianview, 1, NULL, PETSC_FUNCTION));
869566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobianload, 1, NULL, PETSC_FUNCTION));
87ad6bc421SBarry Smith   if (kdm->ops->ijacobianload) {
88800f99ffSJeremy L Thompson     void *ctx;
89800f99ffSJeremy L Thompson 
90800f99ffSJeremy L Thompson     PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer, &ctx));
91800f99ffSJeremy L Thompson     PetscCall((*kdm->ops->ijacobianload)(&ctx, viewer));
92ad6bc421SBarry Smith   }
933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
942d53ad75SBarry Smith }
952d53ad75SBarry Smith 
96d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSView(DMTS kdm, PetscViewer viewer)
97d71ae5a4SJacob Faibussowitsch {
982d53ad75SBarry Smith   PetscBool isascii, isbinary;
992d53ad75SBarry Smith 
1002d53ad75SBarry Smith   PetscFunctionBegin;
1019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1032d53ad75SBarry Smith   if (isascii) {
104c7a10e08SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
1052d53ad75SBarry Smith     const char *fname;
1062d53ad75SBarry Smith 
1079566063dSJacob Faibussowitsch     PetscCall(PetscFPTFind(kdm->ops->ifunction, &fname));
10848a46eb9SPierre Jolivet     if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, "  IFunction used by TS: %s\n", fname));
1099566063dSJacob Faibussowitsch     PetscCall(PetscFPTFind(kdm->ops->ijacobian, &fname));
11048a46eb9SPierre Jolivet     if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, "  IJacobian function used by TS: %s\n", fname));
111c7a10e08SBarry Smith #endif
1122d53ad75SBarry Smith   } else if (isbinary) {
1133964eb88SJed Brown     struct {
1143964eb88SJed Brown       TSIFunction ifunction;
1159200755eSBarry Smith     } funcstruct;
1169200755eSBarry Smith     struct {
1173964eb88SJed Brown       PetscErrorCode (*ifunctionview)(void *, PetscViewer);
1189200755eSBarry Smith     } funcviewstruct;
1199200755eSBarry Smith     struct {
1203964eb88SJed Brown       PetscErrorCode (*ifunctionload)(void **, PetscViewer);
1219200755eSBarry Smith     } funcloadstruct;
1223964eb88SJed Brown     struct {
1233964eb88SJed Brown       TSIJacobian ijacobian;
1249200755eSBarry Smith     } jacstruct;
1259200755eSBarry Smith     struct {
1263964eb88SJed Brown       PetscErrorCode (*ijacobianview)(void *, PetscViewer);
1279200755eSBarry Smith     } jacviewstruct;
1289200755eSBarry Smith     struct {
1293964eb88SJed Brown       PetscErrorCode (*ijacobianload)(void **, PetscViewer);
1309200755eSBarry Smith     } jacloadstruct;
1313964eb88SJed Brown 
1329200755eSBarry Smith     funcstruct.ifunction         = kdm->ops->ifunction;
1339200755eSBarry Smith     funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
1349200755eSBarry Smith     funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
1359566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer, &funcstruct, 1, PETSC_FUNCTION));
1369566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer, &funcviewstruct, 1, PETSC_FUNCTION));
1379566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer, &funcloadstruct, 1, PETSC_FUNCTION));
138800f99ffSJeremy L Thompson     if (kdm->ops->ifunctionview) {
139800f99ffSJeremy L Thompson       void *ctx;
140800f99ffSJeremy L Thompson 
141800f99ffSJeremy L Thompson       PetscCall(PetscContainerGetPointer(kdm->ifunctionctxcontainer, &ctx));
142800f99ffSJeremy L Thompson       PetscCall((*kdm->ops->ifunctionview)(ctx, viewer));
143800f99ffSJeremy L Thompson     }
1449200755eSBarry Smith     jacstruct.ijacobian         = kdm->ops->ijacobian;
1459200755eSBarry Smith     jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
1469200755eSBarry Smith     jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
1479566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer, &jacstruct, 1, PETSC_FUNCTION));
1489566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer, &jacviewstruct, 1, PETSC_FUNCTION));
1499566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer, &jacloadstruct, 1, PETSC_FUNCTION));
150800f99ffSJeremy L Thompson     if (kdm->ops->ijacobianview) {
151800f99ffSJeremy L Thompson       void *ctx;
152800f99ffSJeremy L Thompson 
153800f99ffSJeremy L Thompson       PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer, &ctx));
154800f99ffSJeremy L Thompson       PetscCall((*kdm->ops->ijacobianview)(ctx, viewer));
155800f99ffSJeremy L Thompson     }
1562d53ad75SBarry Smith   }
1573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1582d53ad75SBarry Smith }
1592d53ad75SBarry Smith 
160d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSCreate(MPI_Comm comm, DMTS *kdm)
161d71ae5a4SJacob Faibussowitsch {
162d74926cbSBarry Smith   PetscFunctionBegin;
1639566063dSJacob Faibussowitsch   PetscCall(TSInitializePackage());
1649566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView));
1653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
166d74926cbSBarry Smith }
16724989b8cSPeter Brune 
1682a34c10cSBarry Smith /* Attaches the DMTS to the coarse level.
16924989b8cSPeter Brune  * Under what conditions should we copy versus duplicate?
17024989b8cSPeter Brune  */
171d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_DMTS(DM dm, DM dmc, void *ctx)
172d71ae5a4SJacob Faibussowitsch {
17324989b8cSPeter Brune   PetscFunctionBegin;
1749566063dSJacob Faibussowitsch   PetscCall(DMCopyDMTS(dm, dmc));
1753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17624989b8cSPeter Brune }
17724989b8cSPeter Brune 
17824989b8cSPeter Brune /* This could restrict auxiliary information to the coarse level.
17924989b8cSPeter Brune  */
180d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_DMTS(DM dm, Mat Restrict, Vec rscale, Mat Inject, DM dmc, void *ctx)
181d71ae5a4SJacob Faibussowitsch {
18224989b8cSPeter Brune   PetscFunctionBegin;
1833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18424989b8cSPeter Brune }
18524989b8cSPeter Brune 
186d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_DMTS(DM dm, DM subdm, void *ctx)
187d71ae5a4SJacob Faibussowitsch {
188258e1594SPeter Brune   PetscFunctionBegin;
1899566063dSJacob Faibussowitsch   PetscCall(DMCopyDMTS(dm, subdm));
1903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
191258e1594SPeter Brune }
192258e1594SPeter Brune 
193258e1594SPeter Brune /* This could restrict auxiliary information to the coarse level.
194258e1594SPeter Brune  */
195d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
196d71ae5a4SJacob Faibussowitsch {
197258e1594SPeter Brune   PetscFunctionBegin;
1983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
199258e1594SPeter Brune }
200258e1594SPeter Brune 
201d74926cbSBarry Smith /*@C
202bcf0153eSBarry Smith   DMTSCopy - copies the information in a `DMTS` to another `DMTS`
203d74926cbSBarry Smith 
204d74926cbSBarry Smith   Not Collective
205d74926cbSBarry Smith 
2064165533cSJose E. Roman   Input Parameters:
207bcf0153eSBarry Smith + kdm  - Original `DMTS`
208bcf0153eSBarry Smith - nkdm - `DMTS` to receive the data, should have been created with `DMTSCreate()`
209d74926cbSBarry Smith 
210d74926cbSBarry Smith   Level: developer
211d74926cbSBarry Smith 
2121cc06b55SBarry Smith .seealso: [](ch_ts), `DMTSCreate()`, `DMTSDestroy()`
213d74926cbSBarry Smith @*/
214d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCopy(DMTS kdm, DMTS nkdm)
215d71ae5a4SJacob Faibussowitsch {
21624989b8cSPeter Brune   PetscFunctionBegin;
217d74926cbSBarry Smith   PetscValidHeaderSpecific(kdm, DMTS_CLASSID, 1);
218d74926cbSBarry Smith   PetscValidHeaderSpecific(nkdm, DMTS_CLASSID, 2);
219d74926cbSBarry Smith   nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
220d74926cbSBarry Smith   nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
221d74926cbSBarry Smith   nkdm->ops->ifunction   = kdm->ops->ifunction;
222d74926cbSBarry Smith   nkdm->ops->ijacobian   = kdm->ops->ijacobian;
223efe9872eSLisandro Dalcin   nkdm->ops->i2function  = kdm->ops->i2function;
224efe9872eSLisandro Dalcin   nkdm->ops->i2jacobian  = kdm->ops->i2jacobian;
225d74926cbSBarry Smith   nkdm->ops->solution    = kdm->ops->solution;
226d74926cbSBarry Smith   nkdm->ops->destroy     = kdm->ops->destroy;
227d74926cbSBarry Smith   nkdm->ops->duplicate   = kdm->ops->duplicate;
228d74926cbSBarry Smith 
229d74926cbSBarry Smith   nkdm->solutionctx             = kdm->solutionctx;
230800f99ffSJeremy L Thompson   nkdm->rhsfunctionctxcontainer = kdm->rhsfunctionctxcontainer;
231800f99ffSJeremy L Thompson   nkdm->rhsjacobianctxcontainer = kdm->rhsjacobianctxcontainer;
232800f99ffSJeremy L Thompson   nkdm->ifunctionctxcontainer   = kdm->ifunctionctxcontainer;
233800f99ffSJeremy L Thompson   nkdm->ijacobianctxcontainer   = kdm->ijacobianctxcontainer;
234800f99ffSJeremy L Thompson   nkdm->i2functionctxcontainer  = kdm->i2functionctxcontainer;
235800f99ffSJeremy L Thompson   nkdm->i2jacobianctxcontainer  = kdm->i2jacobianctxcontainer;
236800f99ffSJeremy L Thompson   if (nkdm->rhsfunctionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "rhs function ctx", (PetscObject)nkdm->rhsfunctionctxcontainer));
237800f99ffSJeremy L Thompson   if (nkdm->rhsjacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "rhs jacobian ctx", (PetscObject)nkdm->rhsjacobianctxcontainer));
238800f99ffSJeremy L Thompson   if (nkdm->ifunctionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "ifunction ctx", (PetscObject)nkdm->ifunctionctxcontainer));
239800f99ffSJeremy L Thompson   if (nkdm->ijacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "ijacobian ctx", (PetscObject)nkdm->ijacobianctxcontainer));
240800f99ffSJeremy L Thompson   if (nkdm->i2functionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "i2function ctx", (PetscObject)nkdm->i2functionctxcontainer));
241800f99ffSJeremy L Thompson   if (nkdm->i2jacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "i2jacobian ctx", (PetscObject)nkdm->i2jacobianctxcontainer));
242d74926cbSBarry Smith 
243d74926cbSBarry Smith   nkdm->data = kdm->data;
244d74926cbSBarry Smith 
245d74926cbSBarry Smith   /*
246d74926cbSBarry Smith   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
247d74926cbSBarry Smith   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
248d74926cbSBarry Smith   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
249d74926cbSBarry Smith   */
250d74926cbSBarry Smith 
251d74926cbSBarry Smith   /* implementation specific copy hooks */
252dbbe0bcdSBarry Smith   PetscTryTypeMethod(kdm, duplicate, nkdm);
2533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25424989b8cSPeter Brune }
25524989b8cSPeter Brune 
25624989b8cSPeter Brune /*@C
257bcf0153eSBarry Smith   DMGetDMTS - get read-only private `DMTS` context from a `DM`
25824989b8cSPeter Brune 
25924989b8cSPeter Brune   Not Collective
26024989b8cSPeter Brune 
2614165533cSJose E. Roman   Input Parameter:
262bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
26324989b8cSPeter Brune 
2644165533cSJose E. Roman   Output Parameter:
265bcf0153eSBarry Smith . tsdm - private `DMTS` context
26624989b8cSPeter Brune 
26724989b8cSPeter Brune   Level: developer
26824989b8cSPeter Brune 
26924989b8cSPeter Brune   Notes:
270bcf0153eSBarry Smith   Use `DMGetDMTSWrite()` if write access is needed. The `DMTSSetXXX()` API should be used wherever possible.
27124989b8cSPeter Brune 
2721cc06b55SBarry Smith .seealso: [](ch_ts), `DMTS`, `DMGetDMTSWrite()`
27324989b8cSPeter Brune @*/
274d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetDMTS(DM dm, DMTS *tsdm)
275d71ae5a4SJacob Faibussowitsch {
27624989b8cSPeter Brune   PetscFunctionBegin;
27724989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2782a34c10cSBarry Smith   *tsdm = (DMTS)dm->dmts;
279d74926cbSBarry Smith   if (!*tsdm) {
2809566063dSJacob Faibussowitsch     PetscCall(PetscInfo(dm, "Creating new DMTS\n"));
2819566063dSJacob Faibussowitsch     PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm), tsdm));
2822a34c10cSBarry Smith     dm->dmts            = (PetscObject)*tsdm;
2835c87d4f4SJunchao Zhang     (*tsdm)->originaldm = dm;
2849566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_DMTS, DMRestrictHook_DMTS, NULL));
2859566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_DMTS, DMSubDomainRestrictHook_DMTS, NULL));
28624989b8cSPeter Brune   }
2873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28824989b8cSPeter Brune }
28924989b8cSPeter Brune 
29024989b8cSPeter Brune /*@C
291bcf0153eSBarry Smith   DMGetDMTSWrite - get write access to private `DMTS` context from a `DM`
29224989b8cSPeter Brune 
29324989b8cSPeter Brune   Not Collective
29424989b8cSPeter Brune 
2954165533cSJose E. Roman   Input Parameter:
296bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
29724989b8cSPeter Brune 
2984165533cSJose E. Roman   Output Parameter:
299bcf0153eSBarry Smith . tsdm - private `DMTS` context
30024989b8cSPeter Brune 
30124989b8cSPeter Brune   Level: developer
30224989b8cSPeter Brune 
3031cc06b55SBarry Smith .seealso: [](ch_ts), `DMTS`, `DMGetDMTS()`
30424989b8cSPeter Brune @*/
305d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetDMTSWrite(DM dm, DMTS *tsdm)
306d71ae5a4SJacob Faibussowitsch {
307942e3340SBarry Smith   DMTS sdm;
30824989b8cSPeter Brune 
30924989b8cSPeter Brune   PetscFunctionBegin;
31024989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3119566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &sdm));
3123c633725SBarry Smith   PetscCheck(sdm->originaldm, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DMTS has a NULL originaldm");
31324989b8cSPeter Brune   if (sdm->originaldm != dm) { /* Copy on write */
3142a34c10cSBarry Smith     DMTS oldsdm = sdm;
3159566063dSJacob Faibussowitsch     PetscCall(PetscInfo(dm, "Copying DMTS due to write\n"));
3169566063dSJacob Faibussowitsch     PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm), &sdm));
3179566063dSJacob Faibussowitsch     PetscCall(DMTSCopy(oldsdm, sdm));
3189566063dSJacob Faibussowitsch     PetscCall(DMTSDestroy((DMTS *)&dm->dmts));
3192a34c10cSBarry Smith     dm->dmts        = (PetscObject)sdm;
3205c87d4f4SJunchao Zhang     sdm->originaldm = dm;
32124989b8cSPeter Brune   }
32224989b8cSPeter Brune   *tsdm = sdm;
3233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32424989b8cSPeter Brune }
32524989b8cSPeter Brune 
32624989b8cSPeter Brune /*@C
327bcf0153eSBarry Smith   DMCopyDMTS - copies a `DM` context to a new `DM`
32824989b8cSPeter Brune 
32924989b8cSPeter Brune   Logically Collective
33024989b8cSPeter Brune 
3314165533cSJose E. Roman   Input Parameters:
332bcf0153eSBarry Smith + dmsrc  - `DM` to obtain context from
333bcf0153eSBarry Smith - dmdest - `DM` to add context to
33424989b8cSPeter Brune 
33524989b8cSPeter Brune   Level: developer
33624989b8cSPeter Brune 
33724989b8cSPeter Brune   Note:
33824989b8cSPeter Brune   The context is copied by reference. This function does not ensure that a context exists.
33924989b8cSPeter Brune 
3401cc06b55SBarry Smith .seealso: [](ch_ts), `DMTS`, `DMGetDMTS()`, `TSSetDM()`
34124989b8cSPeter Brune @*/
342d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCopyDMTS(DM dmsrc, DM dmdest)
343d71ae5a4SJacob Faibussowitsch {
34424989b8cSPeter Brune   PetscFunctionBegin;
34524989b8cSPeter Brune   PetscValidHeaderSpecific(dmsrc, DM_CLASSID, 1);
34624989b8cSPeter Brune   PetscValidHeaderSpecific(dmdest, DM_CLASSID, 2);
3479566063dSJacob Faibussowitsch   PetscCall(DMTSDestroy((DMTS *)&dmdest->dmts));
3482a34c10cSBarry Smith   dmdest->dmts = dmsrc->dmts;
3499566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference(dmdest->dmts));
3509566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dmdest, DMCoarsenHook_DMTS, DMRestrictHook_DMTS, NULL));
3519566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dmdest, DMSubDomainHook_DMTS, DMSubDomainRestrictHook_DMTS, NULL));
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35324989b8cSPeter Brune }
35424989b8cSPeter Brune 
35524989b8cSPeter Brune /*@C
356bcf0153eSBarry Smith   DMTSSetIFunction - set `TS` implicit function evaluation function
35724989b8cSPeter Brune 
35824989b8cSPeter Brune   Not Collective
35924989b8cSPeter Brune 
3604165533cSJose E. Roman   Input Parameters:
361bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
362a96d6ef6SBarry Smith . func - function evaluating f(t,u,u_t)
36324989b8cSPeter Brune - ctx  - context for residual evaluation
36424989b8cSPeter Brune 
36524989b8cSPeter Brune   Level: advanced
36624989b8cSPeter Brune 
36724989b8cSPeter Brune   Note:
368bcf0153eSBarry Smith   `TSSetFunction()` is normally used, but it calls this function internally because the user context is actually
369bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
370*14d0ab18SJacob Faibussowitsch   not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
37124989b8cSPeter Brune 
372*14d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `DM`, `TSIFunction`, `DMTSSetContext()`, `TSSetIFunction()`,
373*14d0ab18SJacob Faibussowitsch `DMTSSetJacobian()`
37424989b8cSPeter Brune @*/
375d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIFunction(DM dm, TSIFunction func, void *ctx)
376d71ae5a4SJacob Faibussowitsch {
377942e3340SBarry Smith   DMTS tsdm;
37824989b8cSPeter Brune 
37924989b8cSPeter Brune   PetscFunctionBegin;
38024989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3819566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
382d74926cbSBarry Smith   if (func) tsdm->ops->ifunction = func;
383800f99ffSJeremy L Thompson   if (ctx) {
384800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
385800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
386800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
387800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "ifunction ctx", (PetscObject)ctxcontainer));
388800f99ffSJeremy L Thompson     tsdm->ifunctionctxcontainer = ctxcontainer;
389800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
390800f99ffSJeremy L Thompson   }
3913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
392800f99ffSJeremy L Thompson }
393800f99ffSJeremy L Thompson 
394800f99ffSJeremy L Thompson /*@C
3955cb80ecdSBarry Smith   DMTSSetIFunctionContextDestroy - set `TS` implicit evaluation context destroy function
396800f99ffSJeremy L Thompson 
397800f99ffSJeremy L Thompson   Not Collective
398800f99ffSJeremy L Thompson 
399800f99ffSJeremy L Thompson   Input Parameters:
4005cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
4015cb80ecdSBarry Smith - f  - implicit evaluation context destroy function
402800f99ffSJeremy L Thompson 
403800f99ffSJeremy L Thompson   Level: advanced
404800f99ffSJeremy L Thompson 
4051cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetIFunction()`, `TSSetIFunction()`
406800f99ffSJeremy L Thompson @*/
407d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
408d71ae5a4SJacob Faibussowitsch {
409800f99ffSJeremy L Thompson   DMTS tsdm;
410800f99ffSJeremy L Thompson 
411800f99ffSJeremy L Thompson   PetscFunctionBegin;
412800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
413800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
414800f99ffSJeremy L Thompson   if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ifunctionctxcontainer, f));
4153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
416800f99ffSJeremy L Thompson }
417800f99ffSJeremy L Thompson 
418d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetIFunctionContext_Internal(DM dm)
419d71ae5a4SJacob Faibussowitsch {
420800f99ffSJeremy L Thompson   DMTS tsdm;
421800f99ffSJeremy L Thompson 
422800f99ffSJeremy L Thompson   PetscFunctionBegin;
423800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
424800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
425800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIFunctionContext_DMTS(tsdm));
4263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42724989b8cSPeter Brune }
42824989b8cSPeter Brune 
42924989b8cSPeter Brune /*@C
430bcf0153eSBarry Smith   DMTSGetIFunction - get `TS` implicit residual evaluation function
43124989b8cSPeter Brune 
43224989b8cSPeter Brune   Not Collective
43324989b8cSPeter Brune 
4344165533cSJose E. Roman   Input Parameter:
435bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
43624989b8cSPeter Brune 
4374165533cSJose E. Roman   Output Parameters:
43820f4b53cSBarry Smith + func - function evaluation function, for calling sequence see `TSSetIFunction()`
43924989b8cSPeter Brune - ctx  - context for residual evaluation
44024989b8cSPeter Brune 
44124989b8cSPeter Brune   Level: advanced
44224989b8cSPeter Brune 
44324989b8cSPeter Brune   Note:
444bcf0153eSBarry Smith   `TSGetFunction()` is normally used, but it calls this function internally because the user context is actually
445bcf0153eSBarry Smith   associated with the `DM`.
44624989b8cSPeter Brune 
4471cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `DMTSSetFunction()`, `TSSetFunction()`
44824989b8cSPeter Brune @*/
449d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetIFunction(DM dm, TSIFunction *func, void **ctx)
450d71ae5a4SJacob Faibussowitsch {
451942e3340SBarry Smith   DMTS tsdm;
45224989b8cSPeter Brune 
45324989b8cSPeter Brune   PetscFunctionBegin;
45424989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4559566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
456d74926cbSBarry Smith   if (func) *func = tsdm->ops->ifunction;
457800f99ffSJeremy L Thompson   if (ctx) {
458800f99ffSJeremy L Thompson     if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ifunctionctxcontainer, ctx));
459800f99ffSJeremy L Thompson     else *ctx = NULL;
460800f99ffSJeremy L Thompson   }
4613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46224989b8cSPeter Brune }
46324989b8cSPeter Brune 
464efe9872eSLisandro Dalcin /*@C
465bcf0153eSBarry Smith   DMTSSetI2Function - set `TS` implicit function evaluation function for 2nd order systems
466efe9872eSLisandro Dalcin 
467efe9872eSLisandro Dalcin   Not Collective
468efe9872eSLisandro Dalcin 
4694165533cSJose E. Roman   Input Parameters:
470bcf0153eSBarry Smith + dm  - `DM` to be used with `TS`
471a96d6ef6SBarry Smith . fun - function evaluation routine
472efe9872eSLisandro Dalcin - ctx - context for residual evaluation
473efe9872eSLisandro Dalcin 
474efe9872eSLisandro Dalcin   Level: advanced
475efe9872eSLisandro Dalcin 
476efe9872eSLisandro Dalcin   Note:
477bcf0153eSBarry Smith   `TSSetI2Function()` is normally used, but it calls this function internally because the user context is actually
478bcf0153eSBarry Smith   associated with the `DM`.
479efe9872eSLisandro Dalcin 
4801cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `TSSetI2Function()`
481efe9872eSLisandro Dalcin @*/
482d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetI2Function(DM dm, TSI2Function fun, void *ctx)
483d71ae5a4SJacob Faibussowitsch {
484efe9872eSLisandro Dalcin   DMTS tsdm;
485efe9872eSLisandro Dalcin 
486efe9872eSLisandro Dalcin   PetscFunctionBegin;
487efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4889566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
489efe9872eSLisandro Dalcin   if (fun) tsdm->ops->i2function = fun;
490800f99ffSJeremy L Thompson   if (ctx) {
491800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
492800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
493800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
494800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2function ctx", (PetscObject)ctxcontainer));
495800f99ffSJeremy L Thompson     tsdm->i2functionctxcontainer = ctxcontainer;
496800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
497800f99ffSJeremy L Thompson   }
4983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
499800f99ffSJeremy L Thompson }
500800f99ffSJeremy L Thompson 
501800f99ffSJeremy L Thompson /*@C
5025cb80ecdSBarry Smith   DMTSSetI2FunctionContextDestroy - set `TS` implicit evaluation for 2nd order systems context destroy
503800f99ffSJeremy L Thompson 
504800f99ffSJeremy L Thompson   Not Collective
505800f99ffSJeremy L Thompson 
506800f99ffSJeremy L Thompson   Input Parameters:
5075cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
5085cb80ecdSBarry Smith - f  - implicit evaluation context destroy function
509800f99ffSJeremy L Thompson 
510800f99ffSJeremy L Thompson   Level: advanced
511800f99ffSJeremy L Thompson 
512800f99ffSJeremy L Thompson   Note:
5135cb80ecdSBarry Smith   `TSSetI2FunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
5145cb80ecdSBarry Smith   associated with the `DM`.
515800f99ffSJeremy L Thompson 
5161cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetI2FunctionContextDestroy()`, `DMTSSetI2Function()`, `TSSetI2Function()`
517800f99ffSJeremy L Thompson @*/
518d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetI2FunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
519d71ae5a4SJacob Faibussowitsch {
520800f99ffSJeremy L Thompson   DMTS tsdm;
521800f99ffSJeremy L Thompson 
522800f99ffSJeremy L Thompson   PetscFunctionBegin;
523800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
524800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
525800f99ffSJeremy L Thompson   if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2functionctxcontainer, f));
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
527800f99ffSJeremy L Thompson }
528800f99ffSJeremy L Thompson 
529d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetI2FunctionContext_Internal(DM dm)
530d71ae5a4SJacob Faibussowitsch {
531800f99ffSJeremy L Thompson   DMTS tsdm;
532800f99ffSJeremy L Thompson 
533800f99ffSJeremy L Thompson   PetscFunctionBegin;
534800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
535800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
536800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2FunctionContext_DMTS(tsdm));
5373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
538efe9872eSLisandro Dalcin }
539efe9872eSLisandro Dalcin 
540efe9872eSLisandro Dalcin /*@C
541bcf0153eSBarry Smith   DMTSGetI2Function - get `TS` implicit residual evaluation function for 2nd order systems
542efe9872eSLisandro Dalcin 
543efe9872eSLisandro Dalcin   Not Collective
544efe9872eSLisandro Dalcin 
5454165533cSJose E. Roman   Input Parameter:
546bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
547efe9872eSLisandro Dalcin 
5484165533cSJose E. Roman   Output Parameters:
54920f4b53cSBarry Smith + fun - function evaluation function, for calling sequence see `TSSetI2Function()`
550efe9872eSLisandro Dalcin - ctx - context for residual evaluation
551efe9872eSLisandro Dalcin 
552efe9872eSLisandro Dalcin   Level: advanced
553efe9872eSLisandro Dalcin 
554efe9872eSLisandro Dalcin   Note:
555bcf0153eSBarry Smith   `TSGetI2Function()` is normally used, but it calls this function internally because the user context is actually
556bcf0153eSBarry Smith   associated with the `DM`.
557efe9872eSLisandro Dalcin 
5581cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetI2Function()`, `TSGetI2Function()`
559efe9872eSLisandro Dalcin @*/
560d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetI2Function(DM dm, TSI2Function *fun, void **ctx)
561d71ae5a4SJacob Faibussowitsch {
562efe9872eSLisandro Dalcin   DMTS tsdm;
563efe9872eSLisandro Dalcin 
564efe9872eSLisandro Dalcin   PetscFunctionBegin;
565efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5669566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
567efe9872eSLisandro Dalcin   if (fun) *fun = tsdm->ops->i2function;
568800f99ffSJeremy L Thompson   if (ctx) {
569800f99ffSJeremy L Thompson     if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2functionctxcontainer, ctx));
570800f99ffSJeremy L Thompson     else *ctx = NULL;
571800f99ffSJeremy L Thompson   }
5723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
573efe9872eSLisandro Dalcin }
574efe9872eSLisandro Dalcin 
575efe9872eSLisandro Dalcin /*@C
576bcf0153eSBarry Smith   DMTSSetI2Jacobian - set `TS` implicit Jacobian evaluation function for 2nd order systems
577efe9872eSLisandro Dalcin 
578efe9872eSLisandro Dalcin   Not Collective
579efe9872eSLisandro Dalcin 
5804165533cSJose E. Roman   Input Parameters:
581bcf0153eSBarry Smith + dm  - `DM` to be used with `TS`
582*14d0ab18SJacob Faibussowitsch . jac - Jacobian evaluation routine
583efe9872eSLisandro Dalcin - ctx - context for Jacobian evaluation
584efe9872eSLisandro Dalcin 
585efe9872eSLisandro Dalcin   Level: advanced
586efe9872eSLisandro Dalcin 
587efe9872eSLisandro Dalcin   Note:
588bcf0153eSBarry Smith   `TSSetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
589bcf0153eSBarry Smith   associated with the `DM`.
590efe9872eSLisandro Dalcin 
591*14d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `DM`, `TS`, `TSI2Jacobian`, `TSSetI2Jacobian()`
592efe9872eSLisandro Dalcin @*/
593d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetI2Jacobian(DM dm, TSI2Jacobian jac, void *ctx)
594d71ae5a4SJacob Faibussowitsch {
595efe9872eSLisandro Dalcin   DMTS tsdm;
596efe9872eSLisandro Dalcin 
597efe9872eSLisandro Dalcin   PetscFunctionBegin;
598efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5999566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
600efe9872eSLisandro Dalcin   if (jac) tsdm->ops->i2jacobian = jac;
601800f99ffSJeremy L Thompson   if (ctx) {
602800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
603800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
604800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
605800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2jacobian ctx", (PetscObject)ctxcontainer));
606800f99ffSJeremy L Thompson     tsdm->i2jacobianctxcontainer = ctxcontainer;
607800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
608800f99ffSJeremy L Thompson   }
6093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
610800f99ffSJeremy L Thompson }
611800f99ffSJeremy L Thompson 
612800f99ffSJeremy L Thompson /*@C
6135cb80ecdSBarry Smith   DMTSSetI2JacobianContextDestroy - set `TS` implicit Jacobian evaluation for 2nd order systems context destroy function
614800f99ffSJeremy L Thompson 
615800f99ffSJeremy L Thompson   Not Collective
616800f99ffSJeremy L Thompson 
617800f99ffSJeremy L Thompson   Input Parameters:
6185cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
6195cb80ecdSBarry Smith - f  - implicit Jacobian evaluation context destroy function
620800f99ffSJeremy L Thompson 
62187497f52SBarry Smith   Level: advanced
62287497f52SBarry Smith 
6231cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `TSSetI2JacobianContextDestroy()`, `DMTSSetI2Jacobian()`, `TSSetI2Jacobian()`
624800f99ffSJeremy L Thompson @*/
625d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetI2JacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
626d71ae5a4SJacob Faibussowitsch {
627800f99ffSJeremy L Thompson   DMTS tsdm;
628800f99ffSJeremy L Thompson 
629800f99ffSJeremy L Thompson   PetscFunctionBegin;
630800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
631800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
632800f99ffSJeremy L Thompson   if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2jacobianctxcontainer, f));
6333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
634800f99ffSJeremy L Thompson }
635800f99ffSJeremy L Thompson 
636d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM dm)
637d71ae5a4SJacob Faibussowitsch {
638800f99ffSJeremy L Thompson   DMTS tsdm;
639800f99ffSJeremy L Thompson 
640800f99ffSJeremy L Thompson   PetscFunctionBegin;
641800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
642800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
643800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2JacobianContext_DMTS(tsdm));
6443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
645efe9872eSLisandro Dalcin }
646efe9872eSLisandro Dalcin 
647efe9872eSLisandro Dalcin /*@C
648bcf0153eSBarry Smith   DMTSGetI2Jacobian - get `TS` implicit Jacobian evaluation function for 2nd order systems
649efe9872eSLisandro Dalcin 
650efe9872eSLisandro Dalcin   Not Collective
651efe9872eSLisandro Dalcin 
6524165533cSJose E. Roman   Input Parameter:
653bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
654efe9872eSLisandro Dalcin 
6554165533cSJose E. Roman   Output Parameters:
65620f4b53cSBarry Smith + jac - Jacobian evaluation function,  for calling sequence see `TSSetI2Jacobian()`
657efe9872eSLisandro Dalcin - ctx - context for Jacobian evaluation
658efe9872eSLisandro Dalcin 
659efe9872eSLisandro Dalcin   Level: advanced
660efe9872eSLisandro Dalcin 
661efe9872eSLisandro Dalcin   Note:
662bcf0153eSBarry Smith   `TSGetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
663bcf0153eSBarry Smith   associated with the `DM`.
664efe9872eSLisandro Dalcin 
6651cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetI2Jacobian()`, `TSGetI2Jacobian()`
666efe9872eSLisandro Dalcin @*/
667d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetI2Jacobian(DM dm, TSI2Jacobian *jac, void **ctx)
668d71ae5a4SJacob Faibussowitsch {
669efe9872eSLisandro Dalcin   DMTS tsdm;
670efe9872eSLisandro Dalcin 
671efe9872eSLisandro Dalcin   PetscFunctionBegin;
672efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6739566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
674efe9872eSLisandro Dalcin   if (jac) *jac = tsdm->ops->i2jacobian;
675800f99ffSJeremy L Thompson   if (ctx) {
676800f99ffSJeremy L Thompson     if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2jacobianctxcontainer, ctx));
677800f99ffSJeremy L Thompson     else *ctx = NULL;
678800f99ffSJeremy L Thompson   }
6793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
680efe9872eSLisandro Dalcin }
68124989b8cSPeter Brune 
68224989b8cSPeter Brune /*@C
683bcf0153eSBarry Smith   DMTSSetRHSFunction - set `TS` explicit residual evaluation function
68424989b8cSPeter Brune 
68524989b8cSPeter Brune   Not Collective
68624989b8cSPeter Brune 
6874165533cSJose E. Roman   Input Parameters:
688bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
689a96d6ef6SBarry Smith . func - RHS function evaluation routine
69024989b8cSPeter Brune - ctx  - context for residual evaluation
69124989b8cSPeter Brune 
69224989b8cSPeter Brune   Level: advanced
69324989b8cSPeter Brune 
69424989b8cSPeter Brune   Note:
695bcf0153eSBarry Smith   `TSSetRHSFunction()` is normally used, but it calls this function internally because the user context is actually
696bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
697bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
69824989b8cSPeter Brune 
699*14d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `DM`, `TS`, `TSRHSFunction`, `DMTSSetContext()`, `TSSetRHSFunction()`,
700*14d0ab18SJacob Faibussowitsch `DMTSSetJacobian()`
70124989b8cSPeter Brune @*/
702d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetRHSFunction(DM dm, TSRHSFunction func, void *ctx)
703d71ae5a4SJacob Faibussowitsch {
704942e3340SBarry Smith   DMTS tsdm;
70524989b8cSPeter Brune 
70624989b8cSPeter Brune   PetscFunctionBegin;
70724989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7089566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
709d74926cbSBarry Smith   if (func) tsdm->ops->rhsfunction = func;
710800f99ffSJeremy L Thompson   if (ctx) {
711800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
712800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
713800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
714800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs function ctx", (PetscObject)ctxcontainer));
715800f99ffSJeremy L Thompson     tsdm->rhsfunctionctxcontainer = ctxcontainer;
716800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
717800f99ffSJeremy L Thompson   }
7183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
719800f99ffSJeremy L Thompson }
720800f99ffSJeremy L Thompson 
721800f99ffSJeremy L Thompson /*@C
7225cb80ecdSBarry Smith   DMTSSetRHSFunctionContextDestroy - set `TS` explicit residual evaluation context destroy function
723800f99ffSJeremy L Thompson 
724800f99ffSJeremy L Thompson   Not Collective
725800f99ffSJeremy L Thompson 
726800f99ffSJeremy L Thompson   Input Parameters:
7275cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
7285cb80ecdSBarry Smith - f  - explicit evaluation context destroy function
729800f99ffSJeremy L Thompson 
730800f99ffSJeremy L Thompson   Level: advanced
731800f99ffSJeremy L Thompson 
732800f99ffSJeremy L Thompson   Note:
7335cb80ecdSBarry Smith   `TSSetRHSFunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
7345cb80ecdSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
7355cb80ecdSBarry Smith   not.
7365cb80ecdSBarry Smith 
737b43aa488SJacob Faibussowitsch   Developer Notes:
7385cb80ecdSBarry Smith   If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
739800f99ffSJeremy L Thompson 
7401cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetRHSFunctionContextDestroy()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
741800f99ffSJeremy L Thompson @*/
742d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
743d71ae5a4SJacob Faibussowitsch {
744800f99ffSJeremy L Thompson   DMTS tsdm;
745800f99ffSJeremy L Thompson 
746800f99ffSJeremy L Thompson   PetscFunctionBegin;
747800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
748800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
749800f99ffSJeremy L Thompson   if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsfunctionctxcontainer, f));
7503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
751800f99ffSJeremy L Thompson }
752800f99ffSJeremy L Thompson 
753d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM dm)
754d71ae5a4SJacob Faibussowitsch {
755800f99ffSJeremy L Thompson   DMTS tsdm;
756800f99ffSJeremy L Thompson 
757800f99ffSJeremy L Thompson   PetscFunctionBegin;
758800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
759800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
760800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetRHSFunctionContext_DMTS(tsdm));
761800f99ffSJeremy L Thompson   tsdm->rhsfunctionctxcontainer = NULL;
7623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76324989b8cSPeter Brune }
76424989b8cSPeter Brune 
765ef20d060SBarry Smith /*@C
766e3c11fc1SJed Brown   DMTSSetTransientVariable - sets function to transform from state to transient variables
767e3c11fc1SJed Brown 
768e3c11fc1SJed Brown   Logically Collective
769e3c11fc1SJed Brown 
7704165533cSJose E. Roman   Input Parameters:
771bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
772a96d6ef6SBarry Smith . tvar - a function that transforms to transient variables
773e3c11fc1SJed Brown - ctx  - a context for tvar
774e3c11fc1SJed Brown 
775e3c11fc1SJed Brown   Level: advanced
776e3c11fc1SJed Brown 
777e3c11fc1SJed Brown   Notes:
778bcf0153eSBarry Smith   This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., `TSBDF`)
779e3c11fc1SJed Brown   can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
780e3c11fc1SJed Brown   well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
781e3c11fc1SJed Brown   C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
782e3c11fc1SJed Brown   evaluated via the chain rule, as in
783e3c11fc1SJed Brown 
784e3c11fc1SJed Brown   dF/dP + shift * dF/dCdot dC/dP.
785e3c11fc1SJed Brown 
7861cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSBDF`, `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()`
787e3c11fc1SJed Brown @*/
788d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetTransientVariable(DM dm, TSTransientVariable tvar, void *ctx)
789d71ae5a4SJacob Faibussowitsch {
790e3c11fc1SJed Brown   DMTS dmts;
791e3c11fc1SJed Brown 
792e3c11fc1SJed Brown   PetscFunctionBegin;
793e3c11fc1SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7949566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &dmts));
795e3c11fc1SJed Brown   dmts->ops->transientvar = tvar;
796e3c11fc1SJed Brown   dmts->transientvarctx   = ctx;
7973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
798e3c11fc1SJed Brown }
799e3c11fc1SJed Brown 
800e3c11fc1SJed Brown /*@C
801bcf0153eSBarry Smith   DMTSGetTransientVariable - gets function to transform from state to transient variables set with `DMTSSetTransientVariable()`
802e3c11fc1SJed Brown 
803e3c11fc1SJed Brown   Logically Collective
804e3c11fc1SJed Brown 
8054165533cSJose E. Roman   Input Parameter:
806bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
807e3c11fc1SJed Brown 
8084165533cSJose E. Roman   Output Parameters:
809a96d6ef6SBarry Smith + tvar - a function that transforms to transient variables
810e3c11fc1SJed Brown - ctx  - a context for tvar
811e3c11fc1SJed Brown 
812e3c11fc1SJed Brown   Level: advanced
813e3c11fc1SJed Brown 
8141cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()`
815e3c11fc1SJed Brown @*/
816d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetTransientVariable(DM dm, TSTransientVariable *tvar, void *ctx)
817d71ae5a4SJacob Faibussowitsch {
818e3c11fc1SJed Brown   DMTS dmts;
819e3c11fc1SJed Brown 
820e3c11fc1SJed Brown   PetscFunctionBegin;
821e3c11fc1SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8229566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &dmts));
823e3c11fc1SJed Brown   if (tvar) *tvar = dmts->ops->transientvar;
824e3c11fc1SJed Brown   if (ctx) *(void **)ctx = dmts->transientvarctx;
8253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
826e3c11fc1SJed Brown }
827e3c11fc1SJed Brown 
828e3c11fc1SJed Brown /*@C
829bcf0153eSBarry Smith   DMTSGetSolutionFunction - gets the `TS` solution evaluation function
830ef20d060SBarry Smith 
831ef20d060SBarry Smith   Not Collective
832ef20d060SBarry Smith 
8334165533cSJose E. Roman   Input Parameter:
834bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
835ef20d060SBarry Smith 
836ef20d060SBarry Smith   Output Parameters:
83720f4b53cSBarry Smith + func - solution function evaluation function, for calling sequence see `TSSetSolution()`
838ef20d060SBarry Smith - ctx  - context for solution evaluation
839ef20d060SBarry Smith 
840ef20d060SBarry Smith   Level: advanced
841ef20d060SBarry Smith 
8421cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSSetSolutionFunction()`
843ef20d060SBarry Smith @*/
844d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetSolutionFunction(DM dm, TSSolutionFunction *func, void **ctx)
845d71ae5a4SJacob Faibussowitsch {
846942e3340SBarry Smith   DMTS tsdm;
847ef20d060SBarry Smith 
848ef20d060SBarry Smith   PetscFunctionBegin;
849ef20d060SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8509566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
851d74926cbSBarry Smith   if (func) *func = tsdm->ops->solution;
852ef20d060SBarry Smith   if (ctx) *ctx = tsdm->solutionctx;
8533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
854ef20d060SBarry Smith }
855ef20d060SBarry Smith 
856ef20d060SBarry Smith /*@C
857bcf0153eSBarry Smith   DMTSSetSolutionFunction - set `TS` solution evaluation function
858ef20d060SBarry Smith 
859ef20d060SBarry Smith   Not Collective
860ef20d060SBarry Smith 
8614165533cSJose E. Roman   Input Parameters:
862bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
863a96d6ef6SBarry Smith . func - solution function evaluation routine
864ef20d060SBarry Smith - ctx  - context for solution evaluation
865ef20d060SBarry Smith 
866ef20d060SBarry Smith   Level: advanced
867ef20d060SBarry Smith 
868ef20d060SBarry Smith   Note:
869bcf0153eSBarry Smith   `TSSetSolutionFunction()` is normally used, but it calls this function internally because the user context is actually
870bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
871bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
872ef20d060SBarry Smith 
8731cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSGetSolutionFunction()`
874ef20d060SBarry Smith @*/
875d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetSolutionFunction(DM dm, TSSolutionFunction func, void *ctx)
876d71ae5a4SJacob Faibussowitsch {
877942e3340SBarry Smith   DMTS tsdm;
878ef20d060SBarry Smith 
879ef20d060SBarry Smith   PetscFunctionBegin;
880ef20d060SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8819566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
882d74926cbSBarry Smith   if (func) tsdm->ops->solution = func;
883ef20d060SBarry Smith   if (ctx) tsdm->solutionctx = ctx;
8843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
885ef20d060SBarry Smith }
886ef20d060SBarry Smith 
8879b7cd975SBarry Smith /*@C
888bcf0153eSBarry Smith   DMTSSetForcingFunction - set `TS` forcing function evaluation function
8899b7cd975SBarry Smith 
8909b7cd975SBarry Smith   Not Collective
8919b7cd975SBarry Smith 
8924165533cSJose E. Roman   Input Parameters:
893bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
8942fe279fdSBarry Smith . func - forcing function evaluation routine
8959b7cd975SBarry Smith - ctx  - context for solution evaluation
8969b7cd975SBarry Smith 
8979b7cd975SBarry Smith   Level: advanced
8989b7cd975SBarry Smith 
8999b7cd975SBarry Smith   Note:
900bcf0153eSBarry Smith   `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually
901bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
902bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
9039b7cd975SBarry Smith 
904*14d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `DM`, `TS`, `TSForcingFunction`, `DMTSSetContext()`, `TSSetFunction()`,
905*14d0ab18SJacob Faibussowitsch `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
9069b7cd975SBarry Smith @*/
9072fe279fdSBarry Smith PetscErrorCode DMTSSetForcingFunction(DM dm, TSForcingFunction func, void *ctx)
908d71ae5a4SJacob Faibussowitsch {
9099b7cd975SBarry Smith   DMTS tsdm;
9109b7cd975SBarry Smith 
9119b7cd975SBarry Smith   PetscFunctionBegin;
9129b7cd975SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9139566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
9142fe279fdSBarry Smith   if (func) tsdm->ops->forcing = func;
9159b7cd975SBarry Smith   if (ctx) tsdm->forcingctx = ctx;
9163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9179b7cd975SBarry Smith }
9189b7cd975SBarry Smith 
9199b7cd975SBarry Smith /*@C
920bcf0153eSBarry Smith   DMTSGetForcingFunction - get `TS` forcing function evaluation function
9219b7cd975SBarry Smith 
9229b7cd975SBarry Smith   Not Collective
9239b7cd975SBarry Smith 
9244165533cSJose E. Roman   Input Parameter:
925bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
9269b7cd975SBarry Smith 
9274165533cSJose E. Roman   Output Parameters:
928bcf0153eSBarry Smith + f   - forcing function evaluation function; see `TSForcingFunction` for details
9299b7cd975SBarry Smith - ctx - context for solution evaluation
9309b7cd975SBarry Smith 
9319b7cd975SBarry Smith   Level: advanced
9329b7cd975SBarry Smith 
9339b7cd975SBarry Smith   Note:
934bcf0153eSBarry Smith   `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually
935bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
936bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
9379b7cd975SBarry Smith 
9381cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
9399b7cd975SBarry Smith @*/
940d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetForcingFunction(DM dm, TSForcingFunction *f, void **ctx)
941d71ae5a4SJacob Faibussowitsch {
9429b7cd975SBarry Smith   DMTS tsdm;
9439b7cd975SBarry Smith 
9449b7cd975SBarry Smith   PetscFunctionBegin;
9459b7cd975SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9469566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
947f8b49ee9SBarry Smith   if (f) *f = tsdm->ops->forcing;
9489b7cd975SBarry Smith   if (ctx) *ctx = tsdm->forcingctx;
9493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9509b7cd975SBarry Smith }
9519b7cd975SBarry Smith 
95224989b8cSPeter Brune /*@C
953bcf0153eSBarry Smith   DMTSGetRHSFunction - get `TS` explicit residual evaluation function
95424989b8cSPeter Brune 
95524989b8cSPeter Brune   Not Collective
95624989b8cSPeter Brune 
9574165533cSJose E. Roman   Input Parameter:
958bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
95924989b8cSPeter Brune 
9604165533cSJose E. Roman   Output Parameters:
96120f4b53cSBarry Smith + func - residual evaluation function, for calling sequence see `TSSetRHSFunction()`
96224989b8cSPeter Brune - ctx  - context for residual evaluation
96324989b8cSPeter Brune 
96424989b8cSPeter Brune   Level: advanced
96524989b8cSPeter Brune 
96624989b8cSPeter Brune   Note:
967bcf0153eSBarry Smith   `TSGetFunction()` is normally used, but it calls this function internally because the user context is actually
96824989b8cSPeter Brune   associated with the DM.
96924989b8cSPeter Brune 
9701cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
97124989b8cSPeter Brune @*/
972d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetRHSFunction(DM dm, TSRHSFunction *func, void **ctx)
973d71ae5a4SJacob Faibussowitsch {
974942e3340SBarry Smith   DMTS tsdm;
97524989b8cSPeter Brune 
97624989b8cSPeter Brune   PetscFunctionBegin;
97724989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9789566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
979d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsfunction;
980800f99ffSJeremy L Thompson   if (ctx) {
981800f99ffSJeremy L Thompson     if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsfunctionctxcontainer, ctx));
982800f99ffSJeremy L Thompson     else *ctx = NULL;
983800f99ffSJeremy L Thompson   }
9843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98524989b8cSPeter Brune }
98624989b8cSPeter Brune 
98724989b8cSPeter Brune /*@C
988bcf0153eSBarry Smith   DMTSSetIJacobian - set `TS` Jacobian evaluation function
98924989b8cSPeter Brune 
99024989b8cSPeter Brune   Not Collective
99124989b8cSPeter Brune 
9924165533cSJose E. Roman   Input Parameters:
993bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
994a96d6ef6SBarry Smith . func - Jacobian evaluation routine
99524989b8cSPeter Brune - ctx  - context for residual evaluation
99624989b8cSPeter Brune 
99724989b8cSPeter Brune   Level: advanced
99824989b8cSPeter Brune 
99924989b8cSPeter Brune   Note:
1000bcf0153eSBarry Smith   `TSSetJacobian()` is normally used, but it calls this function internally because the user context is actually
1001bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1002bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
100324989b8cSPeter Brune 
1004*14d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `DM`, `TSIJacobian`, `DMTSSetContext()`, `TSSetRHSFunction()`,
1005*14d0ab18SJacob Faibussowitsch `DMTSGetJacobian()`, `TSSetIJacobian()`, `TSSetIFunction()`
100624989b8cSPeter Brune @*/
1007d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIJacobian(DM dm, TSIJacobian func, void *ctx)
1008d71ae5a4SJacob Faibussowitsch {
1009800f99ffSJeremy L Thompson   DMTS tsdm;
101024989b8cSPeter Brune 
101124989b8cSPeter Brune   PetscFunctionBegin;
101224989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1013800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1014800f99ffSJeremy L Thompson   if (func) tsdm->ops->ijacobian = func;
1015800f99ffSJeremy L Thompson   if (ctx) {
1016800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
1017800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
1018800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
1019800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "ijacobian ctx", (PetscObject)ctxcontainer));
1020800f99ffSJeremy L Thompson     tsdm->ijacobianctxcontainer = ctxcontainer;
1021800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
1022800f99ffSJeremy L Thompson   }
10233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1024800f99ffSJeremy L Thompson }
1025800f99ffSJeremy L Thompson 
1026800f99ffSJeremy L Thompson /*@C
10275cb80ecdSBarry Smith   DMTSSetIJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function
1028800f99ffSJeremy L Thompson 
1029800f99ffSJeremy L Thompson   Not Collective
1030800f99ffSJeremy L Thompson 
1031800f99ffSJeremy L Thompson   Input Parameters:
10325cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
10335cb80ecdSBarry Smith - f  - Jacobian evaluation context destroy function
1034800f99ffSJeremy L Thompson 
1035800f99ffSJeremy L Thompson   Level: advanced
1036800f99ffSJeremy L Thompson 
1037800f99ffSJeremy L Thompson   Note:
10385cb80ecdSBarry Smith   `TSSetIJacobianContextDestroy()` is normally used, but it calls this function internally because the user context is actually
10395cb80ecdSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
10405cb80ecdSBarry Smith   not.
1041800f99ffSJeremy L Thompson 
1042b43aa488SJacob Faibussowitsch   Developer Notes:
10435cb80ecdSBarry Smith   If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
10445cb80ecdSBarry Smith 
10451cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetIJacobianContextDestroy()`, `TSSetI2JacobianContextDestroy()`, `DMTSSetIJacobian()`, `TSSetIJacobian()`
1046800f99ffSJeremy L Thompson @*/
1047d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
1048d71ae5a4SJacob Faibussowitsch {
1049800f99ffSJeremy L Thompson   DMTS tsdm;
1050800f99ffSJeremy L Thompson 
1051800f99ffSJeremy L Thompson   PetscFunctionBegin;
1052800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1053800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1054800f99ffSJeremy L Thompson   if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ijacobianctxcontainer, f));
10553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1056800f99ffSJeremy L Thompson }
1057800f99ffSJeremy L Thompson 
1058d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM dm)
1059d71ae5a4SJacob Faibussowitsch {
1060800f99ffSJeremy L Thompson   DMTS tsdm;
1061800f99ffSJeremy L Thompson 
1062800f99ffSJeremy L Thompson   PetscFunctionBegin;
1063800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1064800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1065800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIJacobianContext_DMTS(tsdm));
10663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106724989b8cSPeter Brune }
106824989b8cSPeter Brune 
106924989b8cSPeter Brune /*@C
1070bcf0153eSBarry Smith   DMTSGetIJacobian - get `TS` Jacobian evaluation function
107124989b8cSPeter Brune 
107224989b8cSPeter Brune   Not Collective
107324989b8cSPeter Brune 
10744165533cSJose E. Roman   Input Parameter:
1075bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
107624989b8cSPeter Brune 
10774165533cSJose E. Roman   Output Parameters:
107820f4b53cSBarry Smith + func - Jacobian evaluation function, for calling sequence see `TSSetIJacobian()`
107924989b8cSPeter Brune - ctx  - context for residual evaluation
108024989b8cSPeter Brune 
108124989b8cSPeter Brune   Level: advanced
108224989b8cSPeter Brune 
108324989b8cSPeter Brune   Note:
1084bcf0153eSBarry Smith   `TSGetJacobian()` is normally used, but it calls this function internally because the user context is actually
1085bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1086bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
108724989b8cSPeter Brune 
10881cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
108924989b8cSPeter Brune @*/
1090d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetIJacobian(DM dm, TSIJacobian *func, void **ctx)
1091d71ae5a4SJacob Faibussowitsch {
1092942e3340SBarry Smith   DMTS tsdm;
109324989b8cSPeter Brune 
109424989b8cSPeter Brune   PetscFunctionBegin;
109524989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10969566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
1097d74926cbSBarry Smith   if (func) *func = tsdm->ops->ijacobian;
1098800f99ffSJeremy L Thompson   if (ctx) {
1099800f99ffSJeremy L Thompson     if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ijacobianctxcontainer, ctx));
1100800f99ffSJeremy L Thompson     else *ctx = NULL;
1101800f99ffSJeremy L Thompson   }
11023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110324989b8cSPeter Brune }
110424989b8cSPeter Brune 
110524989b8cSPeter Brune /*@C
1106bcf0153eSBarry Smith   DMTSSetRHSJacobian - set `TS` Jacobian evaluation function
110724989b8cSPeter Brune 
110824989b8cSPeter Brune   Not Collective
110924989b8cSPeter Brune 
11104165533cSJose E. Roman   Input Parameters:
1111bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
1112a96d6ef6SBarry Smith . func - Jacobian evaluation routine
111324989b8cSPeter Brune - ctx  - context for residual evaluation
111424989b8cSPeter Brune 
111524989b8cSPeter Brune   Level: advanced
111624989b8cSPeter Brune 
111724989b8cSPeter Brune   Note:
1118bcf0153eSBarry Smith   `TSSetJacobian()` is normally used, but it calls this function internally because the user context is actually
1119bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1120bcf0153eSBarry Smith   not.
112124989b8cSPeter Brune 
1122b43aa488SJacob Faibussowitsch   Developer Notes:
1123bcf0153eSBarry Smith   If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1124bcf0153eSBarry Smith 
1125*14d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TSRHSJacobian`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSGetJacobian()`, `TSSetRHSJacobian()`
112624989b8cSPeter Brune @*/
1127d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetRHSJacobian(DM dm, TSRHSJacobian func, void *ctx)
1128d71ae5a4SJacob Faibussowitsch {
1129942e3340SBarry Smith   DMTS tsdm;
113024989b8cSPeter Brune 
113124989b8cSPeter Brune   PetscFunctionBegin;
113224989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11339566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1134d74926cbSBarry Smith   if (func) tsdm->ops->rhsjacobian = func;
1135800f99ffSJeremy L Thompson   if (ctx) {
1136800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
1137800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
1138800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
1139800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs jacobian ctx", (PetscObject)ctxcontainer));
1140800f99ffSJeremy L Thompson     tsdm->rhsjacobianctxcontainer = ctxcontainer;
1141800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
1142800f99ffSJeremy L Thompson   }
11433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1144800f99ffSJeremy L Thompson }
1145800f99ffSJeremy L Thompson 
1146800f99ffSJeremy L Thompson /*@C
11475cb80ecdSBarry Smith   DMTSSetRHSJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function
1148800f99ffSJeremy L Thompson 
1149800f99ffSJeremy L Thompson   Not Collective
1150800f99ffSJeremy L Thompson 
1151800f99ffSJeremy L Thompson   Input Parameters:
11525cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
11535cb80ecdSBarry Smith - f  - Jacobian evaluation context destroy function
11545cb80ecdSBarry Smith 
11555cb80ecdSBarry Smith   Level: advanced
11565cb80ecdSBarry Smith 
11575cb80ecdSBarry Smith   Note:
11585cb80ecdSBarry Smith   The user usually calls `TSSetRHSJacobianContextDestroy()` which calls this routine
1159800f99ffSJeremy L Thompson 
11601cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianContextDestroy()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
1161800f99ffSJeremy L Thompson @*/
1162d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
1163d71ae5a4SJacob Faibussowitsch {
1164800f99ffSJeremy L Thompson   DMTS tsdm;
1165800f99ffSJeremy L Thompson 
1166800f99ffSJeremy L Thompson   PetscFunctionBegin;
1167800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1168800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1169800f99ffSJeremy L Thompson   if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsjacobianctxcontainer, f));
11703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1171800f99ffSJeremy L Thompson }
1172800f99ffSJeremy L Thompson 
1173d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM dm)
1174d71ae5a4SJacob Faibussowitsch {
1175800f99ffSJeremy L Thompson   DMTS tsdm;
1176800f99ffSJeremy L Thompson 
1177800f99ffSJeremy L Thompson   PetscFunctionBegin;
1178800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1179800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1180800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetRHSJacobianContext_DMTS(tsdm));
11813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118224989b8cSPeter Brune }
118324989b8cSPeter Brune 
118424989b8cSPeter Brune /*@C
1185bcf0153eSBarry Smith   DMTSGetRHSJacobian - get `TS` Jacobian evaluation function
118624989b8cSPeter Brune 
118724989b8cSPeter Brune   Not Collective
118824989b8cSPeter Brune 
11894165533cSJose E. Roman   Input Parameter:
1190bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
119124989b8cSPeter Brune 
11924165533cSJose E. Roman   Output Parameters:
119320f4b53cSBarry Smith + func - Jacobian evaluation function, for calling sequence see `TSSetRHSJacobian()`
119424989b8cSPeter Brune - ctx  - context for residual evaluation
119524989b8cSPeter Brune 
119624989b8cSPeter Brune   Level: advanced
119724989b8cSPeter Brune 
119824989b8cSPeter Brune   Note:
1199bcf0153eSBarry Smith   `TSGetJacobian()` is normally used, but it calls this function internally because the user context is actually
1200bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1201bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
120224989b8cSPeter Brune 
12031cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
120424989b8cSPeter Brune @*/
1205d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetRHSJacobian(DM dm, TSRHSJacobian *func, void **ctx)
1206d71ae5a4SJacob Faibussowitsch {
1207942e3340SBarry Smith   DMTS tsdm;
120824989b8cSPeter Brune 
120924989b8cSPeter Brune   PetscFunctionBegin;
121024989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12119566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
1212d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsjacobian;
1213800f99ffSJeremy L Thompson   if (ctx) {
1214800f99ffSJeremy L Thompson     if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsjacobianctxcontainer, ctx));
1215800f99ffSJeremy L Thompson     else *ctx = NULL;
1216800f99ffSJeremy L Thompson   }
12173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121824989b8cSPeter Brune }
1219ad6bc421SBarry Smith 
1220ad6bc421SBarry Smith /*@C
1221ad6bc421SBarry Smith   DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
1222ad6bc421SBarry Smith 
1223ad6bc421SBarry Smith   Not Collective
1224ad6bc421SBarry Smith 
12254165533cSJose E. Roman   Input Parameters:
1226bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
1227ad6bc421SBarry Smith . view - viewer function
1228ad6bc421SBarry Smith - load - loading function
1229ad6bc421SBarry Smith 
1230ad6bc421SBarry Smith   Level: advanced
1231ad6bc421SBarry Smith 
12321cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1233ad6bc421SBarry Smith @*/
1234d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIFunctionSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1235d71ae5a4SJacob Faibussowitsch {
1236ad6bc421SBarry Smith   DMTS tsdm;
1237ad6bc421SBarry Smith 
1238ad6bc421SBarry Smith   PetscFunctionBegin;
1239ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12409566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1241ad6bc421SBarry Smith   tsdm->ops->ifunctionview = view;
1242ad6bc421SBarry Smith   tsdm->ops->ifunctionload = load;
12433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1244ad6bc421SBarry Smith }
1245ad6bc421SBarry Smith 
1246ad6bc421SBarry Smith /*@C
1247ad6bc421SBarry Smith   DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
1248ad6bc421SBarry Smith 
1249ad6bc421SBarry Smith   Not Collective
1250ad6bc421SBarry Smith 
12514165533cSJose E. Roman   Input Parameters:
1252bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
1253ad6bc421SBarry Smith . view - viewer function
1254ad6bc421SBarry Smith - load - loading function
1255ad6bc421SBarry Smith 
1256ad6bc421SBarry Smith   Level: advanced
1257ad6bc421SBarry Smith 
12581cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1259ad6bc421SBarry Smith @*/
1260d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIJacobianSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1261d71ae5a4SJacob Faibussowitsch {
1262ad6bc421SBarry Smith   DMTS tsdm;
1263ad6bc421SBarry Smith 
1264ad6bc421SBarry Smith   PetscFunctionBegin;
1265ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12669566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1267ad6bc421SBarry Smith   tsdm->ops->ijacobianview = view;
1268ad6bc421SBarry Smith   tsdm->ops->ijacobianload = load;
12693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1270ad6bc421SBarry Smith }
1271