xref: /petsc/src/ts/utils/dmts.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
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 
212*1cc06b55SBarry 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 
272*1cc06b55SBarry 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 
303*1cc06b55SBarry 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 
340*1cc06b55SBarry 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 
36520f4b53cSBarry Smith    Calling sequence of `func`:
36620f4b53cSBarry Smith $     PetscErrorCode func(TS ts, PetscReal t, Vec u, Vec u_t, Vec F,void *ctx);
36720f4b53cSBarry Smith +  ts  - the `TS` context obtained from `TSCreate()`
36820f4b53cSBarry Smith .  t   - time at step/stage being solved
369a96d6ef6SBarry Smith .  u   - state vector
370a96d6ef6SBarry Smith .  u_t - time derivative of state vector
371a96d6ef6SBarry Smith .  F   - function vector
372a96d6ef6SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
373a96d6ef6SBarry Smith 
37424989b8cSPeter Brune    Level: advanced
37524989b8cSPeter Brune 
37624989b8cSPeter Brune    Note:
377bcf0153eSBarry Smith    `TSSetFunction()` is normally used, but it calls this function internally because the user context is actually
378bcf0153eSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
37924989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
38024989b8cSPeter Brune 
381*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `TSSetIFunction()`, `DMTSSetJacobian()`
38224989b8cSPeter Brune @*/
383d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIFunction(DM dm, TSIFunction func, void *ctx)
384d71ae5a4SJacob Faibussowitsch {
385942e3340SBarry Smith   DMTS tsdm;
38624989b8cSPeter Brune 
38724989b8cSPeter Brune   PetscFunctionBegin;
38824989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3899566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
390d74926cbSBarry Smith   if (func) tsdm->ops->ifunction = func;
391800f99ffSJeremy L Thompson   if (ctx) {
392800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
393800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
394800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
395800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "ifunction ctx", (PetscObject)ctxcontainer));
396800f99ffSJeremy L Thompson     tsdm->ifunctionctxcontainer = ctxcontainer;
397800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
398800f99ffSJeremy L Thompson   }
3993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
400800f99ffSJeremy L Thompson }
401800f99ffSJeremy L Thompson 
402800f99ffSJeremy L Thompson /*@C
4035cb80ecdSBarry Smith    DMTSSetIFunctionContextDestroy - set `TS` implicit evaluation context destroy function
404800f99ffSJeremy L Thompson 
405800f99ffSJeremy L Thompson    Not Collective
406800f99ffSJeremy L Thompson 
407800f99ffSJeremy L Thompson    Input Parameters:
4085cb80ecdSBarry Smith +  dm - `DM` to be used with `TS`
4095cb80ecdSBarry Smith -  f - implicit evaluation context destroy function
410800f99ffSJeremy L Thompson 
411800f99ffSJeremy L Thompson    Level: advanced
412800f99ffSJeremy L Thompson 
413*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetIFunction()`, `TSSetIFunction()`
414800f99ffSJeremy L Thompson @*/
415d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
416d71ae5a4SJacob Faibussowitsch {
417800f99ffSJeremy L Thompson   DMTS tsdm;
418800f99ffSJeremy L Thompson 
419800f99ffSJeremy L Thompson   PetscFunctionBegin;
420800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
421800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
422800f99ffSJeremy L Thompson   if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ifunctionctxcontainer, f));
4233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
424800f99ffSJeremy L Thompson }
425800f99ffSJeremy L Thompson 
426d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetIFunctionContext_Internal(DM dm)
427d71ae5a4SJacob Faibussowitsch {
428800f99ffSJeremy L Thompson   DMTS tsdm;
429800f99ffSJeremy L Thompson 
430800f99ffSJeremy L Thompson   PetscFunctionBegin;
431800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
432800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
433800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIFunctionContext_DMTS(tsdm));
4343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43524989b8cSPeter Brune }
43624989b8cSPeter Brune 
43724989b8cSPeter Brune /*@C
438bcf0153eSBarry Smith    DMTSGetIFunction - get `TS` implicit residual evaluation function
43924989b8cSPeter Brune 
44024989b8cSPeter Brune    Not Collective
44124989b8cSPeter Brune 
4424165533cSJose E. Roman    Input Parameter:
443bcf0153eSBarry Smith .  dm - `DM` to be used with `TS`
44424989b8cSPeter Brune 
4454165533cSJose E. Roman    Output Parameters:
44620f4b53cSBarry Smith +  func - function evaluation function, for calling sequence see `TSSetIFunction()`
44724989b8cSPeter Brune -  ctx - context for residual evaluation
44824989b8cSPeter Brune 
44924989b8cSPeter Brune    Level: advanced
45024989b8cSPeter Brune 
45124989b8cSPeter Brune    Note:
452bcf0153eSBarry Smith    `TSGetFunction()` is normally used, but it calls this function internally because the user context is actually
453bcf0153eSBarry Smith    associated with the `DM`.
45424989b8cSPeter Brune 
455*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `DMTSSetFunction()`, `TSSetFunction()`
45624989b8cSPeter Brune @*/
457d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetIFunction(DM dm, TSIFunction *func, void **ctx)
458d71ae5a4SJacob Faibussowitsch {
459942e3340SBarry Smith   DMTS tsdm;
46024989b8cSPeter Brune 
46124989b8cSPeter Brune   PetscFunctionBegin;
46224989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4639566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
464d74926cbSBarry Smith   if (func) *func = tsdm->ops->ifunction;
465800f99ffSJeremy L Thompson   if (ctx) {
466800f99ffSJeremy L Thompson     if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ifunctionctxcontainer, ctx));
467800f99ffSJeremy L Thompson     else *ctx = NULL;
468800f99ffSJeremy L Thompson   }
4693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47024989b8cSPeter Brune }
47124989b8cSPeter Brune 
472efe9872eSLisandro Dalcin /*@C
473bcf0153eSBarry Smith    DMTSSetI2Function - set `TS` implicit function evaluation function for 2nd order systems
474efe9872eSLisandro Dalcin 
475efe9872eSLisandro Dalcin    Not Collective
476efe9872eSLisandro Dalcin 
4774165533cSJose E. Roman    Input Parameters:
478bcf0153eSBarry Smith +  dm - `DM` to be used with `TS`
479a96d6ef6SBarry Smith .  fun - function evaluation routine
480efe9872eSLisandro Dalcin -  ctx - context for residual evaluation
481efe9872eSLisandro Dalcin 
48220f4b53cSBarry Smith    Calling sequence of `fun`:
48320f4b53cSBarry Smith $  PetscErrorCode fun(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F,void *ctx);
48420f4b53cSBarry Smith +  ts  - the `TS` context obtained from `TSCreate()`
48520f4b53cSBarry Smith .  t    - time at step/stage being solved
486a96d6ef6SBarry Smith .  U    - state vector
487a96d6ef6SBarry Smith .  U_t  - time derivative of state vector
488a96d6ef6SBarry Smith .  U_tt - second time derivative of state vector
489a96d6ef6SBarry Smith .  F    - function vector
49020f4b53cSBarry Smith -  ctx  - [optional] user-defined context for matrix evaluation routine (may be `NULL`)
491a96d6ef6SBarry Smith 
492efe9872eSLisandro Dalcin    Level: advanced
493efe9872eSLisandro Dalcin 
494efe9872eSLisandro Dalcin    Note:
495bcf0153eSBarry Smith    `TSSetI2Function()` is normally used, but it calls this function internally because the user context is actually
496bcf0153eSBarry Smith    associated with the `DM`.
497efe9872eSLisandro Dalcin 
498*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `TSSetI2Function()`
499efe9872eSLisandro Dalcin @*/
500d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetI2Function(DM dm, TSI2Function fun, void *ctx)
501d71ae5a4SJacob Faibussowitsch {
502efe9872eSLisandro Dalcin   DMTS tsdm;
503efe9872eSLisandro Dalcin 
504efe9872eSLisandro Dalcin   PetscFunctionBegin;
505efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5069566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
507efe9872eSLisandro Dalcin   if (fun) tsdm->ops->i2function = fun;
508800f99ffSJeremy L Thompson   if (ctx) {
509800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
510800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
511800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
512800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2function ctx", (PetscObject)ctxcontainer));
513800f99ffSJeremy L Thompson     tsdm->i2functionctxcontainer = ctxcontainer;
514800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
515800f99ffSJeremy L Thompson   }
5163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
517800f99ffSJeremy L Thompson }
518800f99ffSJeremy L Thompson 
519800f99ffSJeremy L Thompson /*@C
5205cb80ecdSBarry Smith    DMTSSetI2FunctionContextDestroy - set `TS` implicit evaluation for 2nd order systems context destroy
521800f99ffSJeremy L Thompson 
522800f99ffSJeremy L Thompson    Not Collective
523800f99ffSJeremy L Thompson 
524800f99ffSJeremy L Thompson    Input Parameters:
5255cb80ecdSBarry Smith +  dm - `DM` to be used with `TS`
5265cb80ecdSBarry Smith -  f - implicit evaluation context destroy function
527800f99ffSJeremy L Thompson 
528800f99ffSJeremy L Thompson    Level: advanced
529800f99ffSJeremy L Thompson 
530800f99ffSJeremy L Thompson    Note:
5315cb80ecdSBarry Smith    `TSSetI2FunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
5325cb80ecdSBarry Smith    associated with the `DM`.
533800f99ffSJeremy L Thompson 
534*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetI2FunctionContextDestroy()`, `DMTSSetI2Function()`, `TSSetI2Function()`
535800f99ffSJeremy L Thompson @*/
536d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetI2FunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
537d71ae5a4SJacob Faibussowitsch {
538800f99ffSJeremy L Thompson   DMTS tsdm;
539800f99ffSJeremy L Thompson 
540800f99ffSJeremy L Thompson   PetscFunctionBegin;
541800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
542800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
543800f99ffSJeremy L Thompson   if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2functionctxcontainer, f));
5443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
545800f99ffSJeremy L Thompson }
546800f99ffSJeremy L Thompson 
547d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetI2FunctionContext_Internal(DM dm)
548d71ae5a4SJacob Faibussowitsch {
549800f99ffSJeremy L Thompson   DMTS tsdm;
550800f99ffSJeremy L Thompson 
551800f99ffSJeremy L Thompson   PetscFunctionBegin;
552800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
553800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
554800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2FunctionContext_DMTS(tsdm));
5553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
556efe9872eSLisandro Dalcin }
557efe9872eSLisandro Dalcin 
558efe9872eSLisandro Dalcin /*@C
559bcf0153eSBarry Smith    DMTSGetI2Function - get `TS` implicit residual evaluation function for 2nd order systems
560efe9872eSLisandro Dalcin 
561efe9872eSLisandro Dalcin    Not Collective
562efe9872eSLisandro Dalcin 
5634165533cSJose E. Roman    Input Parameter:
564bcf0153eSBarry Smith .  dm - `DM` to be used with `TS`
565efe9872eSLisandro Dalcin 
5664165533cSJose E. Roman    Output Parameters:
56720f4b53cSBarry Smith +  fun - function evaluation function, for calling sequence see `TSSetI2Function()`
568efe9872eSLisandro Dalcin -  ctx - context for residual evaluation
569efe9872eSLisandro Dalcin 
570efe9872eSLisandro Dalcin    Level: advanced
571efe9872eSLisandro Dalcin 
572efe9872eSLisandro Dalcin    Note:
573bcf0153eSBarry Smith    `TSGetI2Function()` is normally used, but it calls this function internally because the user context is actually
574bcf0153eSBarry Smith    associated with the `DM`.
575efe9872eSLisandro Dalcin 
576*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetI2Function()`, `TSGetI2Function()`
577efe9872eSLisandro Dalcin @*/
578d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetI2Function(DM dm, TSI2Function *fun, void **ctx)
579d71ae5a4SJacob Faibussowitsch {
580efe9872eSLisandro Dalcin   DMTS tsdm;
581efe9872eSLisandro Dalcin 
582efe9872eSLisandro Dalcin   PetscFunctionBegin;
583efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5849566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
585efe9872eSLisandro Dalcin   if (fun) *fun = tsdm->ops->i2function;
586800f99ffSJeremy L Thompson   if (ctx) {
587800f99ffSJeremy L Thompson     if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2functionctxcontainer, ctx));
588800f99ffSJeremy L Thompson     else *ctx = NULL;
589800f99ffSJeremy L Thompson   }
5903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
591efe9872eSLisandro Dalcin }
592efe9872eSLisandro Dalcin 
593efe9872eSLisandro Dalcin /*@C
594bcf0153eSBarry Smith    DMTSSetI2Jacobian - set `TS` implicit Jacobian evaluation function for 2nd order systems
595efe9872eSLisandro Dalcin 
596efe9872eSLisandro Dalcin    Not Collective
597efe9872eSLisandro Dalcin 
5984165533cSJose E. Roman    Input Parameters:
599bcf0153eSBarry Smith +  dm - `DM` to be used with `TS`
600a96d6ef6SBarry Smith .  fun - Jacobian evaluation routine
601efe9872eSLisandro Dalcin -  ctx - context for Jacobian evaluation
602efe9872eSLisandro Dalcin 
60320f4b53cSBarry Smith    Calling sequence of `jac`:
60420f4b53cSBarry Smith $    PetscErrorCode jac(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, PetscReal v, PetscReal a, Mat J, Mat P, void *ctx)
60520f4b53cSBarry Smith +  ts  - the `TS` context obtained from `TSCreate()`
60620f4b53cSBarry Smith .  t    - time at step/stage being solved
607a96d6ef6SBarry Smith .  U    - state vector
608a96d6ef6SBarry Smith .  U_t  - time derivative of state vector
609a96d6ef6SBarry Smith .  U_tt - second time derivative of state vector
610a96d6ef6SBarry Smith .  v    - shift for U_t
611a96d6ef6SBarry Smith .  a    - shift for U_tt
612a96d6ef6SBarry Smith .  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
61320f4b53cSBarry Smith .  P    - preconditioning matrix for J, may be same as `J`
614a96d6ef6SBarry Smith -  ctx  - [optional] user-defined context for matrix evaluation routine
615a96d6ef6SBarry Smith 
616efe9872eSLisandro Dalcin    Level: advanced
617efe9872eSLisandro Dalcin 
618efe9872eSLisandro Dalcin    Note:
619bcf0153eSBarry Smith    `TSSetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
620bcf0153eSBarry Smith    associated with the `DM`.
621efe9872eSLisandro Dalcin 
622*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `TSSetI2Jacobian()`
623efe9872eSLisandro Dalcin @*/
624d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetI2Jacobian(DM dm, TSI2Jacobian jac, void *ctx)
625d71ae5a4SJacob Faibussowitsch {
626efe9872eSLisandro Dalcin   DMTS tsdm;
627efe9872eSLisandro Dalcin 
628efe9872eSLisandro Dalcin   PetscFunctionBegin;
629efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6309566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
631efe9872eSLisandro Dalcin   if (jac) tsdm->ops->i2jacobian = jac;
632800f99ffSJeremy L Thompson   if (ctx) {
633800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
634800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
635800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
636800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2jacobian ctx", (PetscObject)ctxcontainer));
637800f99ffSJeremy L Thompson     tsdm->i2jacobianctxcontainer = ctxcontainer;
638800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
639800f99ffSJeremy L Thompson   }
6403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
641800f99ffSJeremy L Thompson }
642800f99ffSJeremy L Thompson 
643800f99ffSJeremy L Thompson /*@C
6445cb80ecdSBarry Smith    DMTSSetI2JacobianContextDestroy - set `TS` implicit Jacobian evaluation for 2nd order systems context destroy function
645800f99ffSJeremy L Thompson 
646800f99ffSJeremy L Thompson    Not Collective
647800f99ffSJeremy L Thompson 
648800f99ffSJeremy L Thompson    Input Parameters:
6495cb80ecdSBarry Smith +  dm - `DM` to be used with `TS`
6505cb80ecdSBarry Smith -  f - implicit Jacobian evaluation context destroy function
651800f99ffSJeremy L Thompson 
65287497f52SBarry Smith    Level: advanced
65387497f52SBarry Smith 
654*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `TSSetI2JacobianContextDestroy()`, `DMTSSetI2Jacobian()`, `TSSetI2Jacobian()`
655800f99ffSJeremy L Thompson @*/
656d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetI2JacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
657d71ae5a4SJacob Faibussowitsch {
658800f99ffSJeremy L Thompson   DMTS tsdm;
659800f99ffSJeremy L Thompson 
660800f99ffSJeremy L Thompson   PetscFunctionBegin;
661800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
662800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
663800f99ffSJeremy L Thompson   if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2jacobianctxcontainer, f));
6643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
665800f99ffSJeremy L Thompson }
666800f99ffSJeremy L Thompson 
667d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM dm)
668d71ae5a4SJacob Faibussowitsch {
669800f99ffSJeremy L Thompson   DMTS tsdm;
670800f99ffSJeremy L Thompson 
671800f99ffSJeremy L Thompson   PetscFunctionBegin;
672800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
673800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
674800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2JacobianContext_DMTS(tsdm));
6753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
676efe9872eSLisandro Dalcin }
677efe9872eSLisandro Dalcin 
678efe9872eSLisandro Dalcin /*@C
679bcf0153eSBarry Smith    DMTSGetI2Jacobian - get `TS` implicit Jacobian evaluation function for 2nd order systems
680efe9872eSLisandro Dalcin 
681efe9872eSLisandro Dalcin    Not Collective
682efe9872eSLisandro Dalcin 
6834165533cSJose E. Roman    Input Parameter:
684bcf0153eSBarry Smith .  dm - `DM` to be used with `TS`
685efe9872eSLisandro Dalcin 
6864165533cSJose E. Roman    Output Parameters:
68720f4b53cSBarry Smith +  jac - Jacobian evaluation function,  for calling sequence see `TSSetI2Jacobian()`
688efe9872eSLisandro Dalcin -  ctx - context for Jacobian evaluation
689efe9872eSLisandro Dalcin 
690efe9872eSLisandro Dalcin    Level: advanced
691efe9872eSLisandro Dalcin 
692efe9872eSLisandro Dalcin    Note:
693bcf0153eSBarry Smith    `TSGetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
694bcf0153eSBarry Smith    associated with the `DM`.
695efe9872eSLisandro Dalcin 
696*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetI2Jacobian()`, `TSGetI2Jacobian()`
697efe9872eSLisandro Dalcin @*/
698d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetI2Jacobian(DM dm, TSI2Jacobian *jac, void **ctx)
699d71ae5a4SJacob Faibussowitsch {
700efe9872eSLisandro Dalcin   DMTS tsdm;
701efe9872eSLisandro Dalcin 
702efe9872eSLisandro Dalcin   PetscFunctionBegin;
703efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7049566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
705efe9872eSLisandro Dalcin   if (jac) *jac = tsdm->ops->i2jacobian;
706800f99ffSJeremy L Thompson   if (ctx) {
707800f99ffSJeremy L Thompson     if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2jacobianctxcontainer, ctx));
708800f99ffSJeremy L Thompson     else *ctx = NULL;
709800f99ffSJeremy L Thompson   }
7103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
711efe9872eSLisandro Dalcin }
71224989b8cSPeter Brune 
71324989b8cSPeter Brune /*@C
714bcf0153eSBarry Smith    DMTSSetRHSFunction - set `TS` explicit residual evaluation function
71524989b8cSPeter Brune 
71624989b8cSPeter Brune    Not Collective
71724989b8cSPeter Brune 
7184165533cSJose E. Roman    Input Parameters:
719bcf0153eSBarry Smith +  dm - `DM` to be used with `TS`
720a96d6ef6SBarry Smith .  func - RHS function evaluation routine
72124989b8cSPeter Brune -  ctx - context for residual evaluation
72224989b8cSPeter Brune 
72320f4b53cSBarry Smith     Calling sequence of `func`:
724a96d6ef6SBarry Smith $   PetscErrorCode func(TS ts, PetscReal t, Vec u, Vec F, void *ctx);
725a96d6ef6SBarry Smith +   ts - timestep context
726a96d6ef6SBarry Smith .   t - current timestep
727a96d6ef6SBarry Smith .   u - input vector
728a96d6ef6SBarry Smith .   F - function vector
729a96d6ef6SBarry Smith -   ctx - [optional] user-defined function context
730a96d6ef6SBarry Smith 
73124989b8cSPeter Brune    Level: advanced
73224989b8cSPeter Brune 
73324989b8cSPeter Brune    Note:
734bcf0153eSBarry Smith    `TSSetRHSFunction()` is normally used, but it calls this function internally because the user context is actually
735bcf0153eSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
736bcf0153eSBarry Smith    not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
73724989b8cSPeter Brune 
738*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetJacobian()`
73924989b8cSPeter Brune @*/
740d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetRHSFunction(DM dm, TSRHSFunction func, void *ctx)
741d71ae5a4SJacob Faibussowitsch {
742942e3340SBarry Smith   DMTS tsdm;
74324989b8cSPeter Brune 
74424989b8cSPeter Brune   PetscFunctionBegin;
74524989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7469566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
747d74926cbSBarry Smith   if (func) tsdm->ops->rhsfunction = func;
748800f99ffSJeremy L Thompson   if (ctx) {
749800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
750800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
751800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
752800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs function ctx", (PetscObject)ctxcontainer));
753800f99ffSJeremy L Thompson     tsdm->rhsfunctionctxcontainer = ctxcontainer;
754800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
755800f99ffSJeremy L Thompson   }
7563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
757800f99ffSJeremy L Thompson }
758800f99ffSJeremy L Thompson 
759800f99ffSJeremy L Thompson /*@C
7605cb80ecdSBarry Smith    DMTSSetRHSFunctionContextDestroy - set `TS` explicit residual evaluation context destroy function
761800f99ffSJeremy L Thompson 
762800f99ffSJeremy L Thompson    Not Collective
763800f99ffSJeremy L Thompson 
764800f99ffSJeremy L Thompson    Input Parameters:
7655cb80ecdSBarry Smith +  dm - `DM` to be used with `TS`
7665cb80ecdSBarry Smith -  f - explicit evaluation context destroy function
767800f99ffSJeremy L Thompson 
768800f99ffSJeremy L Thompson    Level: advanced
769800f99ffSJeremy L Thompson 
770800f99ffSJeremy L Thompson    Note:
7715cb80ecdSBarry Smith    `TSSetRHSFunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
7725cb80ecdSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
7735cb80ecdSBarry Smith    not.
7745cb80ecdSBarry Smith 
7755cb80ecdSBarry Smith    Developer Note:
7765cb80ecdSBarry Smith    If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
777800f99ffSJeremy L Thompson 
778*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetRHSFunctionContextDestroy()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
779800f99ffSJeremy L Thompson @*/
780d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
781d71ae5a4SJacob Faibussowitsch {
782800f99ffSJeremy L Thompson   DMTS tsdm;
783800f99ffSJeremy L Thompson 
784800f99ffSJeremy L Thompson   PetscFunctionBegin;
785800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
786800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
787800f99ffSJeremy L Thompson   if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsfunctionctxcontainer, f));
7883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
789800f99ffSJeremy L Thompson }
790800f99ffSJeremy L Thompson 
791d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM dm)
792d71ae5a4SJacob Faibussowitsch {
793800f99ffSJeremy L Thompson   DMTS tsdm;
794800f99ffSJeremy L Thompson 
795800f99ffSJeremy L Thompson   PetscFunctionBegin;
796800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
797800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
798800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetRHSFunctionContext_DMTS(tsdm));
799800f99ffSJeremy L Thompson   tsdm->rhsfunctionctxcontainer = NULL;
8003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80124989b8cSPeter Brune }
80224989b8cSPeter Brune 
803ef20d060SBarry Smith /*@C
804e3c11fc1SJed Brown    DMTSSetTransientVariable - sets function to transform from state to transient variables
805e3c11fc1SJed Brown 
806e3c11fc1SJed Brown    Logically Collective
807e3c11fc1SJed Brown 
8084165533cSJose E. Roman    Input Parameters:
809bcf0153eSBarry Smith +  dm - `DM` to be used with `TS`
810a96d6ef6SBarry Smith .  tvar - a function that transforms to transient variables
811e3c11fc1SJed Brown -  ctx - a context for tvar
812e3c11fc1SJed Brown 
81320f4b53cSBarry Smith     Calling sequence of `tvar`:
814a96d6ef6SBarry Smith $   PetscErrorCode tvar(TS ts, Vec p, Vec c, void *ctx);
815a96d6ef6SBarry Smith +   ts - timestep context
8166aad120cSJose E. Roman .   p - input vector (primitive form)
817a96d6ef6SBarry Smith .   c - output vector, transient variables (conservative form)
818a96d6ef6SBarry Smith -   ctx - [optional] user-defined function context
819a96d6ef6SBarry Smith 
820e3c11fc1SJed Brown    Level: advanced
821e3c11fc1SJed Brown 
822e3c11fc1SJed Brown    Notes:
823bcf0153eSBarry Smith    This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., `TSBDF`)
824e3c11fc1SJed Brown    can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
825e3c11fc1SJed Brown    well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
826e3c11fc1SJed Brown    C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
827e3c11fc1SJed Brown    evaluated via the chain rule, as in
828e3c11fc1SJed Brown 
829e3c11fc1SJed Brown      dF/dP + shift * dF/dCdot dC/dP.
830e3c11fc1SJed Brown 
831*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSBDF`, `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()`
832e3c11fc1SJed Brown @*/
833d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetTransientVariable(DM dm, TSTransientVariable tvar, void *ctx)
834d71ae5a4SJacob Faibussowitsch {
835e3c11fc1SJed Brown   DMTS dmts;
836e3c11fc1SJed Brown 
837e3c11fc1SJed Brown   PetscFunctionBegin;
838e3c11fc1SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8399566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &dmts));
840e3c11fc1SJed Brown   dmts->ops->transientvar = tvar;
841e3c11fc1SJed Brown   dmts->transientvarctx   = ctx;
8423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
843e3c11fc1SJed Brown }
844e3c11fc1SJed Brown 
845e3c11fc1SJed Brown /*@C
846bcf0153eSBarry Smith    DMTSGetTransientVariable - gets function to transform from state to transient variables set with `DMTSSetTransientVariable()`
847e3c11fc1SJed Brown 
848e3c11fc1SJed Brown    Logically Collective
849e3c11fc1SJed Brown 
8504165533cSJose E. Roman    Input Parameter:
851bcf0153eSBarry Smith .  dm - `DM` to be used with `TS`
852e3c11fc1SJed Brown 
8534165533cSJose E. Roman    Output Parameters:
854a96d6ef6SBarry Smith +  tvar - a function that transforms to transient variables
855e3c11fc1SJed Brown -  ctx - a context for tvar
856e3c11fc1SJed Brown 
857e3c11fc1SJed Brown    Level: advanced
858e3c11fc1SJed Brown 
859*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()`
860e3c11fc1SJed Brown @*/
861d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetTransientVariable(DM dm, TSTransientVariable *tvar, void *ctx)
862d71ae5a4SJacob Faibussowitsch {
863e3c11fc1SJed Brown   DMTS dmts;
864e3c11fc1SJed Brown 
865e3c11fc1SJed Brown   PetscFunctionBegin;
866e3c11fc1SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8679566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &dmts));
868e3c11fc1SJed Brown   if (tvar) *tvar = dmts->ops->transientvar;
869e3c11fc1SJed Brown   if (ctx) *(void **)ctx = dmts->transientvarctx;
8703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
871e3c11fc1SJed Brown }
872e3c11fc1SJed Brown 
873e3c11fc1SJed Brown /*@C
874bcf0153eSBarry Smith    DMTSGetSolutionFunction - gets the `TS` solution evaluation function
875ef20d060SBarry Smith 
876ef20d060SBarry Smith    Not Collective
877ef20d060SBarry Smith 
8784165533cSJose E. Roman    Input Parameter:
879bcf0153eSBarry Smith .  dm - `DM` to be used with `TS`
880ef20d060SBarry Smith 
881ef20d060SBarry Smith    Output Parameters:
88220f4b53cSBarry Smith +  func - solution function evaluation function, for calling sequence see `TSSetSolution()`
883ef20d060SBarry Smith -  ctx - context for solution evaluation
884ef20d060SBarry Smith 
885ef20d060SBarry Smith    Level: advanced
886ef20d060SBarry Smith 
887*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSSetSolutionFunction()`
888ef20d060SBarry Smith @*/
889d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetSolutionFunction(DM dm, TSSolutionFunction *func, void **ctx)
890d71ae5a4SJacob Faibussowitsch {
891942e3340SBarry Smith   DMTS tsdm;
892ef20d060SBarry Smith 
893ef20d060SBarry Smith   PetscFunctionBegin;
894ef20d060SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8959566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
896d74926cbSBarry Smith   if (func) *func = tsdm->ops->solution;
897ef20d060SBarry Smith   if (ctx) *ctx = tsdm->solutionctx;
8983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
899ef20d060SBarry Smith }
900ef20d060SBarry Smith 
901ef20d060SBarry Smith /*@C
902bcf0153eSBarry Smith    DMTSSetSolutionFunction - set `TS` solution evaluation function
903ef20d060SBarry Smith 
904ef20d060SBarry Smith    Not Collective
905ef20d060SBarry Smith 
9064165533cSJose E. Roman    Input Parameters:
907bcf0153eSBarry Smith +  dm - `DM` to be used with `TS`
908a96d6ef6SBarry Smith .  func - solution function evaluation routine
909ef20d060SBarry Smith -  ctx - context for solution evaluation
910ef20d060SBarry Smith 
91120f4b53cSBarry Smith     Calling sequence of `f`:
912a96d6ef6SBarry Smith $   PetscErrorCode f(TS ts, PetscReal t, Vec u, void *ctx);
913a96d6ef6SBarry Smith +   ts - timestep context
914a96d6ef6SBarry Smith .   t - current timestep
915a96d6ef6SBarry Smith .   u - output vector
916a96d6ef6SBarry Smith -   ctx - [optional] user-defined function context
917a96d6ef6SBarry Smith 
918ef20d060SBarry Smith    Level: advanced
919ef20d060SBarry Smith 
920ef20d060SBarry Smith    Note:
921bcf0153eSBarry Smith    `TSSetSolutionFunction()` is normally used, but it calls this function internally because the user context is actually
922bcf0153eSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
923bcf0153eSBarry Smith    not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
924ef20d060SBarry Smith 
925*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSGetSolutionFunction()`
926ef20d060SBarry Smith @*/
927d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetSolutionFunction(DM dm, TSSolutionFunction func, void *ctx)
928d71ae5a4SJacob Faibussowitsch {
929942e3340SBarry Smith   DMTS tsdm;
930ef20d060SBarry Smith 
931ef20d060SBarry Smith   PetscFunctionBegin;
932ef20d060SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9339566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
934d74926cbSBarry Smith   if (func) tsdm->ops->solution = func;
935ef20d060SBarry Smith   if (ctx) tsdm->solutionctx = ctx;
9363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
937ef20d060SBarry Smith }
938ef20d060SBarry Smith 
9399b7cd975SBarry Smith /*@C
940bcf0153eSBarry Smith    DMTSSetForcingFunction - set `TS` forcing function evaluation function
9419b7cd975SBarry Smith 
9429b7cd975SBarry Smith    Not Collective
9439b7cd975SBarry Smith 
9444165533cSJose E. Roman    Input Parameters:
945bcf0153eSBarry Smith +  dm - `DM` to be used with `TS`
9462fe279fdSBarry Smith .  func - forcing function evaluation routine
9479b7cd975SBarry Smith -  ctx - context for solution evaluation
9489b7cd975SBarry Smith 
94920f4b53cSBarry Smith     Calling sequence of `func`:
95020f4b53cSBarry Smith $     PetscErrorCode func (TS ts, PetscReal t, Vec f,void *ctx)
951a96d6ef6SBarry Smith +   ts - timestep context
952a96d6ef6SBarry Smith .   t - current timestep
953a96d6ef6SBarry Smith .   f - output vector
954a96d6ef6SBarry Smith -   ctx - [optional] user-defined function context
955a96d6ef6SBarry Smith 
9569b7cd975SBarry Smith    Level: advanced
9579b7cd975SBarry Smith 
9589b7cd975SBarry Smith    Note:
959bcf0153eSBarry Smith    `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually
960bcf0153eSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
961bcf0153eSBarry Smith    not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
9629b7cd975SBarry Smith 
963*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
9649b7cd975SBarry Smith @*/
9652fe279fdSBarry Smith PetscErrorCode DMTSSetForcingFunction(DM dm, TSForcingFunction func, void *ctx)
966d71ae5a4SJacob Faibussowitsch {
9679b7cd975SBarry Smith   DMTS tsdm;
9689b7cd975SBarry Smith 
9699b7cd975SBarry Smith   PetscFunctionBegin;
9709b7cd975SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9719566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
9722fe279fdSBarry Smith   if (func) tsdm->ops->forcing = func;
9739b7cd975SBarry Smith   if (ctx) tsdm->forcingctx = ctx;
9743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9759b7cd975SBarry Smith }
9769b7cd975SBarry Smith 
9779b7cd975SBarry Smith /*@C
978bcf0153eSBarry Smith    DMTSGetForcingFunction - get `TS` forcing function evaluation function
9799b7cd975SBarry Smith 
9809b7cd975SBarry Smith    Not Collective
9819b7cd975SBarry Smith 
9824165533cSJose E. Roman    Input Parameter:
983bcf0153eSBarry Smith .   dm - `DM` to be used with `TS`
9849b7cd975SBarry Smith 
9854165533cSJose E. Roman    Output Parameters:
986bcf0153eSBarry Smith +  f - forcing function evaluation function; see `TSForcingFunction` for details
9879b7cd975SBarry Smith -  ctx - context for solution evaluation
9889b7cd975SBarry Smith 
9899b7cd975SBarry Smith    Level: advanced
9909b7cd975SBarry Smith 
9919b7cd975SBarry Smith    Note:
992bcf0153eSBarry Smith    `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually
993bcf0153eSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
994bcf0153eSBarry Smith    not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
9959b7cd975SBarry Smith 
996*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
9979b7cd975SBarry Smith @*/
998d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetForcingFunction(DM dm, TSForcingFunction *f, void **ctx)
999d71ae5a4SJacob Faibussowitsch {
10009b7cd975SBarry Smith   DMTS tsdm;
10019b7cd975SBarry Smith 
10029b7cd975SBarry Smith   PetscFunctionBegin;
10039b7cd975SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10049566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1005f8b49ee9SBarry Smith   if (f) *f = tsdm->ops->forcing;
10069b7cd975SBarry Smith   if (ctx) *ctx = tsdm->forcingctx;
10073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10089b7cd975SBarry Smith }
10099b7cd975SBarry Smith 
101024989b8cSPeter Brune /*@C
1011bcf0153eSBarry Smith    DMTSGetRHSFunction - get `TS` explicit residual evaluation function
101224989b8cSPeter Brune 
101324989b8cSPeter Brune    Not Collective
101424989b8cSPeter Brune 
10154165533cSJose E. Roman    Input Parameter:
1016bcf0153eSBarry Smith .  dm - `DM` to be used with `TS`
101724989b8cSPeter Brune 
10184165533cSJose E. Roman    Output Parameters:
101920f4b53cSBarry Smith +  func - residual evaluation function, for calling sequence see `TSSetRHSFunction()`
102024989b8cSPeter Brune -  ctx - context for residual evaluation
102124989b8cSPeter Brune 
102224989b8cSPeter Brune    Level: advanced
102324989b8cSPeter Brune 
102424989b8cSPeter Brune    Note:
1025bcf0153eSBarry Smith    `TSGetFunction()` is normally used, but it calls this function internally because the user context is actually
102624989b8cSPeter Brune    associated with the DM.
102724989b8cSPeter Brune 
1028*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
102924989b8cSPeter Brune @*/
1030d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetRHSFunction(DM dm, TSRHSFunction *func, void **ctx)
1031d71ae5a4SJacob Faibussowitsch {
1032942e3340SBarry Smith   DMTS tsdm;
103324989b8cSPeter Brune 
103424989b8cSPeter Brune   PetscFunctionBegin;
103524989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10369566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
1037d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsfunction;
1038800f99ffSJeremy L Thompson   if (ctx) {
1039800f99ffSJeremy L Thompson     if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsfunctionctxcontainer, ctx));
1040800f99ffSJeremy L Thompson     else *ctx = NULL;
1041800f99ffSJeremy L Thompson   }
10423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104324989b8cSPeter Brune }
104424989b8cSPeter Brune 
104524989b8cSPeter Brune /*@C
1046bcf0153eSBarry Smith    DMTSSetIJacobian - set `TS` Jacobian evaluation function
104724989b8cSPeter Brune 
104824989b8cSPeter Brune    Not Collective
104924989b8cSPeter Brune 
10504165533cSJose E. Roman    Input Parameters:
1051bcf0153eSBarry Smith +  dm - `DM` to be used with `TS`
1052a96d6ef6SBarry Smith .  func - Jacobian evaluation routine
105324989b8cSPeter Brune -  ctx - context for residual evaluation
105424989b8cSPeter Brune 
105520f4b53cSBarry Smith    Calling sequence of `f`:
1056a96d6ef6SBarry Smith $    PetscErrorCode f(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx);
105720f4b53cSBarry Smith +  ts  - the `TS` context obtained from `TSCreate()`
105820f4b53cSBarry Smith .  t    - time at step/stage being solved
1059a96d6ef6SBarry Smith .  U    - state vector
1060a96d6ef6SBarry Smith .  U_t  - time derivative of state vector
1061a96d6ef6SBarry Smith .  a    - shift
1062a96d6ef6SBarry Smith .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
1063a96d6ef6SBarry Smith .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
1064a96d6ef6SBarry Smith -  ctx  - [optional] user-defined context for matrix evaluation routine
1065a96d6ef6SBarry Smith 
106624989b8cSPeter Brune    Level: advanced
106724989b8cSPeter Brune 
106824989b8cSPeter Brune    Note:
1069bcf0153eSBarry Smith    `TSSetJacobian()` is normally used, but it calls this function internally because the user context is actually
1070bcf0153eSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1071bcf0153eSBarry Smith    not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
107224989b8cSPeter Brune 
1073*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSGetJacobian()`, `TSSetIJacobian()`, `TSSetIFunction()`
107424989b8cSPeter Brune @*/
1075d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIJacobian(DM dm, TSIJacobian func, void *ctx)
1076d71ae5a4SJacob Faibussowitsch {
1077800f99ffSJeremy L Thompson   DMTS tsdm;
107824989b8cSPeter Brune 
107924989b8cSPeter Brune   PetscFunctionBegin;
108024989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1081800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1082800f99ffSJeremy L Thompson   if (func) tsdm->ops->ijacobian = func;
1083800f99ffSJeremy L Thompson   if (ctx) {
1084800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
1085800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
1086800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
1087800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "ijacobian ctx", (PetscObject)ctxcontainer));
1088800f99ffSJeremy L Thompson     tsdm->ijacobianctxcontainer = ctxcontainer;
1089800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
1090800f99ffSJeremy L Thompson   }
10913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1092800f99ffSJeremy L Thompson }
1093800f99ffSJeremy L Thompson 
1094800f99ffSJeremy L Thompson /*@C
10955cb80ecdSBarry Smith    DMTSSetIJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function
1096800f99ffSJeremy L Thompson 
1097800f99ffSJeremy L Thompson    Not Collective
1098800f99ffSJeremy L Thompson 
1099800f99ffSJeremy L Thompson    Input Parameters:
11005cb80ecdSBarry Smith +  dm - `DM` to be used with `TS`
11015cb80ecdSBarry Smith -  f - Jacobian evaluation context destroy function
1102800f99ffSJeremy L Thompson 
1103800f99ffSJeremy L Thompson    Level: advanced
1104800f99ffSJeremy L Thompson 
1105800f99ffSJeremy L Thompson    Note:
11065cb80ecdSBarry Smith    `TSSetIJacobianContextDestroy()` is normally used, but it calls this function internally because the user context is actually
11075cb80ecdSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
11085cb80ecdSBarry Smith    not.
1109800f99ffSJeremy L Thompson 
11105cb80ecdSBarry Smith    Developer Note:
11115cb80ecdSBarry Smith    If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
11125cb80ecdSBarry Smith 
1113*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetIJacobianContextDestroy()`, `TSSetI2JacobianContextDestroy()`, `DMTSSetIJacobian()`, `TSSetIJacobian()`
1114800f99ffSJeremy L Thompson @*/
1115d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
1116d71ae5a4SJacob Faibussowitsch {
1117800f99ffSJeremy L Thompson   DMTS tsdm;
1118800f99ffSJeremy L Thompson 
1119800f99ffSJeremy L Thompson   PetscFunctionBegin;
1120800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1121800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1122800f99ffSJeremy L Thompson   if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ijacobianctxcontainer, f));
11233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1124800f99ffSJeremy L Thompson }
1125800f99ffSJeremy L Thompson 
1126d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM dm)
1127d71ae5a4SJacob Faibussowitsch {
1128800f99ffSJeremy L Thompson   DMTS tsdm;
1129800f99ffSJeremy L Thompson 
1130800f99ffSJeremy L Thompson   PetscFunctionBegin;
1131800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1132800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1133800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIJacobianContext_DMTS(tsdm));
11343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113524989b8cSPeter Brune }
113624989b8cSPeter Brune 
113724989b8cSPeter Brune /*@C
1138bcf0153eSBarry Smith    DMTSGetIJacobian - get `TS` Jacobian evaluation function
113924989b8cSPeter Brune 
114024989b8cSPeter Brune    Not Collective
114124989b8cSPeter Brune 
11424165533cSJose E. Roman    Input Parameter:
1143bcf0153eSBarry Smith .  dm - `DM` to be used with `TS`
114424989b8cSPeter Brune 
11454165533cSJose E. Roman    Output Parameters:
114620f4b53cSBarry Smith +  func - Jacobian evaluation function, for calling sequence see `TSSetIJacobian()`
114724989b8cSPeter Brune -  ctx - context for residual evaluation
114824989b8cSPeter Brune 
114924989b8cSPeter Brune    Level: advanced
115024989b8cSPeter Brune 
115124989b8cSPeter Brune    Note:
1152bcf0153eSBarry Smith    `TSGetJacobian()` is normally used, but it calls this function internally because the user context is actually
1153bcf0153eSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1154bcf0153eSBarry Smith    not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
115524989b8cSPeter Brune 
1156*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
115724989b8cSPeter Brune @*/
1158d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetIJacobian(DM dm, TSIJacobian *func, void **ctx)
1159d71ae5a4SJacob Faibussowitsch {
1160942e3340SBarry Smith   DMTS tsdm;
116124989b8cSPeter Brune 
116224989b8cSPeter Brune   PetscFunctionBegin;
116324989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11649566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
1165d74926cbSBarry Smith   if (func) *func = tsdm->ops->ijacobian;
1166800f99ffSJeremy L Thompson   if (ctx) {
1167800f99ffSJeremy L Thompson     if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ijacobianctxcontainer, ctx));
1168800f99ffSJeremy L Thompson     else *ctx = NULL;
1169800f99ffSJeremy L Thompson   }
11703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
117124989b8cSPeter Brune }
117224989b8cSPeter Brune 
117324989b8cSPeter Brune /*@C
1174bcf0153eSBarry Smith    DMTSSetRHSJacobian - set `TS` Jacobian evaluation function
117524989b8cSPeter Brune 
117624989b8cSPeter Brune    Not Collective
117724989b8cSPeter Brune 
11784165533cSJose E. Roman    Input Parameters:
1179bcf0153eSBarry Smith +  dm - `DM` to be used with `TS`
1180a96d6ef6SBarry Smith .  func - Jacobian evaluation routine
118124989b8cSPeter Brune -  ctx - context for residual evaluation
118224989b8cSPeter Brune 
118320f4b53cSBarry Smith    Calling sequence of `func`:
1184a96d6ef6SBarry Smith $     PetscErrorCode func(TS ts, PetscReal t, Vec u, Mat A, Mat B, void *ctx);
118520f4b53cSBarry Smith +  ts  - the `TS` context obtained from `TSCreate()`
118620f4b53cSBarry Smith .  t - current timestep
1187a96d6ef6SBarry Smith .  u - input vector
1188a96d6ef6SBarry Smith .  Amat - (approximate) Jacobian matrix
1189a96d6ef6SBarry Smith .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1190a96d6ef6SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
1191a96d6ef6SBarry Smith 
119224989b8cSPeter Brune    Level: advanced
119324989b8cSPeter Brune 
119424989b8cSPeter Brune    Note:
1195bcf0153eSBarry Smith    `TSSetJacobian()` is normally used, but it calls this function internally because the user context is actually
1196bcf0153eSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1197bcf0153eSBarry Smith    not.
119824989b8cSPeter Brune 
1199bcf0153eSBarry Smith    Developer Note:
1200bcf0153eSBarry Smith    If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1201bcf0153eSBarry Smith 
1202*1cc06b55SBarry Smith .seealso: [](ch_ts), `DMTSSetContext()`, `TSSetFunction()`, `DMTSGetJacobian()`, `TSSetRHSJacobian()`
120324989b8cSPeter Brune @*/
1204d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetRHSJacobian(DM dm, TSRHSJacobian func, void *ctx)
1205d71ae5a4SJacob Faibussowitsch {
1206942e3340SBarry Smith   DMTS tsdm;
120724989b8cSPeter Brune 
120824989b8cSPeter Brune   PetscFunctionBegin;
120924989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12109566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1211d74926cbSBarry Smith   if (func) tsdm->ops->rhsjacobian = func;
1212800f99ffSJeremy L Thompson   if (ctx) {
1213800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
1214800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
1215800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
1216800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs jacobian ctx", (PetscObject)ctxcontainer));
1217800f99ffSJeremy L Thompson     tsdm->rhsjacobianctxcontainer = ctxcontainer;
1218800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
1219800f99ffSJeremy L Thompson   }
12203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1221800f99ffSJeremy L Thompson }
1222800f99ffSJeremy L Thompson 
1223800f99ffSJeremy L Thompson /*@C
12245cb80ecdSBarry Smith    DMTSSetRHSJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function
1225800f99ffSJeremy L Thompson 
1226800f99ffSJeremy L Thompson    Not Collective
1227800f99ffSJeremy L Thompson 
1228800f99ffSJeremy L Thompson    Input Parameters:
12295cb80ecdSBarry Smith +  dm - `DM` to be used with `TS`
12305cb80ecdSBarry Smith -  f - Jacobian evaluation context destroy function
12315cb80ecdSBarry Smith 
12325cb80ecdSBarry Smith    Level: advanced
12335cb80ecdSBarry Smith 
12345cb80ecdSBarry Smith    Note:
12355cb80ecdSBarry Smith    The user usually calls `TSSetRHSJacobianContextDestroy()` which calls this routine
1236800f99ffSJeremy L Thompson 
1237*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianContextDestroy()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
1238800f99ffSJeremy L Thompson @*/
1239d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
1240d71ae5a4SJacob Faibussowitsch {
1241800f99ffSJeremy L Thompson   DMTS tsdm;
1242800f99ffSJeremy L Thompson 
1243800f99ffSJeremy L Thompson   PetscFunctionBegin;
1244800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1245800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1246800f99ffSJeremy L Thompson   if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsjacobianctxcontainer, f));
12473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1248800f99ffSJeremy L Thompson }
1249800f99ffSJeremy L Thompson 
1250d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM dm)
1251d71ae5a4SJacob Faibussowitsch {
1252800f99ffSJeremy L Thompson   DMTS tsdm;
1253800f99ffSJeremy L Thompson 
1254800f99ffSJeremy L Thompson   PetscFunctionBegin;
1255800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1256800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1257800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetRHSJacobianContext_DMTS(tsdm));
12583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125924989b8cSPeter Brune }
126024989b8cSPeter Brune 
126124989b8cSPeter Brune /*@C
1262bcf0153eSBarry Smith    DMTSGetRHSJacobian - get `TS` Jacobian evaluation function
126324989b8cSPeter Brune 
126424989b8cSPeter Brune    Not Collective
126524989b8cSPeter Brune 
12664165533cSJose E. Roman    Input Parameter:
1267bcf0153eSBarry Smith .  dm - `DM` to be used with `TS`
126824989b8cSPeter Brune 
12694165533cSJose E. Roman    Output Parameters:
127020f4b53cSBarry Smith +  func - Jacobian evaluation function, for calling sequence see `TSSetRHSJacobian()`
127124989b8cSPeter Brune -  ctx - context for residual evaluation
127224989b8cSPeter Brune 
127324989b8cSPeter Brune    Level: advanced
127424989b8cSPeter Brune 
127524989b8cSPeter Brune    Note:
1276bcf0153eSBarry Smith    `TSGetJacobian()` is normally used, but it calls this function internally because the user context is actually
1277bcf0153eSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1278bcf0153eSBarry Smith    not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
127924989b8cSPeter Brune 
1280*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
128124989b8cSPeter Brune @*/
1282d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetRHSJacobian(DM dm, TSRHSJacobian *func, void **ctx)
1283d71ae5a4SJacob Faibussowitsch {
1284942e3340SBarry Smith   DMTS tsdm;
128524989b8cSPeter Brune 
128624989b8cSPeter Brune   PetscFunctionBegin;
128724989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12889566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
1289d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsjacobian;
1290800f99ffSJeremy L Thompson   if (ctx) {
1291800f99ffSJeremy L Thompson     if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsjacobianctxcontainer, ctx));
1292800f99ffSJeremy L Thompson     else *ctx = NULL;
1293800f99ffSJeremy L Thompson   }
12943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
129524989b8cSPeter Brune }
1296ad6bc421SBarry Smith 
1297ad6bc421SBarry Smith /*@C
1298ad6bc421SBarry Smith    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
1299ad6bc421SBarry Smith 
1300ad6bc421SBarry Smith    Not Collective
1301ad6bc421SBarry Smith 
13024165533cSJose E. Roman    Input Parameters:
1303bcf0153eSBarry Smith +  dm - `DM` to be used with `TS`
1304ad6bc421SBarry Smith .  view - viewer function
1305ad6bc421SBarry Smith -  load - loading function
1306ad6bc421SBarry Smith 
1307ad6bc421SBarry Smith    Level: advanced
1308ad6bc421SBarry Smith 
1309*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1310ad6bc421SBarry Smith @*/
1311d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIFunctionSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1312d71ae5a4SJacob Faibussowitsch {
1313ad6bc421SBarry Smith   DMTS tsdm;
1314ad6bc421SBarry Smith 
1315ad6bc421SBarry Smith   PetscFunctionBegin;
1316ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13179566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1318ad6bc421SBarry Smith   tsdm->ops->ifunctionview = view;
1319ad6bc421SBarry Smith   tsdm->ops->ifunctionload = load;
13203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1321ad6bc421SBarry Smith }
1322ad6bc421SBarry Smith 
1323ad6bc421SBarry Smith /*@C
1324ad6bc421SBarry Smith    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
1325ad6bc421SBarry Smith 
1326ad6bc421SBarry Smith    Not Collective
1327ad6bc421SBarry Smith 
13284165533cSJose E. Roman    Input Parameters:
1329bcf0153eSBarry Smith +  dm - `DM` to be used with `TS`
1330ad6bc421SBarry Smith .  view - viewer function
1331ad6bc421SBarry Smith -  load - loading function
1332ad6bc421SBarry Smith 
1333ad6bc421SBarry Smith    Level: advanced
1334ad6bc421SBarry Smith 
1335*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1336ad6bc421SBarry Smith @*/
1337d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIJacobianSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1338d71ae5a4SJacob Faibussowitsch {
1339ad6bc421SBarry Smith   DMTS tsdm;
1340ad6bc421SBarry Smith 
1341ad6bc421SBarry Smith   PetscFunctionBegin;
1342ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13439566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1344ad6bc421SBarry Smith   tsdm->ops->ijacobianview = view;
1345ad6bc421SBarry Smith   tsdm->ops->ijacobianload = load;
13463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1347ad6bc421SBarry Smith }
1348