xref: /petsc/src/ts/utils/dmts.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
2af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
324989b8cSPeter Brune 
4800f99ffSJeremy L Thompson static PetscErrorCode DMTSUnsetRHSFunctionContext_DMTS(DMTS tsdm)
5800f99ffSJeremy L Thompson {
6800f99ffSJeremy L Thompson   PetscFunctionBegin;
7800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm,"rhs function ctx",NULL));
8800f99ffSJeremy L Thompson   tsdm->rhsfunctionctxcontainer = NULL;
9800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
10800f99ffSJeremy L Thompson }
11800f99ffSJeremy L Thompson 
12800f99ffSJeremy L Thompson static PetscErrorCode DMTSUnsetRHSJacobianContext_DMTS(DMTS tsdm)
13800f99ffSJeremy L Thompson {
14800f99ffSJeremy L Thompson   PetscFunctionBegin;
15800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm,"rhs jacobian ctx",NULL));
16800f99ffSJeremy L Thompson   tsdm->rhsjacobianctxcontainer = NULL;
17800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
18800f99ffSJeremy L Thompson }
19800f99ffSJeremy L Thompson 
20800f99ffSJeremy L Thompson static PetscErrorCode DMTSUnsetIFunctionContext_DMTS(DMTS tsdm)
21800f99ffSJeremy L Thompson {
22800f99ffSJeremy L Thompson   PetscFunctionBegin;
23800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm,"ifunction ctx",NULL));
24800f99ffSJeremy L Thompson   tsdm->ifunctionctxcontainer = NULL;
25800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
26800f99ffSJeremy L Thompson }
27800f99ffSJeremy L Thompson 
28800f99ffSJeremy L Thompson static PetscErrorCode DMTSUnsetIJacobianContext_DMTS(DMTS tsdm)
29800f99ffSJeremy L Thompson {
30800f99ffSJeremy L Thompson   PetscFunctionBegin;
31800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm,"ijacobian ctx",NULL));
32800f99ffSJeremy L Thompson   tsdm->ijacobianctxcontainer = NULL;
33800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
34800f99ffSJeremy L Thompson }
35800f99ffSJeremy L Thompson 
36800f99ffSJeremy L Thompson static PetscErrorCode DMTSUnsetI2FunctionContext_DMTS(DMTS tsdm)
37800f99ffSJeremy L Thompson {
38800f99ffSJeremy L Thompson   PetscFunctionBegin;
39800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm,"i2function ctx",NULL));
40800f99ffSJeremy L Thompson   tsdm->i2functionctxcontainer = NULL;
41800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
42800f99ffSJeremy L Thompson }
43800f99ffSJeremy L Thompson 
44800f99ffSJeremy L Thompson static PetscErrorCode DMTSUnsetI2JacobianContext_DMTS(DMTS tsdm)
45800f99ffSJeremy L Thompson {
46800f99ffSJeremy L Thompson   PetscFunctionBegin;
47800f99ffSJeremy L Thompson   PetscCall(PetscObjectCompose((PetscObject)tsdm,"i2jacobian ctx",NULL));
48800f99ffSJeremy L Thompson   tsdm->i2jacobianctxcontainer = NULL;
49800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
50800f99ffSJeremy L Thompson }
51800f99ffSJeremy L Thompson 
52d74926cbSBarry Smith static PetscErrorCode DMTSDestroy(DMTS *kdm)
53d74926cbSBarry Smith {
54d74926cbSBarry Smith   PetscFunctionBegin;
55d74926cbSBarry Smith   if (!*kdm) PetscFunctionReturn(0);
56d74926cbSBarry Smith   PetscValidHeaderSpecific((*kdm),DMTS_CLASSID,1);
57c793f718SLisandro Dalcin   if (--((PetscObject)(*kdm))->refct > 0) {*kdm = NULL; PetscFunctionReturn(0);}
58800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetRHSFunctionContext_DMTS(*kdm));
59800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetRHSJacobianContext_DMTS(*kdm));
60800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIFunctionContext_DMTS(*kdm));
61800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIJacobianContext_DMTS(*kdm));
62800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2FunctionContext_DMTS(*kdm));
63800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2JacobianContext_DMTS(*kdm));
64*dbbe0bcdSBarry Smith   PetscTryTypeMethod(*kdm,destroy);
659566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(kdm));
66d74926cbSBarry Smith   PetscFunctionReturn(0);
67d74926cbSBarry Smith }
68d74926cbSBarry Smith 
692d53ad75SBarry Smith PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer)
702d53ad75SBarry Smith {
712d53ad75SBarry Smith   PetscFunctionBegin;
729566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,NULL,PETSC_FUNCTION));
739566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,NULL,PETSC_FUNCTION));
749566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,NULL,PETSC_FUNCTION));
75ad6bc421SBarry Smith   if (kdm->ops->ifunctionload) {
76800f99ffSJeremy L Thompson     void *ctx;
77800f99ffSJeremy L Thompson 
78800f99ffSJeremy L Thompson     PetscCall(PetscContainerGetPointer(kdm->ifunctionctxcontainer,&ctx));
79800f99ffSJeremy L Thompson     PetscCall((*kdm->ops->ifunctionload)(&ctx,viewer));
80ad6bc421SBarry Smith   }
819566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,NULL,PETSC_FUNCTION));
829566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,NULL,PETSC_FUNCTION));
839566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,NULL,PETSC_FUNCTION));
84ad6bc421SBarry Smith   if (kdm->ops->ijacobianload) {
85800f99ffSJeremy L Thompson     void *ctx;
86800f99ffSJeremy L Thompson 
87800f99ffSJeremy L Thompson     PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer,&ctx));
88800f99ffSJeremy L Thompson     PetscCall((*kdm->ops->ijacobianload)(&ctx,viewer));
89ad6bc421SBarry Smith   }
902d53ad75SBarry Smith   PetscFunctionReturn(0);
912d53ad75SBarry Smith }
922d53ad75SBarry Smith 
932d53ad75SBarry Smith PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer)
942d53ad75SBarry Smith {
952d53ad75SBarry Smith   PetscBool      isascii,isbinary;
962d53ad75SBarry Smith 
972d53ad75SBarry Smith   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
999566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
1002d53ad75SBarry Smith   if (isascii) {
101c7a10e08SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
1022d53ad75SBarry Smith     const char *fname;
1032d53ad75SBarry Smith 
1049566063dSJacob Faibussowitsch     PetscCall(PetscFPTFind(kdm->ops->ifunction,&fname));
1052d53ad75SBarry Smith     if (fname) {
1069566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  IFunction used by TS: %s\n",fname));
1072d53ad75SBarry Smith     }
1089566063dSJacob Faibussowitsch     PetscCall(PetscFPTFind(kdm->ops->ijacobian,&fname));
1092d53ad75SBarry Smith     if (fname) {
1109566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  IJacobian function used by TS: %s\n",fname));
1112d53ad75SBarry Smith     }
112c7a10e08SBarry Smith #endif
1132d53ad75SBarry Smith   } else if (isbinary) {
1143964eb88SJed Brown     struct {
1153964eb88SJed Brown       TSIFunction ifunction;
1169200755eSBarry Smith     } funcstruct;
1179200755eSBarry Smith     struct {
1183964eb88SJed Brown       PetscErrorCode (*ifunctionview)(void*,PetscViewer);
1199200755eSBarry Smith     } funcviewstruct;
1209200755eSBarry Smith     struct {
1213964eb88SJed Brown       PetscErrorCode (*ifunctionload)(void**,PetscViewer);
1229200755eSBarry Smith     } funcloadstruct;
1233964eb88SJed Brown     struct {
1243964eb88SJed Brown       TSIJacobian ijacobian;
1259200755eSBarry Smith     } jacstruct;
1269200755eSBarry Smith     struct {
1273964eb88SJed Brown       PetscErrorCode (*ijacobianview)(void*,PetscViewer);
1289200755eSBarry Smith     } jacviewstruct;
1299200755eSBarry Smith     struct {
1303964eb88SJed Brown       PetscErrorCode (*ijacobianload)(void**,PetscViewer);
1319200755eSBarry Smith     } jacloadstruct;
1323964eb88SJed Brown 
1339200755eSBarry Smith     funcstruct.ifunction         = kdm->ops->ifunction;
1349200755eSBarry Smith     funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
1359200755eSBarry Smith     funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
1369566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION));
1379566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION));
1389566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION));
139800f99ffSJeremy L Thompson     if (kdm->ops->ifunctionview) {
140800f99ffSJeremy L Thompson       void *ctx;
141800f99ffSJeremy L Thompson 
142800f99ffSJeremy L Thompson       PetscCall(PetscContainerGetPointer(kdm->ifunctionctxcontainer,&ctx));
143800f99ffSJeremy L Thompson       PetscCall((*kdm->ops->ifunctionview)(ctx,viewer));
144800f99ffSJeremy L Thompson     }
1459200755eSBarry Smith     jacstruct.ijacobian = kdm->ops->ijacobian;
1469200755eSBarry Smith     jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
1479200755eSBarry Smith     jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
1489566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION));
1499566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION));
1509566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION));
151800f99ffSJeremy L Thompson     if (kdm->ops->ijacobianview) {
152800f99ffSJeremy L Thompson       void *ctx;
153800f99ffSJeremy L Thompson 
154800f99ffSJeremy L Thompson       PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer,&ctx));
155800f99ffSJeremy L Thompson       PetscCall((*kdm->ops->ijacobianview)(ctx,viewer));
156800f99ffSJeremy L Thompson     }
1572d53ad75SBarry Smith   }
1582d53ad75SBarry Smith   PetscFunctionReturn(0);
1592d53ad75SBarry Smith }
1602d53ad75SBarry Smith 
161d74926cbSBarry Smith static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
162d74926cbSBarry Smith {
163d74926cbSBarry Smith   PetscFunctionBegin;
1649566063dSJacob Faibussowitsch   PetscCall(TSInitializePackage());
1659566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView));
166d74926cbSBarry Smith   PetscFunctionReturn(0);
167d74926cbSBarry Smith }
16824989b8cSPeter Brune 
1692a34c10cSBarry Smith /* Attaches the DMTS to the coarse level.
17024989b8cSPeter Brune  * Under what conditions should we copy versus duplicate?
17124989b8cSPeter Brune  */
172942e3340SBarry Smith static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
17324989b8cSPeter Brune {
17424989b8cSPeter Brune   PetscFunctionBegin;
1759566063dSJacob Faibussowitsch   PetscCall(DMCopyDMTS(dm,dmc));
17624989b8cSPeter Brune   PetscFunctionReturn(0);
17724989b8cSPeter Brune }
17824989b8cSPeter Brune 
17924989b8cSPeter Brune /* This could restrict auxiliary information to the coarse level.
18024989b8cSPeter Brune  */
181942e3340SBarry Smith static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
18224989b8cSPeter Brune {
18324989b8cSPeter Brune   PetscFunctionBegin;
18424989b8cSPeter Brune   PetscFunctionReturn(0);
18524989b8cSPeter Brune }
18624989b8cSPeter Brune 
187258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
188258e1594SPeter Brune {
189258e1594SPeter Brune   PetscFunctionBegin;
1909566063dSJacob Faibussowitsch   PetscCall(DMCopyDMTS(dm,subdm));
191258e1594SPeter Brune   PetscFunctionReturn(0);
192258e1594SPeter Brune }
193258e1594SPeter Brune 
194258e1594SPeter Brune /* This could restrict auxiliary information to the coarse level.
195258e1594SPeter Brune  */
196258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
197258e1594SPeter Brune {
198258e1594SPeter Brune   PetscFunctionBegin;
199258e1594SPeter Brune   PetscFunctionReturn(0);
200258e1594SPeter Brune }
201258e1594SPeter Brune 
202d74926cbSBarry Smith /*@C
203d74926cbSBarry Smith    DMTSCopy - copies the information in a DMTS to another DMTS
204d74926cbSBarry Smith 
205d74926cbSBarry Smith    Not Collective
206d74926cbSBarry Smith 
2074165533cSJose E. Roman    Input Parameters:
208d74926cbSBarry Smith +  kdm - Original DMTS
209d74926cbSBarry Smith -  nkdm - DMTS to receive the data, should have been created with DMTSCreate()
210d74926cbSBarry Smith 
211d74926cbSBarry Smith    Level: developer
212d74926cbSBarry Smith 
213db781477SPatrick Sanan .seealso: `DMTSCreate()`, `DMTSDestroy()`
214d74926cbSBarry Smith @*/
215d74926cbSBarry Smith PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
21624989b8cSPeter Brune {
21724989b8cSPeter Brune   PetscFunctionBegin;
218d74926cbSBarry Smith   PetscValidHeaderSpecific(kdm,DMTS_CLASSID,1);
219d74926cbSBarry Smith   PetscValidHeaderSpecific(nkdm,DMTS_CLASSID,2);
220d74926cbSBarry Smith   nkdm->ops->rhsfunction        = kdm->ops->rhsfunction;
221d74926cbSBarry Smith   nkdm->ops->rhsjacobian        = kdm->ops->rhsjacobian;
222d74926cbSBarry Smith   nkdm->ops->ifunction          = kdm->ops->ifunction;
223d74926cbSBarry Smith   nkdm->ops->ijacobian          = kdm->ops->ijacobian;
224efe9872eSLisandro Dalcin   nkdm->ops->i2function         = kdm->ops->i2function;
225efe9872eSLisandro Dalcin   nkdm->ops->i2jacobian         = kdm->ops->i2jacobian;
226d74926cbSBarry Smith   nkdm->ops->solution           = kdm->ops->solution;
227d74926cbSBarry Smith   nkdm->ops->destroy            = kdm->ops->destroy;
228d74926cbSBarry Smith   nkdm->ops->duplicate          = kdm->ops->duplicate;
229d74926cbSBarry Smith 
230d74926cbSBarry Smith   nkdm->solutionctx             = kdm->solutionctx;
231800f99ffSJeremy L Thompson   nkdm->rhsfunctionctxcontainer = kdm->rhsfunctionctxcontainer;
232800f99ffSJeremy L Thompson   nkdm->rhsjacobianctxcontainer = kdm->rhsjacobianctxcontainer;
233800f99ffSJeremy L Thompson   nkdm->ifunctionctxcontainer   = kdm->ifunctionctxcontainer;
234800f99ffSJeremy L Thompson   nkdm->ijacobianctxcontainer   = kdm->ijacobianctxcontainer;
235800f99ffSJeremy L Thompson   nkdm->i2functionctxcontainer  = kdm->i2functionctxcontainer;
236800f99ffSJeremy L Thompson   nkdm->i2jacobianctxcontainer  = kdm->i2jacobianctxcontainer;
237800f99ffSJeremy L Thompson   if (nkdm->rhsfunctionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm,"rhs function ctx",(PetscObject)nkdm->rhsfunctionctxcontainer));
238800f99ffSJeremy L Thompson   if (nkdm->rhsjacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm,"rhs jacobian ctx",(PetscObject)nkdm->rhsjacobianctxcontainer));
239800f99ffSJeremy L Thompson   if (nkdm->ifunctionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm,"ifunction ctx",(PetscObject)nkdm->ifunctionctxcontainer));
240800f99ffSJeremy L Thompson   if (nkdm->ijacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm,"ijacobian ctx",(PetscObject)nkdm->ijacobianctxcontainer));
241800f99ffSJeremy L Thompson   if (nkdm->i2functionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm,"i2function ctx",(PetscObject)nkdm->i2functionctxcontainer));
242800f99ffSJeremy L Thompson   if (nkdm->i2jacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm,"i2jacobian ctx",(PetscObject)nkdm->i2jacobianctxcontainer));
243d74926cbSBarry Smith 
244d74926cbSBarry Smith   nkdm->data = kdm->data;
245d74926cbSBarry Smith 
246d74926cbSBarry Smith   /*
247d74926cbSBarry Smith   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
248d74926cbSBarry Smith   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
249d74926cbSBarry Smith   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
250d74926cbSBarry Smith   */
251d74926cbSBarry Smith 
252d74926cbSBarry Smith   /* implementation specific copy hooks */
253*dbbe0bcdSBarry Smith   PetscTryTypeMethod(kdm,duplicate,nkdm);
25424989b8cSPeter Brune   PetscFunctionReturn(0);
25524989b8cSPeter Brune }
25624989b8cSPeter Brune 
25724989b8cSPeter Brune /*@C
258942e3340SBarry Smith    DMGetDMTS - get read-only private DMTS context from a DM
25924989b8cSPeter Brune 
26024989b8cSPeter Brune    Not Collective
26124989b8cSPeter Brune 
2624165533cSJose E. Roman    Input Parameter:
26324989b8cSPeter Brune .  dm - DM to be used with TS
26424989b8cSPeter Brune 
2654165533cSJose E. Roman    Output Parameter:
266942e3340SBarry Smith .  tsdm - private DMTS context
26724989b8cSPeter Brune 
26824989b8cSPeter Brune    Level: developer
26924989b8cSPeter Brune 
27024989b8cSPeter Brune    Notes:
271942e3340SBarry Smith    Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.
27224989b8cSPeter Brune 
273db781477SPatrick Sanan .seealso: `DMGetDMTSWrite()`
27424989b8cSPeter Brune @*/
275942e3340SBarry Smith PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
27624989b8cSPeter Brune {
27724989b8cSPeter Brune   PetscFunctionBegin;
27824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2792a34c10cSBarry Smith   *tsdm = (DMTS) dm->dmts;
280d74926cbSBarry Smith   if (!*tsdm) {
2819566063dSJacob Faibussowitsch     PetscCall(PetscInfo(dm,"Creating new DMTS\n"));
2829566063dSJacob Faibussowitsch     PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm));
2832a34c10cSBarry Smith     dm->dmts = (PetscObject) *tsdm;
2845c87d4f4SJunchao Zhang     (*tsdm)->originaldm = dm;
2859566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL));
2869566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL));
28724989b8cSPeter Brune   }
28824989b8cSPeter Brune   PetscFunctionReturn(0);
28924989b8cSPeter Brune }
29024989b8cSPeter Brune 
29124989b8cSPeter Brune /*@C
292942e3340SBarry Smith    DMGetDMTSWrite - get write access to private DMTS context from a DM
29324989b8cSPeter Brune 
29424989b8cSPeter Brune    Not Collective
29524989b8cSPeter Brune 
2964165533cSJose E. Roman    Input Parameter:
29724989b8cSPeter Brune .  dm - DM to be used with TS
29824989b8cSPeter Brune 
2994165533cSJose E. Roman    Output Parameter:
300942e3340SBarry Smith .  tsdm - private DMTS context
30124989b8cSPeter Brune 
30224989b8cSPeter Brune    Level: developer
30324989b8cSPeter Brune 
304db781477SPatrick Sanan .seealso: `DMGetDMTS()`
30524989b8cSPeter Brune @*/
306942e3340SBarry Smith PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
30724989b8cSPeter Brune {
308942e3340SBarry Smith   DMTS           sdm;
30924989b8cSPeter Brune 
31024989b8cSPeter Brune   PetscFunctionBegin;
31124989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3129566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&sdm));
3133c633725SBarry Smith   PetscCheck(sdm->originaldm,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DMTS has a NULL originaldm");
31424989b8cSPeter Brune   if (sdm->originaldm != dm) {  /* Copy on write */
3152a34c10cSBarry Smith     DMTS oldsdm = sdm;
3169566063dSJacob Faibussowitsch     PetscCall(PetscInfo(dm,"Copying DMTS due to write\n"));
3179566063dSJacob Faibussowitsch     PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm));
3189566063dSJacob Faibussowitsch     PetscCall(DMTSCopy(oldsdm,sdm));
3199566063dSJacob Faibussowitsch     PetscCall(DMTSDestroy((DMTS*)&dm->dmts));
3202a34c10cSBarry Smith     dm->dmts = (PetscObject) sdm;
3215c87d4f4SJunchao Zhang     sdm->originaldm = dm;
32224989b8cSPeter Brune   }
32324989b8cSPeter Brune   *tsdm = sdm;
32424989b8cSPeter Brune   PetscFunctionReturn(0);
32524989b8cSPeter Brune }
32624989b8cSPeter Brune 
32724989b8cSPeter Brune /*@C
328942e3340SBarry Smith    DMCopyDMTS - copies a DM context to a new DM
32924989b8cSPeter Brune 
33024989b8cSPeter Brune    Logically Collective
33124989b8cSPeter Brune 
3324165533cSJose E. Roman    Input Parameters:
33324989b8cSPeter Brune +  dmsrc - DM to obtain context from
33424989b8cSPeter Brune -  dmdest - DM to add context to
33524989b8cSPeter Brune 
33624989b8cSPeter Brune    Level: developer
33724989b8cSPeter Brune 
33824989b8cSPeter Brune    Note:
33924989b8cSPeter Brune    The context is copied by reference. This function does not ensure that a context exists.
34024989b8cSPeter Brune 
341db781477SPatrick Sanan .seealso: `DMGetDMTS()`, `TSSetDM()`
34224989b8cSPeter Brune @*/
343942e3340SBarry Smith PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
34424989b8cSPeter Brune {
34524989b8cSPeter Brune   PetscFunctionBegin;
34624989b8cSPeter Brune   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
34724989b8cSPeter Brune   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
3489566063dSJacob Faibussowitsch   PetscCall(DMTSDestroy((DMTS*)&dmdest->dmts));
3492a34c10cSBarry Smith   dmdest->dmts = dmsrc->dmts;
3509566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference(dmdest->dmts));
3519566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL));
3529566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL));
35324989b8cSPeter Brune   PetscFunctionReturn(0);
35424989b8cSPeter Brune }
35524989b8cSPeter Brune 
35624989b8cSPeter Brune /*@C
35724989b8cSPeter Brune    DMTSSetIFunction - set TS implicit function evaluation function
35824989b8cSPeter Brune 
35924989b8cSPeter Brune    Not Collective
36024989b8cSPeter Brune 
3614165533cSJose E. Roman    Input Parameters:
36224989b8cSPeter Brune +  dm - DM to be used with TS
363a96d6ef6SBarry Smith .  func - function evaluating f(t,u,u_t)
36424989b8cSPeter Brune -  ctx - context for residual evaluation
36524989b8cSPeter Brune 
366a96d6ef6SBarry Smith    Calling sequence of func:
367a96d6ef6SBarry Smith $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
368a96d6ef6SBarry Smith 
369a96d6ef6SBarry Smith +  t   - time at step/stage being solved
370a96d6ef6SBarry Smith .  u   - state vector
371a96d6ef6SBarry Smith .  u_t - time derivative of state vector
372a96d6ef6SBarry Smith .  F   - function vector
373a96d6ef6SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
374a96d6ef6SBarry Smith 
37524989b8cSPeter Brune    Level: advanced
37624989b8cSPeter Brune 
37724989b8cSPeter Brune    Note:
37824989b8cSPeter Brune    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
37924989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
38024989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
38124989b8cSPeter Brune 
382db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetIFunction()`, `DMTSSetJacobian()`
38324989b8cSPeter Brune @*/
38424989b8cSPeter Brune PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
38524989b8cSPeter Brune {
386942e3340SBarry Smith   DMTS           tsdm;
38724989b8cSPeter Brune 
38824989b8cSPeter Brune   PetscFunctionBegin;
38924989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3909566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
391d74926cbSBarry Smith   if (func) tsdm->ops->ifunction = func;
392800f99ffSJeremy L Thompson   if (ctx) {
393800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
394800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm),&ctxcontainer));
395800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer,ctx));
396800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm,"ifunction ctx",(PetscObject)ctxcontainer));
397800f99ffSJeremy L Thompson     tsdm->ifunctionctxcontainer = ctxcontainer;
398800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
399800f99ffSJeremy L Thompson   }
400800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
401800f99ffSJeremy L Thompson }
402800f99ffSJeremy L Thompson 
403800f99ffSJeremy L Thompson /*@C
4045cb80ecdSBarry Smith    DMTSSetIFunctionContextDestroy - set `TS` implicit evaluation context destroy function
405800f99ffSJeremy L Thompson 
406800f99ffSJeremy L Thompson    Not Collective
407800f99ffSJeremy L Thompson 
408800f99ffSJeremy L Thompson    Input Parameters:
4095cb80ecdSBarry Smith +  dm - `DM` to be used with `TS`
4105cb80ecdSBarry Smith -  f - implicit evaluation context destroy function
411800f99ffSJeremy L Thompson 
412800f99ffSJeremy L Thompson    Level: advanced
413800f99ffSJeremy L Thompson 
41490f54644SBarry Smith .seealso: `DMTSSetIFunction()`, `TSSetIFunction()`
415800f99ffSJeremy L Thompson @*/
416800f99ffSJeremy L Thompson PetscErrorCode DMTSSetIFunctionContextDestroy(DM dm,PetscErrorCode (*f)(void*))
417800f99ffSJeremy L Thompson {
418800f99ffSJeremy L Thompson   DMTS           tsdm;
419800f99ffSJeremy L Thompson 
420800f99ffSJeremy L Thompson   PetscFunctionBegin;
421800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
422800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm,&tsdm));
423800f99ffSJeremy L Thompson   if (tsdm->ifunctionctxcontainer)PetscCall(PetscContainerSetUserDestroy(tsdm->ifunctionctxcontainer,f));
424800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
425800f99ffSJeremy L Thompson }
426800f99ffSJeremy L Thompson 
427800f99ffSJeremy L Thompson PetscErrorCode DMTSUnsetIFunctionContext_Internal(DM dm)
428800f99ffSJeremy L Thompson {
429800f99ffSJeremy L Thompson   DMTS           tsdm;
430800f99ffSJeremy L Thompson 
431800f99ffSJeremy L Thompson   PetscFunctionBegin;
432800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
433800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm,&tsdm));
434800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIFunctionContext_DMTS(tsdm));
43524989b8cSPeter Brune   PetscFunctionReturn(0);
43624989b8cSPeter Brune }
43724989b8cSPeter Brune 
43824989b8cSPeter Brune /*@C
43924989b8cSPeter Brune    DMTSGetIFunction - get TS implicit residual evaluation function
44024989b8cSPeter Brune 
44124989b8cSPeter Brune    Not Collective
44224989b8cSPeter Brune 
4434165533cSJose E. Roman    Input Parameter:
44424989b8cSPeter Brune .  dm - DM to be used with TS
44524989b8cSPeter Brune 
4464165533cSJose E. Roman    Output Parameters:
44724989b8cSPeter Brune +  func - function evaluation function, see TSSetIFunction() for calling sequence
44824989b8cSPeter Brune -  ctx - context for residual evaluation
44924989b8cSPeter Brune 
45024989b8cSPeter Brune    Level: advanced
45124989b8cSPeter Brune 
45224989b8cSPeter Brune    Note:
45324989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
45424989b8cSPeter Brune    associated with the DM.
45524989b8cSPeter Brune 
456db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `DMTSSetFunction()`, `TSSetFunction()`
45724989b8cSPeter Brune @*/
45824989b8cSPeter Brune PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
45924989b8cSPeter Brune {
460942e3340SBarry Smith   DMTS           tsdm;
46124989b8cSPeter Brune 
46224989b8cSPeter Brune   PetscFunctionBegin;
46324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4649566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
465d74926cbSBarry Smith   if (func) *func = tsdm->ops->ifunction;
466800f99ffSJeremy L Thompson   if (ctx) {
467800f99ffSJeremy L Thompson     if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ifunctionctxcontainer,ctx));
468800f99ffSJeremy L Thompson     else *ctx = NULL;
469800f99ffSJeremy L Thompson   }
47024989b8cSPeter Brune   PetscFunctionReturn(0);
47124989b8cSPeter Brune }
47224989b8cSPeter Brune 
473efe9872eSLisandro Dalcin /*@C
474efe9872eSLisandro Dalcin    DMTSSetI2Function - set TS implicit function evaluation function for 2nd order systems
475efe9872eSLisandro Dalcin 
476efe9872eSLisandro Dalcin    Not Collective
477efe9872eSLisandro Dalcin 
4784165533cSJose E. Roman    Input Parameters:
479efe9872eSLisandro Dalcin +  dm - DM to be used with TS
480a96d6ef6SBarry Smith .  fun - function evaluation routine
481efe9872eSLisandro Dalcin -  ctx - context for residual evaluation
482efe9872eSLisandro Dalcin 
483a96d6ef6SBarry Smith    Calling sequence of fun:
484a96d6ef6SBarry Smith $     PetscErrorCode fun(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,Vec F,ctx);
485a96d6ef6SBarry Smith 
486a96d6ef6SBarry Smith +  t    - time at step/stage being solved
487a96d6ef6SBarry Smith .  U    - state vector
488a96d6ef6SBarry Smith .  U_t  - time derivative of state vector
489a96d6ef6SBarry Smith .  U_tt - second time derivative of state vector
490a96d6ef6SBarry Smith .  F    - function vector
491a96d6ef6SBarry Smith -  ctx  - [optional] user-defined context for matrix evaluation routine (may be NULL)
492a96d6ef6SBarry Smith 
493efe9872eSLisandro Dalcin    Level: advanced
494efe9872eSLisandro Dalcin 
495efe9872eSLisandro Dalcin    Note:
496efe9872eSLisandro Dalcin    TSSetI2Function() is normally used, but it calls this function internally because the user context is actually
497efe9872eSLisandro Dalcin    associated with the DM.
498efe9872eSLisandro Dalcin 
499db781477SPatrick Sanan .seealso: `TSSetI2Function()`
500efe9872eSLisandro Dalcin @*/
501efe9872eSLisandro Dalcin PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx)
502efe9872eSLisandro Dalcin {
503efe9872eSLisandro Dalcin   DMTS           tsdm;
504efe9872eSLisandro Dalcin 
505efe9872eSLisandro Dalcin   PetscFunctionBegin;
506efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5079566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
508efe9872eSLisandro Dalcin   if (fun) tsdm->ops->i2function = fun;
509800f99ffSJeremy L Thompson   if (ctx) {
510800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
511800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm),&ctxcontainer));
512800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer,ctx));
513800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm,"i2function ctx",(PetscObject)ctxcontainer));
514800f99ffSJeremy L Thompson     tsdm->i2functionctxcontainer = ctxcontainer;
515800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
516800f99ffSJeremy L Thompson   }
517800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
518800f99ffSJeremy L Thompson }
519800f99ffSJeremy L Thompson 
520800f99ffSJeremy L Thompson /*@C
5215cb80ecdSBarry Smith    DMTSSetI2FunctionContextDestroy - set `TS` implicit evaluation for 2nd order systems context destroy
522800f99ffSJeremy L Thompson 
523800f99ffSJeremy L Thompson    Not Collective
524800f99ffSJeremy L Thompson 
525800f99ffSJeremy L Thompson    Input Parameters:
5265cb80ecdSBarry Smith +  dm - `DM` to be used with `TS`
5275cb80ecdSBarry Smith -  f - implicit evaluation context destroy function
528800f99ffSJeremy L Thompson 
529800f99ffSJeremy L Thompson    Level: advanced
530800f99ffSJeremy L Thompson 
531800f99ffSJeremy L Thompson    Note:
5325cb80ecdSBarry Smith    `TSSetI2FunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
5335cb80ecdSBarry Smith    associated with the `DM`.
534800f99ffSJeremy L Thompson 
535800f99ffSJeremy L Thompson .seealso: `TSSetI2FunctionContextDestroy()`, `DMTSSetI2Function()`, `TSSetI2Function()`
536800f99ffSJeremy L Thompson @*/
537800f99ffSJeremy L Thompson PetscErrorCode DMTSSetI2FunctionContextDestroy(DM dm,PetscErrorCode (*f)(void*))
538800f99ffSJeremy L Thompson {
539800f99ffSJeremy L Thompson   DMTS           tsdm;
540800f99ffSJeremy L Thompson 
541800f99ffSJeremy L Thompson   PetscFunctionBegin;
542800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
543800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm,&tsdm));
544800f99ffSJeremy L Thompson   if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2functionctxcontainer,f));
545800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
546800f99ffSJeremy L Thompson }
547800f99ffSJeremy L Thompson 
548800f99ffSJeremy L Thompson PetscErrorCode DMTSUnsetI2FunctionContext_Internal(DM dm)
549800f99ffSJeremy L Thompson {
550800f99ffSJeremy L Thompson   DMTS           tsdm;
551800f99ffSJeremy L Thompson 
552800f99ffSJeremy L Thompson   PetscFunctionBegin;
553800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
554800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm,&tsdm));
555800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2FunctionContext_DMTS(tsdm));
556efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
557efe9872eSLisandro Dalcin }
558efe9872eSLisandro Dalcin 
559efe9872eSLisandro Dalcin /*@C
560efe9872eSLisandro Dalcin    DMTSGetI2Function - get TS implicit residual evaluation function for 2nd order systems
561efe9872eSLisandro Dalcin 
562efe9872eSLisandro Dalcin    Not Collective
563efe9872eSLisandro Dalcin 
5644165533cSJose E. Roman    Input Parameter:
565efe9872eSLisandro Dalcin .  dm - DM to be used with TS
566efe9872eSLisandro Dalcin 
5674165533cSJose E. Roman    Output Parameters:
568efe9872eSLisandro Dalcin +  fun - function evaluation function, see TSSetI2Function() for calling sequence
569efe9872eSLisandro Dalcin -  ctx - context for residual evaluation
570efe9872eSLisandro Dalcin 
571efe9872eSLisandro Dalcin    Level: advanced
572efe9872eSLisandro Dalcin 
573efe9872eSLisandro Dalcin    Note:
574efe9872eSLisandro Dalcin    TSGetI2Function() is normally used, but it calls this function internally because the user context is actually
575efe9872eSLisandro Dalcin    associated with the DM.
576efe9872eSLisandro Dalcin 
577c2e3fba1SPatrick Sanan .seealso: `DMTSSetI2Function()`, `TSGetI2Function()`
578efe9872eSLisandro Dalcin @*/
579efe9872eSLisandro Dalcin PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx)
580efe9872eSLisandro Dalcin {
581efe9872eSLisandro Dalcin   DMTS           tsdm;
582efe9872eSLisandro Dalcin 
583efe9872eSLisandro Dalcin   PetscFunctionBegin;
584efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5859566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
586efe9872eSLisandro Dalcin   if (fun) *fun = tsdm->ops->i2function;
587800f99ffSJeremy L Thompson   if (ctx) {
588800f99ffSJeremy L Thompson     if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2functionctxcontainer,ctx));
589800f99ffSJeremy L Thompson     else *ctx = NULL;
590800f99ffSJeremy L Thompson   }
591efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
592efe9872eSLisandro Dalcin }
593efe9872eSLisandro Dalcin 
594efe9872eSLisandro Dalcin /*@C
595efe9872eSLisandro Dalcin    DMTSSetI2Jacobian - set TS implicit Jacobian evaluation function for 2nd order systems
596efe9872eSLisandro Dalcin 
597efe9872eSLisandro Dalcin    Not Collective
598efe9872eSLisandro Dalcin 
5994165533cSJose E. Roman    Input Parameters:
600efe9872eSLisandro Dalcin +  dm - DM to be used with TS
601a96d6ef6SBarry Smith .  fun - Jacobian evaluation routine
602efe9872eSLisandro Dalcin -  ctx - context for Jacobian evaluation
603efe9872eSLisandro Dalcin 
604a96d6ef6SBarry Smith    Calling sequence of jac:
605a96d6ef6SBarry 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);
606a96d6ef6SBarry Smith 
607a96d6ef6SBarry Smith +  t    - time at step/stage being solved
608a96d6ef6SBarry Smith .  U    - state vector
609a96d6ef6SBarry Smith .  U_t  - time derivative of state vector
610a96d6ef6SBarry Smith .  U_tt - second time derivative of state vector
611a96d6ef6SBarry Smith .  v    - shift for U_t
612a96d6ef6SBarry Smith .  a    - shift for U_tt
613a96d6ef6SBarry 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
614a96d6ef6SBarry Smith .  P    - preconditioning matrix for J, may be same as J
615a96d6ef6SBarry Smith -  ctx  - [optional] user-defined context for matrix evaluation routine
616a96d6ef6SBarry Smith 
617efe9872eSLisandro Dalcin    Level: advanced
618efe9872eSLisandro Dalcin 
619efe9872eSLisandro Dalcin    Note:
620efe9872eSLisandro Dalcin    TSSetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
621efe9872eSLisandro Dalcin    associated with the DM.
622efe9872eSLisandro Dalcin 
623db781477SPatrick Sanan .seealso: `TSSetI2Jacobian()`
624efe9872eSLisandro Dalcin @*/
625efe9872eSLisandro Dalcin PetscErrorCode DMTSSetI2Jacobian(DM dm,TSI2Jacobian jac,void *ctx)
626efe9872eSLisandro Dalcin {
627efe9872eSLisandro Dalcin   DMTS           tsdm;
628efe9872eSLisandro Dalcin 
629efe9872eSLisandro Dalcin   PetscFunctionBegin;
630efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6319566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
632efe9872eSLisandro Dalcin   if (jac) tsdm->ops->i2jacobian = jac;
633800f99ffSJeremy L Thompson   if (ctx) {
634800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
635800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm),&ctxcontainer));
636800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer,ctx));
637800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm,"i2jacobian ctx",(PetscObject)ctxcontainer));
638800f99ffSJeremy L Thompson     tsdm->i2jacobianctxcontainer = ctxcontainer;
639800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
640800f99ffSJeremy L Thompson   }
641800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
642800f99ffSJeremy L Thompson }
643800f99ffSJeremy L Thompson 
644800f99ffSJeremy L Thompson /*@C
6455cb80ecdSBarry Smith    DMTSSetI2JacobianContextDestroy - set `TS` implicit Jacobian evaluation for 2nd order systems context destroy function
646800f99ffSJeremy L Thompson 
647800f99ffSJeremy L Thompson    Not Collective
648800f99ffSJeremy L Thompson 
649800f99ffSJeremy L Thompson    Input Parameters:
6505cb80ecdSBarry Smith +  dm - `DM` to be used with `TS`
6515cb80ecdSBarry Smith -  f - implicit Jacobian evaluation context destroy function
652800f99ffSJeremy L Thompson 
653800f99ffSJeremy L Thompson .seealso: `TSSetI2JacobianContextDestroy()`, `DMTSSetI2Jacobian()`, `TSSetI2Jacobian()`
654800f99ffSJeremy L Thompson @*/
655800f99ffSJeremy L Thompson PetscErrorCode DMTSSetI2JacobianContextDestroy(DM dm,PetscErrorCode (*f)(void*))
656800f99ffSJeremy L Thompson {
657800f99ffSJeremy L Thompson   DMTS           tsdm;
658800f99ffSJeremy L Thompson 
659800f99ffSJeremy L Thompson   PetscFunctionBegin;
660800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
661800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm,&tsdm));
662800f99ffSJeremy L Thompson   if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2jacobianctxcontainer,f));
663800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
664800f99ffSJeremy L Thompson }
665800f99ffSJeremy L Thompson 
666800f99ffSJeremy L Thompson PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM dm)
667800f99ffSJeremy L Thompson {
668800f99ffSJeremy L Thompson   DMTS           tsdm;
669800f99ffSJeremy L Thompson 
670800f99ffSJeremy L Thompson   PetscFunctionBegin;
671800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
672800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm,&tsdm));
673800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetI2JacobianContext_DMTS(tsdm));
674efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
675efe9872eSLisandro Dalcin }
676efe9872eSLisandro Dalcin 
677efe9872eSLisandro Dalcin /*@C
678efe9872eSLisandro Dalcin    DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems
679efe9872eSLisandro Dalcin 
680efe9872eSLisandro Dalcin    Not Collective
681efe9872eSLisandro Dalcin 
6824165533cSJose E. Roman    Input Parameter:
683efe9872eSLisandro Dalcin .  dm - DM to be used with TS
684efe9872eSLisandro Dalcin 
6854165533cSJose E. Roman    Output Parameters:
686efe9872eSLisandro Dalcin +  jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
687efe9872eSLisandro Dalcin -  ctx - context for Jacobian evaluation
688efe9872eSLisandro Dalcin 
689efe9872eSLisandro Dalcin    Level: advanced
690efe9872eSLisandro Dalcin 
691efe9872eSLisandro Dalcin    Note:
692efe9872eSLisandro Dalcin    TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
693efe9872eSLisandro Dalcin    associated with the DM.
694efe9872eSLisandro Dalcin 
695c2e3fba1SPatrick Sanan .seealso: `DMTSSetI2Jacobian()`, `TSGetI2Jacobian()`
696efe9872eSLisandro Dalcin @*/
697efe9872eSLisandro Dalcin PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx)
698efe9872eSLisandro Dalcin {
699efe9872eSLisandro Dalcin   DMTS           tsdm;
700efe9872eSLisandro Dalcin 
701efe9872eSLisandro Dalcin   PetscFunctionBegin;
702efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7039566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
704efe9872eSLisandro Dalcin   if (jac) *jac = tsdm->ops->i2jacobian;
705800f99ffSJeremy L Thompson   if (ctx) {
706800f99ffSJeremy L Thompson     if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2jacobianctxcontainer,ctx));
707800f99ffSJeremy L Thompson     else *ctx = NULL;
708800f99ffSJeremy L Thompson   }
709efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
710efe9872eSLisandro Dalcin }
71124989b8cSPeter Brune 
71224989b8cSPeter Brune /*@C
71324989b8cSPeter Brune    DMTSSetRHSFunction - set TS explicit residual evaluation function
71424989b8cSPeter Brune 
71524989b8cSPeter Brune    Not Collective
71624989b8cSPeter Brune 
7174165533cSJose E. Roman    Input Parameters:
71824989b8cSPeter Brune +  dm - DM to be used with TS
719a96d6ef6SBarry Smith .  func - RHS function evaluation routine
72024989b8cSPeter Brune -  ctx - context for residual evaluation
72124989b8cSPeter Brune 
722a96d6ef6SBarry Smith     Calling sequence of func:
723a96d6ef6SBarry Smith $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec F,void *ctx);
724a96d6ef6SBarry Smith 
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:
734800f99ffSJeremy L Thompson    TSSetRHSFunction() is normally used, but it calls this function internally because the user context is actually
73524989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
73624989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
73724989b8cSPeter Brune 
738db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetJacobian()`
73924989b8cSPeter Brune @*/
74024989b8cSPeter Brune PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
74124989b8cSPeter Brune {
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   }
756800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
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 
778800f99ffSJeremy L Thompson .seealso: `TSSetRHSFunctionContextDestroy()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
779800f99ffSJeremy L Thompson @*/
780800f99ffSJeremy L Thompson PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM dm,PetscErrorCode (*f)(void*))
781800f99ffSJeremy L Thompson {
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));
788800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
789800f99ffSJeremy L Thompson }
790800f99ffSJeremy L Thompson 
791800f99ffSJeremy L Thompson PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM dm)
792800f99ffSJeremy L Thompson {
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;
80024989b8cSPeter Brune   PetscFunctionReturn(0);
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:
809e3c11fc1SJed Brown +  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 
813a96d6ef6SBarry Smith     Calling sequence of tvar:
814a96d6ef6SBarry Smith $     PetscErrorCode tvar(TS ts,Vec p,Vec c,void *ctx);
815a96d6ef6SBarry Smith 
816a96d6ef6SBarry Smith +   ts - timestep context
8176aad120cSJose E. Roman .   p - input vector (primitive form)
818a96d6ef6SBarry Smith .   c - output vector, transient variables (conservative form)
819a96d6ef6SBarry Smith -   ctx - [optional] user-defined function context
820a96d6ef6SBarry Smith 
821e3c11fc1SJed Brown    Level: advanced
822e3c11fc1SJed Brown 
823e3c11fc1SJed Brown    Notes:
824e3c11fc1SJed Brown    This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., TSBDF)
825e3c11fc1SJed Brown    can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
826e3c11fc1SJed Brown    well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
827e3c11fc1SJed Brown    C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
828e3c11fc1SJed Brown    evaluated via the chain rule, as in
829e3c11fc1SJed Brown 
830e3c11fc1SJed Brown      dF/dP + shift * dF/dCdot dC/dP.
831e3c11fc1SJed Brown 
832db781477SPatrick Sanan .seealso: `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()`
833e3c11fc1SJed Brown @*/
834e3c11fc1SJed Brown PetscErrorCode DMTSSetTransientVariable(DM dm,TSTransientVariable tvar,void *ctx)
835e3c11fc1SJed Brown {
836e3c11fc1SJed Brown   DMTS           dmts;
837e3c11fc1SJed Brown 
838e3c11fc1SJed Brown   PetscFunctionBegin;
839e3c11fc1SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
8409566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&dmts));
841e3c11fc1SJed Brown   dmts->ops->transientvar = tvar;
842e3c11fc1SJed Brown   dmts->transientvarctx = ctx;
843e3c11fc1SJed Brown   PetscFunctionReturn(0);
844e3c11fc1SJed Brown }
845e3c11fc1SJed Brown 
846e3c11fc1SJed Brown /*@C
847e3c11fc1SJed Brown    DMTSGetTransientVariable - gets function to transform from state to transient variables
848e3c11fc1SJed Brown 
849e3c11fc1SJed Brown    Logically Collective
850e3c11fc1SJed Brown 
8514165533cSJose E. Roman    Input Parameter:
852e3c11fc1SJed Brown .  dm - DM to be used with TS
853e3c11fc1SJed Brown 
8544165533cSJose E. Roman    Output Parameters:
855a96d6ef6SBarry Smith +  tvar - a function that transforms to transient variables
856e3c11fc1SJed Brown -  ctx - a context for tvar
857e3c11fc1SJed Brown 
858e3c11fc1SJed Brown    Level: advanced
859e3c11fc1SJed Brown 
860db781477SPatrick Sanan .seealso: `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()`
861e3c11fc1SJed Brown @*/
862e3c11fc1SJed Brown PetscErrorCode DMTSGetTransientVariable(DM dm,TSTransientVariable *tvar,void *ctx)
863e3c11fc1SJed Brown {
864e3c11fc1SJed Brown   DMTS           dmts;
865e3c11fc1SJed Brown 
866e3c11fc1SJed Brown   PetscFunctionBegin;
867e3c11fc1SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
8689566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&dmts));
869e3c11fc1SJed Brown   if (tvar) *tvar = dmts->ops->transientvar;
870e3c11fc1SJed Brown   if (ctx)  *(void**)ctx = dmts->transientvarctx;
871e3c11fc1SJed Brown   PetscFunctionReturn(0);
872e3c11fc1SJed Brown }
873e3c11fc1SJed Brown 
874e3c11fc1SJed Brown /*@C
875ef20d060SBarry Smith    DMTSGetSolutionFunction - gets the TS solution evaluation function
876ef20d060SBarry Smith 
877ef20d060SBarry Smith    Not Collective
878ef20d060SBarry Smith 
8794165533cSJose E. Roman    Input Parameter:
880ef20d060SBarry Smith .  dm - DM to be used with TS
881ef20d060SBarry Smith 
882ef20d060SBarry Smith    Output Parameters:
883ef20d060SBarry Smith +  func - solution function evaluation function, see TSSetSolution() for calling sequence
884ef20d060SBarry Smith -  ctx - context for solution evaluation
885ef20d060SBarry Smith 
886ef20d060SBarry Smith    Level: advanced
887ef20d060SBarry Smith 
888db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSSetSolutionFunction()`
889ef20d060SBarry Smith @*/
890ef20d060SBarry Smith PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
891ef20d060SBarry Smith {
892942e3340SBarry Smith   DMTS           tsdm;
893ef20d060SBarry Smith 
894ef20d060SBarry Smith   PetscFunctionBegin;
895ef20d060SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
8969566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
897d74926cbSBarry Smith   if (func) *func = tsdm->ops->solution;
898ef20d060SBarry Smith   if (ctx)  *ctx  = tsdm->solutionctx;
899ef20d060SBarry Smith   PetscFunctionReturn(0);
900ef20d060SBarry Smith }
901ef20d060SBarry Smith 
902ef20d060SBarry Smith /*@C
903ef20d060SBarry Smith    DMTSSetSolutionFunction - set TS solution evaluation function
904ef20d060SBarry Smith 
905ef20d060SBarry Smith    Not Collective
906ef20d060SBarry Smith 
9074165533cSJose E. Roman    Input Parameters:
908ef20d060SBarry Smith +  dm - DM to be used with TS
909a96d6ef6SBarry Smith .  func - solution function evaluation routine
910ef20d060SBarry Smith -  ctx - context for solution evaluation
911ef20d060SBarry Smith 
912a96d6ef6SBarry Smith     Calling sequence of f:
913a96d6ef6SBarry Smith $     PetscErrorCode f(TS ts,PetscReal t,Vec u,void *ctx);
914a96d6ef6SBarry Smith 
915a96d6ef6SBarry Smith +   ts - timestep context
916a96d6ef6SBarry Smith .   t - current timestep
917a96d6ef6SBarry Smith .   u - output vector
918a96d6ef6SBarry Smith -   ctx - [optional] user-defined function context
919a96d6ef6SBarry Smith 
920ef20d060SBarry Smith    Level: advanced
921ef20d060SBarry Smith 
922ef20d060SBarry Smith    Note:
923ef20d060SBarry Smith    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
924ef20d060SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
925ef20d060SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
926ef20d060SBarry Smith 
927db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSGetSolutionFunction()`
928ef20d060SBarry Smith @*/
929ef20d060SBarry Smith PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
930ef20d060SBarry Smith {
931942e3340SBarry Smith   DMTS           tsdm;
932ef20d060SBarry Smith 
933ef20d060SBarry Smith   PetscFunctionBegin;
934ef20d060SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
9359566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
936d74926cbSBarry Smith   if (func) tsdm->ops->solution = func;
937ef20d060SBarry Smith   if (ctx)  tsdm->solutionctx   = ctx;
938ef20d060SBarry Smith   PetscFunctionReturn(0);
939ef20d060SBarry Smith }
940ef20d060SBarry Smith 
9419b7cd975SBarry Smith /*@C
9429b7cd975SBarry Smith    DMTSSetForcingFunction - set TS forcing function evaluation function
9439b7cd975SBarry Smith 
9449b7cd975SBarry Smith    Not Collective
9459b7cd975SBarry Smith 
9464165533cSJose E. Roman    Input Parameters:
9479b7cd975SBarry Smith +  dm - DM to be used with TS
948a96d6ef6SBarry Smith .  f - forcing function evaluation routine
9499b7cd975SBarry Smith -  ctx - context for solution evaluation
9509b7cd975SBarry Smith 
951a96d6ef6SBarry Smith     Calling sequence of func:
952a96d6ef6SBarry Smith $     PetscErrorCode func (TS ts,PetscReal t,Vec f,void *ctx);
953a96d6ef6SBarry Smith 
954a96d6ef6SBarry Smith +   ts - timestep context
955a96d6ef6SBarry Smith .   t - current timestep
956a96d6ef6SBarry Smith .   f - output vector
957a96d6ef6SBarry Smith -   ctx - [optional] user-defined function context
958a96d6ef6SBarry Smith 
9599b7cd975SBarry Smith    Level: advanced
9609b7cd975SBarry Smith 
9619b7cd975SBarry Smith    Note:
9629b7cd975SBarry Smith    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
9639b7cd975SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
9649b7cd975SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
9659b7cd975SBarry Smith 
966db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
9679b7cd975SBarry Smith @*/
968d56366bfSLisandro Dalcin PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx)
9699b7cd975SBarry Smith {
9709b7cd975SBarry Smith   DMTS           tsdm;
9719b7cd975SBarry Smith 
9729b7cd975SBarry Smith   PetscFunctionBegin;
9739b7cd975SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
9749566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
975f8b49ee9SBarry Smith   if (f)    tsdm->ops->forcing = f;
9769b7cd975SBarry Smith   if (ctx)  tsdm->forcingctx   = ctx;
9779b7cd975SBarry Smith   PetscFunctionReturn(0);
9789b7cd975SBarry Smith }
9799b7cd975SBarry Smith 
9809b7cd975SBarry Smith /*@C
9819b7cd975SBarry Smith    DMTSGetForcingFunction - get TS forcing function evaluation function
9829b7cd975SBarry Smith 
9839b7cd975SBarry Smith    Not Collective
9849b7cd975SBarry Smith 
9854165533cSJose E. Roman    Input Parameter:
9869b7cd975SBarry Smith .   dm - DM to be used with TS
9879b7cd975SBarry Smith 
9884165533cSJose E. Roman    Output Parameters:
989f8b49ee9SBarry Smith +  f - forcing function evaluation function; see TSForcingFunction for details
9909b7cd975SBarry Smith -  ctx - context for solution evaluation
9919b7cd975SBarry Smith 
9929b7cd975SBarry Smith    Level: advanced
9939b7cd975SBarry Smith 
9949b7cd975SBarry Smith    Note:
9959b7cd975SBarry Smith    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
9969b7cd975SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
9979b7cd975SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
9989b7cd975SBarry Smith 
999db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
10009b7cd975SBarry Smith @*/
1001d56366bfSLisandro Dalcin PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
10029b7cd975SBarry Smith {
10039b7cd975SBarry Smith   DMTS           tsdm;
10049b7cd975SBarry Smith 
10059b7cd975SBarry Smith   PetscFunctionBegin;
10069b7cd975SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
10079566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1008f8b49ee9SBarry Smith   if (f)   *f   = tsdm->ops->forcing;
10099b7cd975SBarry Smith   if (ctx) *ctx = tsdm->forcingctx;
10109b7cd975SBarry Smith   PetscFunctionReturn(0);
10119b7cd975SBarry Smith }
10129b7cd975SBarry Smith 
101324989b8cSPeter Brune /*@C
101424989b8cSPeter Brune    DMTSGetRHSFunction - get TS explicit residual evaluation function
101524989b8cSPeter Brune 
101624989b8cSPeter Brune    Not Collective
101724989b8cSPeter Brune 
10184165533cSJose E. Roman    Input Parameter:
101924989b8cSPeter Brune .  dm - DM to be used with TS
102024989b8cSPeter Brune 
10214165533cSJose E. Roman    Output Parameters:
102224989b8cSPeter Brune +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
102324989b8cSPeter Brune -  ctx - context for residual evaluation
102424989b8cSPeter Brune 
102524989b8cSPeter Brune    Level: advanced
102624989b8cSPeter Brune 
102724989b8cSPeter Brune    Note:
102824989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
102924989b8cSPeter Brune    associated with the DM.
103024989b8cSPeter Brune 
1031db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
103224989b8cSPeter Brune @*/
103324989b8cSPeter Brune PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
103424989b8cSPeter Brune {
1035942e3340SBarry Smith   DMTS           tsdm;
103624989b8cSPeter Brune 
103724989b8cSPeter Brune   PetscFunctionBegin;
103824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
10399566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
1040d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsfunction;
1041800f99ffSJeremy L Thompson   if (ctx) {
1042800f99ffSJeremy L Thompson     if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsfunctionctxcontainer,ctx));
1043800f99ffSJeremy L Thompson     else *ctx = NULL;
1044800f99ffSJeremy L Thompson   }
104524989b8cSPeter Brune   PetscFunctionReturn(0);
104624989b8cSPeter Brune }
104724989b8cSPeter Brune 
104824989b8cSPeter Brune /*@C
104924989b8cSPeter Brune    DMTSSetIJacobian - set TS Jacobian evaluation function
105024989b8cSPeter Brune 
105124989b8cSPeter Brune    Not Collective
105224989b8cSPeter Brune 
10534165533cSJose E. Roman    Input Parameters:
105424989b8cSPeter Brune +  dm - DM to be used with TS
1055a96d6ef6SBarry Smith .  func - Jacobian evaluation routine
105624989b8cSPeter Brune -  ctx - context for residual evaluation
105724989b8cSPeter Brune 
1058a96d6ef6SBarry Smith    Calling sequence of f:
1059a96d6ef6SBarry Smith $    PetscErrorCode f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx);
1060a96d6ef6SBarry Smith 
1061a96d6ef6SBarry Smith +  t    - time at step/stage being solved
1062a96d6ef6SBarry Smith .  U    - state vector
1063a96d6ef6SBarry Smith .  U_t  - time derivative of state vector
1064a96d6ef6SBarry Smith .  a    - shift
1065a96d6ef6SBarry Smith .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
1066a96d6ef6SBarry Smith .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
1067a96d6ef6SBarry Smith -  ctx  - [optional] user-defined context for matrix evaluation routine
1068a96d6ef6SBarry Smith 
106924989b8cSPeter Brune    Level: advanced
107024989b8cSPeter Brune 
107124989b8cSPeter Brune    Note:
107224989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
107324989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
107424989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
107524989b8cSPeter Brune 
1076db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSGetJacobian()`, `TSSetIJacobian()`, `TSSetIFunction()`
107724989b8cSPeter Brune @*/
107824989b8cSPeter Brune PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
107924989b8cSPeter Brune {
1080800f99ffSJeremy L Thompson   DMTS           tsdm;
108124989b8cSPeter Brune 
108224989b8cSPeter Brune   PetscFunctionBegin;
108324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1084800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1085800f99ffSJeremy L Thompson   if (func) tsdm->ops->ijacobian = func;
1086800f99ffSJeremy L Thompson   if (ctx) {
1087800f99ffSJeremy L Thompson     PetscContainer ctxcontainer;
1088800f99ffSJeremy L Thompson     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm),&ctxcontainer));
1089800f99ffSJeremy L Thompson     PetscCall(PetscContainerSetPointer(ctxcontainer,ctx));
1090800f99ffSJeremy L Thompson     PetscCall(PetscObjectCompose((PetscObject)tsdm,"ijacobian ctx",(PetscObject)ctxcontainer));
1091800f99ffSJeremy L Thompson     tsdm->ijacobianctxcontainer = ctxcontainer;
1092800f99ffSJeremy L Thompson     PetscCall(PetscContainerDestroy(&ctxcontainer));
1093800f99ffSJeremy L Thompson   }
1094800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
1095800f99ffSJeremy L Thompson }
1096800f99ffSJeremy L Thompson 
1097800f99ffSJeremy L Thompson /*@C
10985cb80ecdSBarry Smith    DMTSSetIJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function
1099800f99ffSJeremy L Thompson 
1100800f99ffSJeremy L Thompson    Not Collective
1101800f99ffSJeremy L Thompson 
1102800f99ffSJeremy L Thompson    Input Parameters:
11035cb80ecdSBarry Smith +  dm - `DM` to be used with `TS`
11045cb80ecdSBarry Smith -  f - Jacobian evaluation context destroy function
1105800f99ffSJeremy L Thompson 
1106800f99ffSJeremy L Thompson    Level: advanced
1107800f99ffSJeremy L Thompson 
1108800f99ffSJeremy L Thompson    Note:
11095cb80ecdSBarry Smith    `TSSetIJacobianContextDestroy()` is normally used, but it calls this function internally because the user context is actually
11105cb80ecdSBarry Smith    associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
11115cb80ecdSBarry Smith    not.
1112800f99ffSJeremy L Thompson 
11135cb80ecdSBarry Smith    Developer Note:
11145cb80ecdSBarry Smith    If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
11155cb80ecdSBarry Smith 
11165cb80ecdSBarry Smith .seealso: `TSSetIJacobianContextDestroy()`, `TSSetI2JacobianContextDestroy()`, `DMTSSetIJacobian()`, `TSSetIJacobian()`
1117800f99ffSJeremy L Thompson @*/
1118800f99ffSJeremy L Thompson PetscErrorCode DMTSSetIJacobianContextDestroy(DM dm,PetscErrorCode (*f)(void*))
1119800f99ffSJeremy L Thompson {
1120800f99ffSJeremy L Thompson   DMTS           tsdm;
1121800f99ffSJeremy L Thompson 
1122800f99ffSJeremy L Thompson   PetscFunctionBegin;
1123800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1124800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1125800f99ffSJeremy L Thompson   if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ijacobianctxcontainer,f));
1126800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
1127800f99ffSJeremy L Thompson }
1128800f99ffSJeremy L Thompson 
1129800f99ffSJeremy L Thompson PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM dm)
1130800f99ffSJeremy L Thompson {
1131800f99ffSJeremy L Thompson   DMTS           tsdm;
1132800f99ffSJeremy L Thompson 
1133800f99ffSJeremy L Thompson   PetscFunctionBegin;
1134800f99ffSJeremy L Thompson   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1135800f99ffSJeremy L Thompson   PetscCall(DMGetDMTSWrite(dm,&tsdm));
1136800f99ffSJeremy L Thompson   PetscCall(DMTSUnsetIJacobianContext_DMTS(tsdm));
113724989b8cSPeter Brune   PetscFunctionReturn(0);
113824989b8cSPeter Brune }
113924989b8cSPeter Brune 
114024989b8cSPeter Brune /*@C
114124989b8cSPeter Brune    DMTSGetIJacobian - get TS Jacobian evaluation function
114224989b8cSPeter Brune 
114324989b8cSPeter Brune    Not Collective
114424989b8cSPeter Brune 
11454165533cSJose E. Roman    Input Parameter:
114624989b8cSPeter Brune .  dm - DM to be used with TS
114724989b8cSPeter Brune 
11484165533cSJose E. Roman    Output Parameters:
114924989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
115024989b8cSPeter Brune -  ctx - context for residual evaluation
115124989b8cSPeter Brune 
115224989b8cSPeter Brune    Level: advanced
115324989b8cSPeter Brune 
115424989b8cSPeter Brune    Note:
115524989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
115624989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
115724989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
115824989b8cSPeter Brune 
1159db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
116024989b8cSPeter Brune @*/
116124989b8cSPeter Brune PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
116224989b8cSPeter Brune {
1163942e3340SBarry Smith   DMTS           tsdm;
116424989b8cSPeter Brune 
116524989b8cSPeter Brune   PetscFunctionBegin;
116624989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
11679566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
1168d74926cbSBarry Smith   if (func) *func = tsdm->ops->ijacobian;
1169800f99ffSJeremy L Thompson   if (ctx) {
1170800f99ffSJeremy L Thompson     if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ijacobianctxcontainer,ctx));
1171800f99ffSJeremy L Thompson     else *ctx = NULL;
1172800f99ffSJeremy L Thompson   }
117324989b8cSPeter Brune   PetscFunctionReturn(0);
117424989b8cSPeter Brune }
117524989b8cSPeter Brune 
117624989b8cSPeter Brune /*@C
117724989b8cSPeter Brune    DMTSSetRHSJacobian - set TS Jacobian evaluation function
117824989b8cSPeter Brune 
117924989b8cSPeter Brune    Not Collective
118024989b8cSPeter Brune 
11814165533cSJose E. Roman    Input Parameters:
118224989b8cSPeter Brune +  dm - DM to be used with TS
1183a96d6ef6SBarry Smith .  func - Jacobian evaluation routine
118424989b8cSPeter Brune -  ctx - context for residual evaluation
118524989b8cSPeter Brune 
1186a96d6ef6SBarry Smith    Calling sequence of func:
1187a96d6ef6SBarry Smith $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);
1188a96d6ef6SBarry Smith 
1189a96d6ef6SBarry Smith +  t - current timestep
1190a96d6ef6SBarry Smith .  u - input vector
1191a96d6ef6SBarry Smith .  Amat - (approximate) Jacobian matrix
1192a96d6ef6SBarry Smith .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1193a96d6ef6SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
1194a96d6ef6SBarry Smith 
119524989b8cSPeter Brune    Level: advanced
119624989b8cSPeter Brune 
119724989b8cSPeter Brune    Note:
119824989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
119924989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
120024989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
120124989b8cSPeter Brune 
1202db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSGetJacobian()`, `TSSetRHSJacobian()`
120324989b8cSPeter Brune @*/
120424989b8cSPeter Brune PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
120524989b8cSPeter Brune {
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   }
1220800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
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 
1237800f99ffSJeremy L Thompson .seealso: `TSSetRHSJacobianContextDestroy()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
1238800f99ffSJeremy L Thompson @*/
1239800f99ffSJeremy L Thompson PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM dm,PetscErrorCode (*f)(void*))
1240800f99ffSJeremy L Thompson {
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));
1247800f99ffSJeremy L Thompson   PetscFunctionReturn(0);
1248800f99ffSJeremy L Thompson }
1249800f99ffSJeremy L Thompson 
1250800f99ffSJeremy L Thompson PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM dm)
1251800f99ffSJeremy L Thompson {
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));
125824989b8cSPeter Brune   PetscFunctionReturn(0);
125924989b8cSPeter Brune }
126024989b8cSPeter Brune 
126124989b8cSPeter Brune /*@C
126224989b8cSPeter Brune    DMTSGetRHSJacobian - get TS Jacobian evaluation function
126324989b8cSPeter Brune 
126424989b8cSPeter Brune    Not Collective
126524989b8cSPeter Brune 
12664165533cSJose E. Roman    Input Parameter:
126724989b8cSPeter Brune .  dm - DM to be used with TS
126824989b8cSPeter Brune 
12694165533cSJose E. Roman    Output Parameters:
127024989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
127124989b8cSPeter Brune -  ctx - context for residual evaluation
127224989b8cSPeter Brune 
127324989b8cSPeter Brune    Level: advanced
127424989b8cSPeter Brune 
127524989b8cSPeter Brune    Note:
127624989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
127724989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
127824989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
127924989b8cSPeter Brune 
1280db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
128124989b8cSPeter Brune @*/
128224989b8cSPeter Brune PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
128324989b8cSPeter Brune {
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   }
129424989b8cSPeter Brune   PetscFunctionReturn(0);
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:
1303ad6bc421SBarry 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 
1309db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1310ad6bc421SBarry Smith @*/
1311ad6bc421SBarry Smith PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
1312ad6bc421SBarry Smith {
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;
1320ad6bc421SBarry Smith   PetscFunctionReturn(0);
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:
1329ad6bc421SBarry 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 
1335db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
1336ad6bc421SBarry Smith @*/
1337ad6bc421SBarry Smith PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
1338ad6bc421SBarry Smith {
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;
1346ad6bc421SBarry Smith   PetscFunctionReturn(0);
1347ad6bc421SBarry Smith }
1348