xref: /petsc/src/ts/utils/dmts.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
2af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
324989b8cSPeter Brune 
DMTSUnsetRHSFunctionContext_DMTS(DMTS tsdm)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 
DMTSUnsetRHSJacobianContext_DMTS(DMTS tsdm)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 
DMTSUnsetIFunctionContext_DMTS(DMTS tsdm)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 
DMTSUnsetIJacobianContext_DMTS(DMTS tsdm)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 
DMTSUnsetI2FunctionContext_DMTS(DMTS tsdm)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 
DMTSUnsetI2JacobianContext_DMTS(DMTS tsdm)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 
DMTSDestroy(DMTS * kdm)52d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSDestroy(DMTS *kdm)
53d71ae5a4SJacob Faibussowitsch {
54d74926cbSBarry Smith   PetscFunctionBegin;
553ba16761SJacob Faibussowitsch   if (!*kdm) PetscFunctionReturn(PETSC_SUCCESS);
56f4f49eeaSPierre Jolivet   PetscValidHeaderSpecific(*kdm, DMTS_CLASSID, 1);
57f4f49eeaSPierre Jolivet   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 
DMTSLoad(DMTS kdm,PetscViewer viewer)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) {
79*2a8381b2SBarry Smith     PetscCtx 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) {
88*2a8381b2SBarry Smith     PetscCtx 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 
DMTSView(DMTS kdm,PetscViewer viewer)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 {
1148434afd1SBarry Smith       TSIFunctionFn *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 {
1238434afd1SBarry Smith       TSIJacobianFn *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) {
139*2a8381b2SBarry Smith       PetscCtx 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) {
151*2a8381b2SBarry Smith       PetscCtx 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 
DMTSCreate(MPI_Comm comm,DMTS * kdm)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  */
DMCoarsenHook_DMTS(DM dm,DM dmc,PetscCtx ctx)171*2a8381b2SBarry Smith static PetscErrorCode DMCoarsenHook_DMTS(DM dm, DM dmc, PetscCtx 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  */
DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,PetscCtx ctx)180*2a8381b2SBarry Smith static PetscErrorCode DMRestrictHook_DMTS(DM dm, Mat Restrict, Vec rscale, Mat Inject, DM dmc, PetscCtx ctx)
181d71ae5a4SJacob Faibussowitsch {
18224989b8cSPeter Brune   PetscFunctionBegin;
1833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18424989b8cSPeter Brune }
18524989b8cSPeter Brune 
DMSubDomainHook_DMTS(DM dm,DM subdm,PetscCtx ctx)186*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainHook_DMTS(DM dm, DM subdm, PetscCtx 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  */
DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,PetscCtx ctx)195*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, PetscCtx 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 @*/
DMTSCopy(DMTS kdm,DMTS nkdm)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 @*/
DMGetDMTS(DM dm,DMTS * tsdm)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 @*/
DMGetDMTSWrite(DM dm,DMTS * tsdm)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
327d1c5d1fcSBarry Smith   DMCopyDMTS - copies a `DMTS` 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 @*/
DMCopyDMTS(DM dmsrc,DM dmdest)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
356d1c5d1fcSBarry Smith   DMTSSetIFunction - set `TS` implicit function evaluation function into a `DMTS`
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 
365d1c5d1fcSBarry Smith   Level: developer
36624989b8cSPeter Brune 
36724989b8cSPeter Brune   Note:
368346ce620SStefano Zampini   `TSSetIFunction()` 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
37014d0ab18SJacob 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 
3728434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `TSIFunctionFn`
37324989b8cSPeter Brune @*/
DMTSSetIFunction(DM dm,TSIFunctionFn * func,PetscCtx ctx)374*2a8381b2SBarry Smith PetscErrorCode DMTSSetIFunction(DM dm, TSIFunctionFn *func, PetscCtx ctx)
375d71ae5a4SJacob Faibussowitsch {
376942e3340SBarry Smith   DMTS tsdm;
37724989b8cSPeter Brune 
37824989b8cSPeter Brune   PetscFunctionBegin;
37924989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3809566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
381d74926cbSBarry Smith   if (func) tsdm->ops->ifunction = func;
382800f99ffSJeremy L Thompson   if (ctx) {
383800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
384800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
385800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
386800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "ifunction ctx", (PetscObject)ctxcontainer));
387800f99ffSJeremy L Thompson     tsdm->ifunctionctxcontainer = ctxcontainer;
388800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
389800f99ffSJeremy L Thompson   }
3903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
391800f99ffSJeremy L Thompson }
392800f99ffSJeremy L Thompson 
393800f99ffSJeremy L Thompson /*@C
394d1c5d1fcSBarry Smith   DMTSSetIFunctionContextDestroy - set `TS` implicit evaluation context destroy function into a `DMTS`
395800f99ffSJeremy L Thompson 
396800f99ffSJeremy L Thompson   Not Collective
397800f99ffSJeremy L Thompson 
398800f99ffSJeremy L Thompson   Input Parameters:
3995cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
40049abdd8aSBarry Smith - f  - implicit evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
401800f99ffSJeremy L Thompson 
402d1c5d1fcSBarry Smith   Level: developer
403800f99ffSJeremy L Thompson 
40449abdd8aSBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetIFunction()`, `TSSetIFunction()`, `PetscCtxDestroyFn`
405800f99ffSJeremy L Thompson @*/
DMTSSetIFunctionContextDestroy(DM dm,PetscCtxDestroyFn * f)40649abdd8aSBarry Smith PetscErrorCode DMTSSetIFunctionContextDestroy(DM dm, PetscCtxDestroyFn *f)
407d71ae5a4SJacob Faibussowitsch {
408800f99ffSJeremy L Thompson   DMTS tsdm;
409800f99ffSJeremy L Thompson 
410800f99ffSJeremy L Thompson   PetscFunctionBegin;
411800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
412800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
41349abdd8aSBarry Smith   if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->ifunctionctxcontainer, f));
4143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
415800f99ffSJeremy L Thompson }
416800f99ffSJeremy L Thompson 
DMTSUnsetIFunctionContext_Internal(DM dm)417d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetIFunctionContext_Internal(DM dm)
418d71ae5a4SJacob Faibussowitsch {
419800f99ffSJeremy L Thompson   DMTS tsdm;
420800f99ffSJeremy L Thompson 
421800f99ffSJeremy L Thompson   PetscFunctionBegin;
422800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
423800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
424800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIFunctionContext_DMTS(tsdm));
4253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42624989b8cSPeter Brune }
42724989b8cSPeter Brune 
42824989b8cSPeter Brune /*@C
429d1c5d1fcSBarry Smith   DMTSGetIFunction - get `TS` implicit residual evaluation function from a `DMTS`
43024989b8cSPeter Brune 
43124989b8cSPeter Brune   Not Collective
43224989b8cSPeter Brune 
4334165533cSJose E. Roman   Input Parameter:
434bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
43524989b8cSPeter Brune 
4364165533cSJose E. Roman   Output Parameters:
4378434afd1SBarry Smith + func - function evaluation function, for calling sequence see `TSIFunctionFn`
43824989b8cSPeter Brune - ctx  - context for residual evaluation
43924989b8cSPeter Brune 
440d1c5d1fcSBarry Smith   Level: developer
44124989b8cSPeter Brune 
44224989b8cSPeter Brune   Note:
443346ce620SStefano Zampini   `TSGetIFunction()` is normally used, but it calls this function internally because the user context is actually
444bcf0153eSBarry Smith   associated with the `DM`.
44524989b8cSPeter Brune 
4468434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `DMTSSetIFunction()`, `TSIFunctionFn`
44724989b8cSPeter Brune @*/
DMTSGetIFunction(DM dm,TSIFunctionFn ** func,PetscCtxRt ctx)448*2a8381b2SBarry Smith PetscErrorCode DMTSGetIFunction(DM dm, TSIFunctionFn **func, PetscCtxRt ctx)
449d71ae5a4SJacob Faibussowitsch {
450942e3340SBarry Smith   DMTS tsdm;
45124989b8cSPeter Brune 
45224989b8cSPeter Brune   PetscFunctionBegin;
45324989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4549566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
455d74926cbSBarry Smith   if (func) *func = tsdm->ops->ifunction;
456800f99ffSJeremy L Thompson   if (ctx) {
457800f99ffSJeremy L Thompson     if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ifunctionctxcontainer, ctx));
458*2a8381b2SBarry Smith     else *(void **)ctx = NULL;
459800f99ffSJeremy L Thompson   }
4603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46124989b8cSPeter Brune }
46224989b8cSPeter Brune 
463efe9872eSLisandro Dalcin /*@C
464d1c5d1fcSBarry Smith   DMTSSetI2Function - set `TS` implicit function evaluation function for 2nd order systems into a `TSDM`
465efe9872eSLisandro Dalcin 
466efe9872eSLisandro Dalcin   Not Collective
467efe9872eSLisandro Dalcin 
4684165533cSJose E. Roman   Input Parameters:
469bcf0153eSBarry Smith + dm  - `DM` to be used with `TS`
470a96d6ef6SBarry Smith . fun - function evaluation routine
471efe9872eSLisandro Dalcin - ctx - context for residual evaluation
472efe9872eSLisandro Dalcin 
473d1c5d1fcSBarry Smith   Level: developer
474efe9872eSLisandro Dalcin 
475efe9872eSLisandro Dalcin   Note:
476bcf0153eSBarry Smith   `TSSetI2Function()` is normally used, but it calls this function internally because the user context is actually
477bcf0153eSBarry Smith   associated with the `DM`.
478efe9872eSLisandro Dalcin 
479d1c5d1fcSBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSSetI2Function()`
480efe9872eSLisandro Dalcin @*/
DMTSSetI2Function(DM dm,TSI2FunctionFn * fun,PetscCtx ctx)481*2a8381b2SBarry Smith PetscErrorCode DMTSSetI2Function(DM dm, TSI2FunctionFn *fun, PetscCtx ctx)
482d71ae5a4SJacob Faibussowitsch {
483efe9872eSLisandro Dalcin   DMTS tsdm;
484efe9872eSLisandro Dalcin 
485efe9872eSLisandro Dalcin   PetscFunctionBegin;
486efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4879566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
488efe9872eSLisandro Dalcin   if (fun) tsdm->ops->i2function = fun;
489800f99ffSJeremy L Thompson   if (ctx) {
490800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
491800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
492800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
493800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2function ctx", (PetscObject)ctxcontainer));
494800f99ffSJeremy L Thompson     tsdm->i2functionctxcontainer = ctxcontainer;
495800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
496800f99ffSJeremy L Thompson   }
4973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
498800f99ffSJeremy L Thompson }
499800f99ffSJeremy L Thompson 
500800f99ffSJeremy L Thompson /*@C
501d1c5d1fcSBarry Smith   DMTSSetI2FunctionContextDestroy - set `TS` implicit evaluation for 2nd order systems context destroy into a `DMTS`
502800f99ffSJeremy L Thompson 
503800f99ffSJeremy L Thompson   Not Collective
504800f99ffSJeremy L Thompson 
505800f99ffSJeremy L Thompson   Input Parameters:
5065cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
50749abdd8aSBarry Smith - f  - implicit evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
508800f99ffSJeremy L Thompson 
509d1c5d1fcSBarry Smith   Level: developer
510800f99ffSJeremy L Thompson 
511800f99ffSJeremy L Thompson   Note:
5125cb80ecdSBarry Smith   `TSSetI2FunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
5135cb80ecdSBarry Smith   associated with the `DM`.
514800f99ffSJeremy L Thompson 
515d1c5d1fcSBarry Smith .seealso: [](ch_ts), `DMTS`, `TSSetI2FunctionContextDestroy()`, `DMTSSetI2Function()`, `TSSetI2Function()`
516800f99ffSJeremy L Thompson @*/
DMTSSetI2FunctionContextDestroy(DM dm,PetscCtxDestroyFn * f)51749abdd8aSBarry Smith PetscErrorCode DMTSSetI2FunctionContextDestroy(DM dm, PetscCtxDestroyFn *f)
518d71ae5a4SJacob Faibussowitsch {
519800f99ffSJeremy L Thompson   DMTS tsdm;
520800f99ffSJeremy L Thompson 
521800f99ffSJeremy L Thompson   PetscFunctionBegin;
522800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
523800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
52449abdd8aSBarry Smith   if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->i2functionctxcontainer, f));
5253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
526800f99ffSJeremy L Thompson }
527800f99ffSJeremy L Thompson 
DMTSUnsetI2FunctionContext_Internal(DM dm)528d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetI2FunctionContext_Internal(DM dm)
529d71ae5a4SJacob Faibussowitsch {
530800f99ffSJeremy L Thompson   DMTS tsdm;
531800f99ffSJeremy L Thompson 
532800f99ffSJeremy L Thompson   PetscFunctionBegin;
533800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
534800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
535800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2FunctionContext_DMTS(tsdm));
5363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
537efe9872eSLisandro Dalcin }
538efe9872eSLisandro Dalcin 
539efe9872eSLisandro Dalcin /*@C
540d1c5d1fcSBarry Smith   DMTSGetI2Function - get `TS` implicit residual evaluation function for 2nd order systems from a `DMTS`
541efe9872eSLisandro Dalcin 
542efe9872eSLisandro Dalcin   Not Collective
543efe9872eSLisandro Dalcin 
5444165533cSJose E. Roman   Input Parameter:
545bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
546efe9872eSLisandro Dalcin 
5474165533cSJose E. Roman   Output Parameters:
54820f4b53cSBarry Smith + fun - function evaluation function, for calling sequence see `TSSetI2Function()`
549efe9872eSLisandro Dalcin - ctx - context for residual evaluation
550efe9872eSLisandro Dalcin 
551d1c5d1fcSBarry Smith   Level: developer
552efe9872eSLisandro Dalcin 
553efe9872eSLisandro Dalcin   Note:
554bcf0153eSBarry Smith   `TSGetI2Function()` is normally used, but it calls this function internally because the user context is actually
555bcf0153eSBarry Smith   associated with the `DM`.
556efe9872eSLisandro Dalcin 
557d1c5d1fcSBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetI2Function()`, `TSGetI2Function()`
558efe9872eSLisandro Dalcin @*/
DMTSGetI2Function(DM dm,TSI2FunctionFn ** fun,PetscCtxRt ctx)559*2a8381b2SBarry Smith PetscErrorCode DMTSGetI2Function(DM dm, TSI2FunctionFn **fun, PetscCtxRt ctx)
560d71ae5a4SJacob Faibussowitsch {
561efe9872eSLisandro Dalcin   DMTS tsdm;
562efe9872eSLisandro Dalcin 
563efe9872eSLisandro Dalcin   PetscFunctionBegin;
564efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5659566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
566efe9872eSLisandro Dalcin   if (fun) *fun = tsdm->ops->i2function;
567800f99ffSJeremy L Thompson   if (ctx) {
568800f99ffSJeremy L Thompson     if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2functionctxcontainer, ctx));
569*2a8381b2SBarry Smith     else *(void **)ctx = NULL;
570800f99ffSJeremy L Thompson   }
5713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
572efe9872eSLisandro Dalcin }
573efe9872eSLisandro Dalcin 
574efe9872eSLisandro Dalcin /*@C
575d1c5d1fcSBarry Smith   DMTSSetI2Jacobian - set `TS` implicit Jacobian evaluation function for 2nd order systems from a `DMTS`
576efe9872eSLisandro Dalcin 
577efe9872eSLisandro Dalcin   Not Collective
578efe9872eSLisandro Dalcin 
5794165533cSJose E. Roman   Input Parameters:
580bcf0153eSBarry Smith + dm  - `DM` to be used with `TS`
58114d0ab18SJacob Faibussowitsch . jac - Jacobian evaluation routine
582efe9872eSLisandro Dalcin - ctx - context for Jacobian evaluation
583efe9872eSLisandro Dalcin 
584d1c5d1fcSBarry Smith   Level: developer
585efe9872eSLisandro Dalcin 
586efe9872eSLisandro Dalcin   Note:
587bcf0153eSBarry Smith   `TSSetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
588bcf0153eSBarry Smith   associated with the `DM`.
589efe9872eSLisandro Dalcin 
5908434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSI2JacobianFn`, `TSSetI2Jacobian()`
591efe9872eSLisandro Dalcin @*/
DMTSSetI2Jacobian(DM dm,TSI2JacobianFn * jac,PetscCtx ctx)592*2a8381b2SBarry Smith PetscErrorCode DMTSSetI2Jacobian(DM dm, TSI2JacobianFn *jac, PetscCtx ctx)
593d71ae5a4SJacob Faibussowitsch {
594efe9872eSLisandro Dalcin   DMTS tsdm;
595efe9872eSLisandro Dalcin 
596efe9872eSLisandro Dalcin   PetscFunctionBegin;
597efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5989566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
599efe9872eSLisandro Dalcin   if (jac) tsdm->ops->i2jacobian = jac;
600800f99ffSJeremy L Thompson   if (ctx) {
601800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
602800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
603800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
604800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2jacobian ctx", (PetscObject)ctxcontainer));
605800f99ffSJeremy L Thompson     tsdm->i2jacobianctxcontainer = ctxcontainer;
606800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
607800f99ffSJeremy L Thompson   }
6083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
609800f99ffSJeremy L Thompson }
610800f99ffSJeremy L Thompson 
611800f99ffSJeremy L Thompson /*@C
612d1c5d1fcSBarry Smith   DMTSSetI2JacobianContextDestroy - set `TS` implicit Jacobian evaluation for 2nd order systems context destroy function into a `DMTS`
613800f99ffSJeremy L Thompson 
614800f99ffSJeremy L Thompson   Not Collective
615800f99ffSJeremy L Thompson 
616800f99ffSJeremy L Thompson   Input Parameters:
6175cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
61849abdd8aSBarry Smith - f  - implicit Jacobian evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
619800f99ffSJeremy L Thompson 
620d1c5d1fcSBarry Smith   Level: developer
62187497f52SBarry Smith 
622d1c5d1fcSBarry Smith   Note:
623d1c5d1fcSBarry Smith   Normally `TSSetI2JacobianContextDestroy()` is used
624d1c5d1fcSBarry Smith 
625d1c5d1fcSBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSSetI2JacobianContextDestroy()`, `DMTSSetI2Jacobian()`, `TSSetI2Jacobian()`
626800f99ffSJeremy L Thompson @*/
DMTSSetI2JacobianContextDestroy(DM dm,PetscCtxDestroyFn * f)62749abdd8aSBarry Smith PetscErrorCode DMTSSetI2JacobianContextDestroy(DM dm, PetscCtxDestroyFn *f)
628d71ae5a4SJacob Faibussowitsch {
629800f99ffSJeremy L Thompson   DMTS tsdm;
630800f99ffSJeremy L Thompson 
631800f99ffSJeremy L Thompson   PetscFunctionBegin;
632800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
633800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
63449abdd8aSBarry Smith   if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->i2jacobianctxcontainer, f));
6353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
636800f99ffSJeremy L Thompson }
637800f99ffSJeremy L Thompson 
DMTSUnsetI2JacobianContext_Internal(DM dm)638d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM dm)
639d71ae5a4SJacob Faibussowitsch {
640800f99ffSJeremy L Thompson   DMTS tsdm;
641800f99ffSJeremy L Thompson 
642800f99ffSJeremy L Thompson   PetscFunctionBegin;
643800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
644800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
645800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2JacobianContext_DMTS(tsdm));
6463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
647efe9872eSLisandro Dalcin }
648efe9872eSLisandro Dalcin 
649efe9872eSLisandro Dalcin /*@C
650d1c5d1fcSBarry Smith   DMTSGetI2Jacobian - get `TS` implicit Jacobian evaluation function for 2nd order systems from a `DMTS`
651efe9872eSLisandro Dalcin 
652efe9872eSLisandro Dalcin   Not Collective
653efe9872eSLisandro Dalcin 
6544165533cSJose E. Roman   Input Parameter:
655bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
656efe9872eSLisandro Dalcin 
6574165533cSJose E. Roman   Output Parameters:
6588434afd1SBarry Smith + jac - Jacobian evaluation function,  for calling sequence see `TSI2JacobianFn`
659efe9872eSLisandro Dalcin - ctx - context for Jacobian evaluation
660efe9872eSLisandro Dalcin 
661d1c5d1fcSBarry Smith   Level: developer
662efe9872eSLisandro Dalcin 
663efe9872eSLisandro Dalcin   Note:
664bcf0153eSBarry Smith   `TSGetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
665bcf0153eSBarry Smith   associated with the `DM`.
666efe9872eSLisandro Dalcin 
6678434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetI2Jacobian()`, `TSGetI2Jacobian()`, `TSI2JacobianFn`
668efe9872eSLisandro Dalcin @*/
DMTSGetI2Jacobian(DM dm,TSI2JacobianFn ** jac,PetscCtxRt ctx)669*2a8381b2SBarry Smith PetscErrorCode DMTSGetI2Jacobian(DM dm, TSI2JacobianFn **jac, PetscCtxRt ctx)
670d71ae5a4SJacob Faibussowitsch {
671efe9872eSLisandro Dalcin   DMTS tsdm;
672efe9872eSLisandro Dalcin 
673efe9872eSLisandro Dalcin   PetscFunctionBegin;
674efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6759566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
676efe9872eSLisandro Dalcin   if (jac) *jac = tsdm->ops->i2jacobian;
677800f99ffSJeremy L Thompson   if (ctx) {
678800f99ffSJeremy L Thompson     if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2jacobianctxcontainer, ctx));
679*2a8381b2SBarry Smith     else *(void **)ctx = NULL;
680800f99ffSJeremy L Thompson   }
6813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
682efe9872eSLisandro Dalcin }
68324989b8cSPeter Brune 
68424989b8cSPeter Brune /*@C
685d1c5d1fcSBarry Smith   DMTSSetRHSFunction - set `TS` explicit residual evaluation function into a `DMTS`
68624989b8cSPeter Brune 
68724989b8cSPeter Brune   Not Collective
68824989b8cSPeter Brune 
6894165533cSJose E. Roman   Input Parameters:
690bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
6918434afd1SBarry Smith . func - RHS function evaluation routine, see `TSRHSFunctionFn` for the calling sequence
69224989b8cSPeter Brune - ctx  - context for residual evaluation
69324989b8cSPeter Brune 
694d1c5d1fcSBarry Smith   Level: developer
69524989b8cSPeter Brune 
69624989b8cSPeter Brune   Note:
697bcf0153eSBarry Smith   `TSSetRHSFunction()` is normally used, but it calls this function internally because the user context is actually
698bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
699bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
70024989b8cSPeter Brune 
7018434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSRHSFunctionFn`
70224989b8cSPeter Brune @*/
DMTSSetRHSFunction(DM dm,TSRHSFunctionFn * func,PetscCtx ctx)703*2a8381b2SBarry Smith PetscErrorCode DMTSSetRHSFunction(DM dm, TSRHSFunctionFn *func, PetscCtx ctx)
704d71ae5a4SJacob Faibussowitsch {
705942e3340SBarry Smith   DMTS tsdm;
70624989b8cSPeter Brune 
70724989b8cSPeter Brune   PetscFunctionBegin;
70824989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7099566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
710d74926cbSBarry Smith   if (func) tsdm->ops->rhsfunction = func;
711800f99ffSJeremy L Thompson   if (ctx) {
712800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
713800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
714800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
715800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs function ctx", (PetscObject)ctxcontainer));
716800f99ffSJeremy L Thompson     tsdm->rhsfunctionctxcontainer = ctxcontainer;
717800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
718800f99ffSJeremy L Thompson   }
7193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
720800f99ffSJeremy L Thompson }
721800f99ffSJeremy L Thompson 
722800f99ffSJeremy L Thompson /*@C
723d1c5d1fcSBarry Smith   DMTSSetRHSFunctionContextDestroy - set `TS` explicit residual evaluation context destroy function into a `DMTS`
724800f99ffSJeremy L Thompson 
725800f99ffSJeremy L Thompson   Not Collective
726800f99ffSJeremy L Thompson 
727800f99ffSJeremy L Thompson   Input Parameters:
7285cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
72949abdd8aSBarry Smith - f  - explicit evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
730800f99ffSJeremy L Thompson 
731d1c5d1fcSBarry Smith   Level: developer
732800f99ffSJeremy L Thompson 
733800f99ffSJeremy L Thompson   Note:
7345cb80ecdSBarry Smith   `TSSetRHSFunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
7355cb80ecdSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
7365cb80ecdSBarry Smith   not.
7375cb80ecdSBarry Smith 
738b43aa488SJacob Faibussowitsch   Developer Notes:
7395cb80ecdSBarry Smith   If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
740800f99ffSJeremy L Thompson 
741d1c5d1fcSBarry Smith .seealso: [](ch_ts), `DMTS`, `TSSetRHSFunctionContextDestroy()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
742800f99ffSJeremy L Thompson @*/
DMTSSetRHSFunctionContextDestroy(DM dm,PetscCtxDestroyFn * f)74349abdd8aSBarry Smith PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM dm, PetscCtxDestroyFn *f)
744d71ae5a4SJacob Faibussowitsch {
745800f99ffSJeremy L Thompson   DMTS tsdm;
746800f99ffSJeremy L Thompson 
747800f99ffSJeremy L Thompson   PetscFunctionBegin;
748800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
749800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
75049abdd8aSBarry Smith   if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->rhsfunctionctxcontainer, f));
7513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
752800f99ffSJeremy L Thompson }
753800f99ffSJeremy L Thompson 
DMTSUnsetRHSFunctionContext_Internal(DM dm)754d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM dm)
755d71ae5a4SJacob Faibussowitsch {
756800f99ffSJeremy L Thompson   DMTS tsdm;
757800f99ffSJeremy L Thompson 
758800f99ffSJeremy L Thompson   PetscFunctionBegin;
759800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
760800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
761800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetRHSFunctionContext_DMTS(tsdm));
762800f99ffSJeremy L Thompson   tsdm->rhsfunctionctxcontainer = NULL;
7633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76424989b8cSPeter Brune }
76524989b8cSPeter Brune 
766ef20d060SBarry Smith /*@C
767d1c5d1fcSBarry Smith   DMTSSetTransientVariable - sets function to transform from state to transient variables into a `DMTS`
768e3c11fc1SJed Brown 
769e3c11fc1SJed Brown   Logically Collective
770e3c11fc1SJed Brown 
7714165533cSJose E. Roman   Input Parameters:
772bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
7738434afd1SBarry Smith . tvar - a function that transforms to transient variables, see `TSTransientVariableFn` for the calling sequence
774e3c11fc1SJed Brown - ctx  - a context for tvar
775e3c11fc1SJed Brown 
776d1c5d1fcSBarry Smith   Level: developer
777e3c11fc1SJed Brown 
778e3c11fc1SJed Brown   Notes:
779d1c5d1fcSBarry Smith   Normally `TSSetTransientVariable()` is used
780d1c5d1fcSBarry Smith 
781bcf0153eSBarry Smith   This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., `TSBDF`)
782e3c11fc1SJed Brown   can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
783e3c11fc1SJed Brown   well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
784e3c11fc1SJed Brown   C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
785e3c11fc1SJed Brown   evaluated via the chain rule, as in
786e3c11fc1SJed Brown 
787d1c5d1fcSBarry Smith   $$
788e3c11fc1SJed Brown   dF/dP + shift * dF/dCdot dC/dP.
789d1c5d1fcSBarry Smith   $$
790e3c11fc1SJed Brown 
7918434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `TS`, `TSBDF`, `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()`, `TSTransientVariableFn`
792e3c11fc1SJed Brown @*/
DMTSSetTransientVariable(DM dm,TSTransientVariableFn * tvar,PetscCtx ctx)793*2a8381b2SBarry Smith PetscErrorCode DMTSSetTransientVariable(DM dm, TSTransientVariableFn *tvar, PetscCtx ctx)
794d71ae5a4SJacob Faibussowitsch {
795e3c11fc1SJed Brown   DMTS dmts;
796e3c11fc1SJed Brown 
797e3c11fc1SJed Brown   PetscFunctionBegin;
798e3c11fc1SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7999566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &dmts));
800e3c11fc1SJed Brown   dmts->ops->transientvar = tvar;
801e3c11fc1SJed Brown   dmts->transientvarctx   = ctx;
8023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
803e3c11fc1SJed Brown }
804e3c11fc1SJed Brown 
805e3c11fc1SJed Brown /*@C
806d1c5d1fcSBarry Smith   DMTSGetTransientVariable - gets function to transform from state to transient variables set with `DMTSSetTransientVariable()` from a `TSDM`
807e3c11fc1SJed Brown 
808e3c11fc1SJed Brown   Logically Collective
809e3c11fc1SJed Brown 
8104165533cSJose E. Roman   Input Parameter:
811bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
812e3c11fc1SJed Brown 
8134165533cSJose E. Roman   Output Parameters:
8148434afd1SBarry Smith + tvar - a function that transforms to transient variables, see `TSTransientVariableFn` for the calling sequence
815e3c11fc1SJed Brown - ctx  - a context for tvar
816e3c11fc1SJed Brown 
817d1c5d1fcSBarry Smith   Level: developer
818e3c11fc1SJed Brown 
819d1c5d1fcSBarry Smith   Note:
820d1c5d1fcSBarry Smith   Normally `TSSetTransientVariable()` is used
821d1c5d1fcSBarry Smith 
8228434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()`, `TSTransientVariableFn`
823e3c11fc1SJed Brown @*/
DMTSGetTransientVariable(DM dm,TSTransientVariableFn ** tvar,PetscCtx ctx)824*2a8381b2SBarry Smith PetscErrorCode DMTSGetTransientVariable(DM dm, TSTransientVariableFn **tvar, PetscCtx ctx)
825d71ae5a4SJacob Faibussowitsch {
826e3c11fc1SJed Brown   DMTS dmts;
827e3c11fc1SJed Brown 
828e3c11fc1SJed Brown   PetscFunctionBegin;
829e3c11fc1SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8309566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &dmts));
831e3c11fc1SJed Brown   if (tvar) *tvar = dmts->ops->transientvar;
832e3c11fc1SJed Brown   if (ctx) *(void **)ctx = dmts->transientvarctx;
8333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
834e3c11fc1SJed Brown }
835e3c11fc1SJed Brown 
836e3c11fc1SJed Brown /*@C
837d1c5d1fcSBarry Smith   DMTSGetSolutionFunction - gets the `TS` solution evaluation function from a `DMTS`
838ef20d060SBarry Smith 
839ef20d060SBarry Smith   Not Collective
840ef20d060SBarry Smith 
8414165533cSJose E. Roman   Input Parameter:
842bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
843ef20d060SBarry Smith 
844ef20d060SBarry Smith   Output Parameters:
8458434afd1SBarry Smith + func - solution function evaluation function, for calling sequence see `TSSolutionFn`
846ef20d060SBarry Smith - ctx  - context for solution evaluation
847ef20d060SBarry Smith 
848d1c5d1fcSBarry Smith   Level: developer
849ef20d060SBarry Smith 
8508434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `DMTSSetSolutionFunction()`, `TSSolutionFn`
851ef20d060SBarry Smith @*/
DMTSGetSolutionFunction(DM dm,TSSolutionFn ** func,PetscCtxRt ctx)852*2a8381b2SBarry Smith PetscErrorCode DMTSGetSolutionFunction(DM dm, TSSolutionFn **func, PetscCtxRt ctx)
853d71ae5a4SJacob Faibussowitsch {
854942e3340SBarry Smith   DMTS tsdm;
855ef20d060SBarry Smith 
856ef20d060SBarry Smith   PetscFunctionBegin;
857ef20d060SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8589566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
859d74926cbSBarry Smith   if (func) *func = tsdm->ops->solution;
860*2a8381b2SBarry Smith   if (ctx) *(void **)ctx = tsdm->solutionctx;
8613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
862ef20d060SBarry Smith }
863ef20d060SBarry Smith 
864ef20d060SBarry Smith /*@C
865d1c5d1fcSBarry Smith   DMTSSetSolutionFunction - set `TS` solution evaluation function into a `DMTS`
866ef20d060SBarry Smith 
867ef20d060SBarry Smith   Not Collective
868ef20d060SBarry Smith 
8694165533cSJose E. Roman   Input Parameters:
870bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
8718434afd1SBarry Smith . func - solution function evaluation routine, for calling sequence see `TSSolutionFn`
872ef20d060SBarry Smith - ctx  - context for solution evaluation
873ef20d060SBarry Smith 
874d1c5d1fcSBarry Smith   Level: developer
875ef20d060SBarry Smith 
876ef20d060SBarry Smith   Note:
877bcf0153eSBarry Smith   `TSSetSolutionFunction()` is normally used, but it calls this function internally because the user context is actually
878bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
879bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
880ef20d060SBarry Smith 
8818434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSGetSolutionFunction()`, `TSSolutionFn`
882ef20d060SBarry Smith @*/
DMTSSetSolutionFunction(DM dm,TSSolutionFn * func,PetscCtx ctx)883*2a8381b2SBarry Smith PetscErrorCode DMTSSetSolutionFunction(DM dm, TSSolutionFn *func, PetscCtx ctx)
884d71ae5a4SJacob Faibussowitsch {
885942e3340SBarry Smith   DMTS tsdm;
886ef20d060SBarry Smith 
887ef20d060SBarry Smith   PetscFunctionBegin;
888ef20d060SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8899566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
890d74926cbSBarry Smith   if (func) tsdm->ops->solution = func;
891ef20d060SBarry Smith   if (ctx) tsdm->solutionctx = ctx;
8923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
893ef20d060SBarry Smith }
894ef20d060SBarry Smith 
8959b7cd975SBarry Smith /*@C
896d1c5d1fcSBarry Smith   DMTSSetForcingFunction - set `TS` forcing function evaluation function into a `DMTS`
8979b7cd975SBarry Smith 
8989b7cd975SBarry Smith   Not Collective
8999b7cd975SBarry Smith 
9004165533cSJose E. Roman   Input Parameters:
901bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
9028434afd1SBarry Smith . func - forcing function evaluation routine, for calling sequence see `TSForcingFn`
9039b7cd975SBarry Smith - ctx  - context for solution evaluation
9049b7cd975SBarry Smith 
905d1c5d1fcSBarry Smith   Level: developer
9069b7cd975SBarry Smith 
9079b7cd975SBarry Smith   Note:
908bcf0153eSBarry Smith   `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually
909bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
910bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
9119b7cd975SBarry Smith 
9128434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSForcingFn`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
9139b7cd975SBarry Smith @*/
DMTSSetForcingFunction(DM dm,TSForcingFn * func,PetscCtx ctx)914*2a8381b2SBarry Smith PetscErrorCode DMTSSetForcingFunction(DM dm, TSForcingFn *func, PetscCtx ctx)
915d71ae5a4SJacob Faibussowitsch {
9169b7cd975SBarry Smith   DMTS tsdm;
9179b7cd975SBarry Smith 
9189b7cd975SBarry Smith   PetscFunctionBegin;
9199b7cd975SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9209566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
9212fe279fdSBarry Smith   if (func) tsdm->ops->forcing = func;
9229b7cd975SBarry Smith   if (ctx) tsdm->forcingctx = ctx;
9233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9249b7cd975SBarry Smith }
9259b7cd975SBarry Smith 
9269b7cd975SBarry Smith /*@C
927d1c5d1fcSBarry Smith   DMTSGetForcingFunction - get `TS` forcing function evaluation function from a `DMTS`
9289b7cd975SBarry Smith 
9299b7cd975SBarry Smith   Not Collective
9309b7cd975SBarry Smith 
9314165533cSJose E. Roman   Input Parameter:
932bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
9339b7cd975SBarry Smith 
9344165533cSJose E. Roman   Output Parameters:
9358434afd1SBarry Smith + f   - forcing function evaluation function; see `TSForcingFn` for the calling sequence
9369b7cd975SBarry Smith - ctx - context for solution evaluation
9379b7cd975SBarry Smith 
938d1c5d1fcSBarry Smith   Level: developer
9399b7cd975SBarry Smith 
9409b7cd975SBarry Smith   Note:
941bcf0153eSBarry Smith   `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually
942bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
943bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
9449b7cd975SBarry Smith 
9458434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `TSSetForcingFunction()`, `TSForcingFn`
9469b7cd975SBarry Smith @*/
DMTSGetForcingFunction(DM dm,TSForcingFn ** f,PetscCtxRt ctx)947*2a8381b2SBarry Smith PetscErrorCode DMTSGetForcingFunction(DM dm, TSForcingFn **f, PetscCtxRt ctx)
948d71ae5a4SJacob Faibussowitsch {
9499b7cd975SBarry Smith   DMTS tsdm;
9509b7cd975SBarry Smith 
9519b7cd975SBarry Smith   PetscFunctionBegin;
9529b7cd975SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9539566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
954f8b49ee9SBarry Smith   if (f) *f = tsdm->ops->forcing;
955*2a8381b2SBarry Smith   if (ctx) *(void **)ctx = tsdm->forcingctx;
9563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9579b7cd975SBarry Smith }
9589b7cd975SBarry Smith 
95924989b8cSPeter Brune /*@C
960d1c5d1fcSBarry Smith   DMTSGetRHSFunction - get `TS` explicit residual evaluation function from a `DMTS`
96124989b8cSPeter Brune 
96224989b8cSPeter Brune   Not Collective
96324989b8cSPeter Brune 
9644165533cSJose E. Roman   Input Parameter:
965bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
96624989b8cSPeter Brune 
9674165533cSJose E. Roman   Output Parameters:
9688434afd1SBarry Smith + func - residual evaluation function, for calling sequence see `TSRHSFunctionFn`
96924989b8cSPeter Brune - ctx  - context for residual evaluation
97024989b8cSPeter Brune 
971d1c5d1fcSBarry Smith   Level: developer
97224989b8cSPeter Brune 
97324989b8cSPeter Brune   Note:
974346ce620SStefano Zampini   `TSGetRHSFunction()` is normally used, but it calls this function internally because the user context is actually
97524989b8cSPeter Brune   associated with the DM.
97624989b8cSPeter Brune 
9778434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSRHSFunctionFn`, `TSGetRHSFunction()`
97824989b8cSPeter Brune @*/
DMTSGetRHSFunction(DM dm,TSRHSFunctionFn ** func,PetscCtxRt ctx)979*2a8381b2SBarry Smith PetscErrorCode DMTSGetRHSFunction(DM dm, TSRHSFunctionFn **func, PetscCtxRt ctx)
980d71ae5a4SJacob Faibussowitsch {
981942e3340SBarry Smith   DMTS tsdm;
98224989b8cSPeter Brune 
98324989b8cSPeter Brune   PetscFunctionBegin;
98424989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9859566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
986d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsfunction;
987800f99ffSJeremy L Thompson   if (ctx) {
988800f99ffSJeremy L Thompson     if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsfunctionctxcontainer, ctx));
989*2a8381b2SBarry Smith     else *(void **)ctx = NULL;
990800f99ffSJeremy L Thompson   }
9913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99224989b8cSPeter Brune }
99324989b8cSPeter Brune 
99424989b8cSPeter Brune /*@C
995d1c5d1fcSBarry Smith   DMTSSetIJacobian - set `TS` Jacobian evaluation function into a `DMTS`
99624989b8cSPeter Brune 
99724989b8cSPeter Brune   Not Collective
99824989b8cSPeter Brune 
9994165533cSJose E. Roman   Input Parameters:
1000bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
10018434afd1SBarry Smith . func - Jacobian evaluation routine, see `TSIJacobianFn` for the calling sequence
100224989b8cSPeter Brune - ctx  - context for residual evaluation
100324989b8cSPeter Brune 
1004d1c5d1fcSBarry Smith   Level: developer
100524989b8cSPeter Brune 
100624989b8cSPeter Brune   Note:
1007346ce620SStefano Zampini   `TSSetIJacobian()` is normally used, but it calls this function internally because the user context is actually
1008bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1009bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
101024989b8cSPeter Brune 
10118434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `TSIJacobianFn`, `DMTSGetIJacobian()`, `TSSetIJacobian()`
101224989b8cSPeter Brune @*/
DMTSSetIJacobian(DM dm,TSIJacobianFn * func,PetscCtx ctx)1013*2a8381b2SBarry Smith PetscErrorCode DMTSSetIJacobian(DM dm, TSIJacobianFn *func, PetscCtx ctx)
1014d71ae5a4SJacob Faibussowitsch {
1015800f99ffSJeremy L Thompson   DMTS tsdm;
101624989b8cSPeter Brune 
101724989b8cSPeter Brune   PetscFunctionBegin;
101824989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1019800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1020800f99ffSJeremy L Thompson   if (func) tsdm->ops->ijacobian = func;
1021800f99ffSJeremy L Thompson   if (ctx) {
1022800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
1023800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
1024800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
1025800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "ijacobian ctx", (PetscObject)ctxcontainer));
1026800f99ffSJeremy L Thompson     tsdm->ijacobianctxcontainer = ctxcontainer;
1027800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
1028800f99ffSJeremy L Thompson   }
10293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1030800f99ffSJeremy L Thompson }
1031800f99ffSJeremy L Thompson 
1032800f99ffSJeremy L Thompson /*@C
1033d1c5d1fcSBarry Smith   DMTSSetIJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function into a `DMTS`
1034800f99ffSJeremy L Thompson 
1035800f99ffSJeremy L Thompson   Not Collective
1036800f99ffSJeremy L Thompson 
1037800f99ffSJeremy L Thompson   Input Parameters:
10385cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
103949abdd8aSBarry Smith - f  - Jacobian evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
1040800f99ffSJeremy L Thompson 
1041d1c5d1fcSBarry Smith   Level: developer
1042800f99ffSJeremy L Thompson 
1043800f99ffSJeremy L Thompson   Note:
10445cb80ecdSBarry Smith   `TSSetIJacobianContextDestroy()` is normally used, but it calls this function internally because the user context is actually
10455cb80ecdSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
10465cb80ecdSBarry Smith   not.
1047800f99ffSJeremy L Thompson 
1048b43aa488SJacob Faibussowitsch   Developer Notes:
10495cb80ecdSBarry Smith   If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
10505cb80ecdSBarry Smith 
1051d1c5d1fcSBarry Smith .seealso: [](ch_ts), `DMTS`, `TSSetIJacobianContextDestroy()`, `TSSetI2JacobianContextDestroy()`, `DMTSSetIJacobian()`, `TSSetIJacobian()`
1052800f99ffSJeremy L Thompson @*/
DMTSSetIJacobianContextDestroy(DM dm,PetscCtxDestroyFn * f)105349abdd8aSBarry Smith PetscErrorCode DMTSSetIJacobianContextDestroy(DM dm, PetscCtxDestroyFn *f)
1054d71ae5a4SJacob Faibussowitsch {
1055800f99ffSJeremy L Thompson   DMTS tsdm;
1056800f99ffSJeremy L Thompson 
1057800f99ffSJeremy L Thompson   PetscFunctionBegin;
1058800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1059800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
106049abdd8aSBarry Smith   if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->ijacobianctxcontainer, f));
10613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1062800f99ffSJeremy L Thompson }
1063800f99ffSJeremy L Thompson 
DMTSUnsetIJacobianContext_Internal(DM dm)1064d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM dm)
1065d71ae5a4SJacob Faibussowitsch {
1066800f99ffSJeremy L Thompson   DMTS tsdm;
1067800f99ffSJeremy L Thompson 
1068800f99ffSJeremy L Thompson   PetscFunctionBegin;
1069800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1070800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1071800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIJacobianContext_DMTS(tsdm));
10723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
107324989b8cSPeter Brune }
107424989b8cSPeter Brune 
107524989b8cSPeter Brune /*@C
1076d1c5d1fcSBarry Smith   DMTSGetIJacobian - get `TS` Jacobian evaluation function from a `DMTS`
107724989b8cSPeter Brune 
107824989b8cSPeter Brune   Not Collective
107924989b8cSPeter Brune 
10804165533cSJose E. Roman   Input Parameter:
1081bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
108224989b8cSPeter Brune 
10834165533cSJose E. Roman   Output Parameters:
10848434afd1SBarry Smith + func - Jacobian evaluation function, for calling sequence see `TSIJacobianFn`
108524989b8cSPeter Brune - ctx  - context for residual evaluation
108624989b8cSPeter Brune 
1087d1c5d1fcSBarry Smith   Level: developer
108824989b8cSPeter Brune 
108924989b8cSPeter Brune   Note:
1090346ce620SStefano Zampini   `TSGetIJacobian()` is normally used, but it calls this function internally because the user context is actually
1091bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1092bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
109324989b8cSPeter Brune 
10948434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetIJacobian()`, `TSIJacobianFn`
109524989b8cSPeter Brune @*/
DMTSGetIJacobian(DM dm,TSIJacobianFn ** func,PetscCtxRt ctx)1096*2a8381b2SBarry Smith PetscErrorCode DMTSGetIJacobian(DM dm, TSIJacobianFn **func, PetscCtxRt ctx)
1097d71ae5a4SJacob Faibussowitsch {
1098942e3340SBarry Smith   DMTS tsdm;
109924989b8cSPeter Brune 
110024989b8cSPeter Brune   PetscFunctionBegin;
110124989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11029566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
1103d74926cbSBarry Smith   if (func) *func = tsdm->ops->ijacobian;
1104800f99ffSJeremy L Thompson   if (ctx) {
1105800f99ffSJeremy L Thompson     if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ijacobianctxcontainer, ctx));
1106*2a8381b2SBarry Smith     else *(void **)ctx = NULL;
1107800f99ffSJeremy L Thompson   }
11083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110924989b8cSPeter Brune }
111024989b8cSPeter Brune 
111124989b8cSPeter Brune /*@C
1112d1c5d1fcSBarry Smith   DMTSSetRHSJacobian - set `TS` Jacobian evaluation function into a `DMTS`
111324989b8cSPeter Brune 
111424989b8cSPeter Brune   Not Collective
111524989b8cSPeter Brune 
11164165533cSJose E. Roman   Input Parameters:
1117bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
11188434afd1SBarry Smith . func - Jacobian evaluation routine, for calling sequence see `TSIJacobianFn`
111924989b8cSPeter Brune - ctx  - context for residual evaluation
112024989b8cSPeter Brune 
1121d1c5d1fcSBarry Smith   Level: developer
112224989b8cSPeter Brune 
112324989b8cSPeter Brune   Note:
1124346ce620SStefano Zampini   `TSSetRHSJacobian()` is normally used, but it calls this function internally because the user context is actually
1125bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1126bcf0153eSBarry Smith   not.
112724989b8cSPeter Brune 
1128b43aa488SJacob Faibussowitsch   Developer Notes:
1129bcf0153eSBarry Smith   If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1130bcf0153eSBarry Smith 
11318434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `TSRHSJacobianFn`, `DMTSGetRHSJacobian()`, `TSSetRHSJacobian()`
113224989b8cSPeter Brune @*/
DMTSSetRHSJacobian(DM dm,TSRHSJacobianFn * func,PetscCtx ctx)1133*2a8381b2SBarry Smith PetscErrorCode DMTSSetRHSJacobian(DM dm, TSRHSJacobianFn *func, PetscCtx ctx)
1134d71ae5a4SJacob Faibussowitsch {
1135942e3340SBarry Smith   DMTS tsdm;
113624989b8cSPeter Brune 
113724989b8cSPeter Brune   PetscFunctionBegin;
113824989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11399566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1140d74926cbSBarry Smith   if (func) tsdm->ops->rhsjacobian = func;
1141800f99ffSJeremy L Thompson   if (ctx) {
1142800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
1143800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
1144800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
1145800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs jacobian ctx", (PetscObject)ctxcontainer));
1146800f99ffSJeremy L Thompson     tsdm->rhsjacobianctxcontainer = ctxcontainer;
1147800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
1148800f99ffSJeremy L Thompson   }
11493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1150800f99ffSJeremy L Thompson }
1151800f99ffSJeremy L Thompson 
1152800f99ffSJeremy L Thompson /*@C
1153d1c5d1fcSBarry Smith   DMTSSetRHSJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function from a `DMTS`
1154800f99ffSJeremy L Thompson 
1155800f99ffSJeremy L Thompson   Not Collective
1156800f99ffSJeremy L Thompson 
1157800f99ffSJeremy L Thompson   Input Parameters:
11585cb80ecdSBarry Smith + dm - `DM` to be used with `TS`
115949abdd8aSBarry Smith - f  - Jacobian evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
11605cb80ecdSBarry Smith 
1161d1c5d1fcSBarry Smith   Level: developer
11625cb80ecdSBarry Smith 
11635cb80ecdSBarry Smith   Note:
11645cb80ecdSBarry Smith   The user usually calls `TSSetRHSJacobianContextDestroy()` which calls this routine
1165800f99ffSJeremy L Thompson 
1166d1c5d1fcSBarry Smith .seealso: [](ch_ts), `DMTS`, `TS`, `TSSetRHSJacobianContextDestroy()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
1167800f99ffSJeremy L Thompson @*/
DMTSSetRHSJacobianContextDestroy(DM dm,PetscCtxDestroyFn * f)116849abdd8aSBarry Smith PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM dm, PetscCtxDestroyFn *f)
1169d71ae5a4SJacob Faibussowitsch {
1170800f99ffSJeremy L Thompson   DMTS tsdm;
1171800f99ffSJeremy L Thompson 
1172800f99ffSJeremy L Thompson   PetscFunctionBegin;
1173800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1174800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
117549abdd8aSBarry Smith   if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->rhsjacobianctxcontainer, f));
11763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1177800f99ffSJeremy L Thompson }
1178800f99ffSJeremy L Thompson 
DMTSUnsetRHSJacobianContext_Internal(DM dm)1179d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM dm)
1180d71ae5a4SJacob Faibussowitsch {
1181800f99ffSJeremy L Thompson   DMTS tsdm;
1182800f99ffSJeremy L Thompson 
1183800f99ffSJeremy L Thompson   PetscFunctionBegin;
1184800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1185800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1186800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetRHSJacobianContext_DMTS(tsdm));
11873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118824989b8cSPeter Brune }
118924989b8cSPeter Brune 
119024989b8cSPeter Brune /*@C
1191d1c5d1fcSBarry Smith   DMTSGetRHSJacobian - get `TS` Jacobian evaluation function from a `DMTS`
119224989b8cSPeter Brune 
119324989b8cSPeter Brune   Not Collective
119424989b8cSPeter Brune 
11954165533cSJose E. Roman   Input Parameter:
1196bcf0153eSBarry Smith . dm - `DM` to be used with `TS`
119724989b8cSPeter Brune 
11984165533cSJose E. Roman   Output Parameters:
11998434afd1SBarry Smith + func - Jacobian evaluation function, for calling sequence see `TSRHSJacobianFn`
120024989b8cSPeter Brune - ctx  - context for residual evaluation
120124989b8cSPeter Brune 
1202d1c5d1fcSBarry Smith   Level: developer
120324989b8cSPeter Brune 
120424989b8cSPeter Brune   Note:
1205346ce620SStefano Zampini   `TSGetRHSJacobian()` is normally used, but it calls this function internally because the user context is actually
1206bcf0153eSBarry Smith   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
1207bcf0153eSBarry Smith   not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
120824989b8cSPeter Brune 
12098434afd1SBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetRHSJacobian()`, `TSRHSJacobianFn`
121024989b8cSPeter Brune @*/
DMTSGetRHSJacobian(DM dm,TSRHSJacobianFn ** func,PetscCtxRt ctx)1211*2a8381b2SBarry Smith PetscErrorCode DMTSGetRHSJacobian(DM dm, TSRHSJacobianFn **func, PetscCtxRt ctx)
1212d71ae5a4SJacob Faibussowitsch {
1213942e3340SBarry Smith   DMTS tsdm;
121424989b8cSPeter Brune 
121524989b8cSPeter Brune   PetscFunctionBegin;
121624989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12179566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
1218d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsjacobian;
1219800f99ffSJeremy L Thompson   if (ctx) {
1220800f99ffSJeremy L Thompson     if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsjacobianctxcontainer, ctx));
1221*2a8381b2SBarry Smith     else *(void **)ctx = NULL;
1222800f99ffSJeremy L Thompson   }
12233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122424989b8cSPeter Brune }
1225ad6bc421SBarry Smith 
1226ad6bc421SBarry Smith /*@C
12278434afd1SBarry Smith   DMTSSetIFunctionSerialize - sets functions used to view and load a `TSIFunctionFn` context
1228ad6bc421SBarry Smith 
1229ad6bc421SBarry Smith   Not Collective
1230ad6bc421SBarry Smith 
12314165533cSJose E. Roman   Input Parameters:
1232bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
1233ad6bc421SBarry Smith . view - viewer function
1234ad6bc421SBarry Smith - load - loading function
1235ad6bc421SBarry Smith 
1236d1c5d1fcSBarry Smith   Level: developer
1237ad6bc421SBarry Smith 
1238d1c5d1fcSBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`
1239ad6bc421SBarry Smith @*/
DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (* view)(void *,PetscViewer),PetscErrorCode (* load)(void **,PetscViewer))1240d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIFunctionSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1241d71ae5a4SJacob Faibussowitsch {
1242ad6bc421SBarry Smith   DMTS tsdm;
1243ad6bc421SBarry Smith 
1244ad6bc421SBarry Smith   PetscFunctionBegin;
1245ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12469566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1247ad6bc421SBarry Smith   tsdm->ops->ifunctionview = view;
1248ad6bc421SBarry Smith   tsdm->ops->ifunctionload = load;
12493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1250ad6bc421SBarry Smith }
1251ad6bc421SBarry Smith 
1252ad6bc421SBarry Smith /*@C
12538434afd1SBarry Smith   DMTSSetIJacobianSerialize - sets functions used to view and load a `TSIJacobianFn` context
1254ad6bc421SBarry Smith 
1255ad6bc421SBarry Smith   Not Collective
1256ad6bc421SBarry Smith 
12574165533cSJose E. Roman   Input Parameters:
1258bcf0153eSBarry Smith + dm   - `DM` to be used with `TS`
1259ad6bc421SBarry Smith . view - viewer function
1260ad6bc421SBarry Smith - load - loading function
1261ad6bc421SBarry Smith 
1262d1c5d1fcSBarry Smith   Level: developer
1263ad6bc421SBarry Smith 
1264d1c5d1fcSBarry Smith .seealso: [](ch_ts), `DMTS`, `DM`, `TS`
1265ad6bc421SBarry Smith @*/
DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (* view)(void *,PetscViewer),PetscErrorCode (* load)(void **,PetscViewer))1266d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIJacobianSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1267d71ae5a4SJacob Faibussowitsch {
1268ad6bc421SBarry Smith   DMTS tsdm;
1269ad6bc421SBarry Smith 
1270ad6bc421SBarry Smith   PetscFunctionBegin;
1271ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12729566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tsdm));
1273ad6bc421SBarry Smith   tsdm->ops->ijacobianview = view;
1274ad6bc421SBarry Smith   tsdm->ops->ijacobianload = load;
12753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1276ad6bc421SBarry Smith }
1277