xref: /petsc/src/ts/utils/dmts.c (revision c2e3fba1fe1cda7e6350bbca19c4ed35ce95940a)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
2af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
324989b8cSPeter Brune 
4d74926cbSBarry Smith static PetscErrorCode DMTSDestroy(DMTS *kdm)
5d74926cbSBarry Smith {
6d74926cbSBarry Smith   PetscFunctionBegin;
7d74926cbSBarry Smith   if (!*kdm) PetscFunctionReturn(0);
8d74926cbSBarry Smith   PetscValidHeaderSpecific((*kdm),DMTS_CLASSID,1);
9c793f718SLisandro Dalcin   if (--((PetscObject)(*kdm))->refct > 0) {*kdm = NULL; PetscFunctionReturn(0);}
109566063dSJacob Faibussowitsch   if ((*kdm)->ops->destroy) PetscCall(((*kdm)->ops->destroy)(*kdm));
119566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(kdm));
12d74926cbSBarry Smith   PetscFunctionReturn(0);
13d74926cbSBarry Smith }
14d74926cbSBarry Smith 
152d53ad75SBarry Smith PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer)
162d53ad75SBarry Smith {
172d53ad75SBarry Smith   PetscFunctionBegin;
189566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,NULL,PETSC_FUNCTION));
199566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,NULL,PETSC_FUNCTION));
209566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,NULL,PETSC_FUNCTION));
21ad6bc421SBarry Smith   if (kdm->ops->ifunctionload) {
229566063dSJacob Faibussowitsch     PetscCall((*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer));
23ad6bc421SBarry Smith   }
249566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,NULL,PETSC_FUNCTION));
259566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,NULL,PETSC_FUNCTION));
269566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,NULL,PETSC_FUNCTION));
27ad6bc421SBarry Smith   if (kdm->ops->ijacobianload) {
289566063dSJacob Faibussowitsch     PetscCall((*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer));
29ad6bc421SBarry Smith   }
302d53ad75SBarry Smith   PetscFunctionReturn(0);
312d53ad75SBarry Smith }
322d53ad75SBarry Smith 
332d53ad75SBarry Smith PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer)
342d53ad75SBarry Smith {
352d53ad75SBarry Smith   PetscBool      isascii,isbinary;
362d53ad75SBarry Smith 
372d53ad75SBarry Smith   PetscFunctionBegin;
389566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
399566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
402d53ad75SBarry Smith   if (isascii) {
41c7a10e08SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
422d53ad75SBarry Smith     const char *fname;
432d53ad75SBarry Smith 
449566063dSJacob Faibussowitsch     PetscCall(PetscFPTFind(kdm->ops->ifunction,&fname));
452d53ad75SBarry Smith     if (fname) {
469566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  IFunction used by TS: %s\n",fname));
472d53ad75SBarry Smith     }
489566063dSJacob Faibussowitsch     PetscCall(PetscFPTFind(kdm->ops->ijacobian,&fname));
492d53ad75SBarry Smith     if (fname) {
509566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  IJacobian function used by TS: %s\n",fname));
512d53ad75SBarry Smith     }
52c7a10e08SBarry Smith #endif
532d53ad75SBarry Smith   } else if (isbinary) {
543964eb88SJed Brown     struct {
553964eb88SJed Brown       TSIFunction ifunction;
569200755eSBarry Smith     } funcstruct;
579200755eSBarry Smith     struct {
583964eb88SJed Brown       PetscErrorCode (*ifunctionview)(void*,PetscViewer);
599200755eSBarry Smith     } funcviewstruct;
609200755eSBarry Smith     struct {
613964eb88SJed Brown       PetscErrorCode (*ifunctionload)(void**,PetscViewer);
629200755eSBarry Smith     } funcloadstruct;
633964eb88SJed Brown     struct {
643964eb88SJed Brown       TSIJacobian ijacobian;
659200755eSBarry Smith     } jacstruct;
669200755eSBarry Smith     struct {
673964eb88SJed Brown       PetscErrorCode (*ijacobianview)(void*,PetscViewer);
689200755eSBarry Smith     } jacviewstruct;
699200755eSBarry Smith     struct {
703964eb88SJed Brown       PetscErrorCode (*ijacobianload)(void**,PetscViewer);
719200755eSBarry Smith     } jacloadstruct;
723964eb88SJed Brown 
739200755eSBarry Smith     funcstruct.ifunction         = kdm->ops->ifunction;
749200755eSBarry Smith     funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
759200755eSBarry Smith     funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
769566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION));
779566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION));
789566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION));
79ad6bc421SBarry Smith     if (kdm->ops->ifunctionview) {
809566063dSJacob Faibussowitsch       PetscCall((*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer));
81ad6bc421SBarry Smith     }
829200755eSBarry Smith     jacstruct.ijacobian = kdm->ops->ijacobian;
839200755eSBarry Smith     jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
849200755eSBarry Smith     jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
859566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION));
869566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION));
879566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION));
88ad6bc421SBarry Smith     if (kdm->ops->ijacobianview) {
899566063dSJacob Faibussowitsch       PetscCall((*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer));
90ad6bc421SBarry Smith     }
912d53ad75SBarry Smith   }
922d53ad75SBarry Smith   PetscFunctionReturn(0);
932d53ad75SBarry Smith }
942d53ad75SBarry Smith 
95d74926cbSBarry Smith static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
96d74926cbSBarry Smith {
97d74926cbSBarry Smith   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCall(TSInitializePackage());
999566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView));
100d74926cbSBarry Smith   PetscFunctionReturn(0);
101d74926cbSBarry Smith }
10224989b8cSPeter Brune 
1032a34c10cSBarry Smith /* Attaches the DMTS to the coarse level.
10424989b8cSPeter Brune  * Under what conditions should we copy versus duplicate?
10524989b8cSPeter Brune  */
106942e3340SBarry Smith static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
10724989b8cSPeter Brune {
10824989b8cSPeter Brune   PetscFunctionBegin;
1099566063dSJacob Faibussowitsch   PetscCall(DMCopyDMTS(dm,dmc));
11024989b8cSPeter Brune   PetscFunctionReturn(0);
11124989b8cSPeter Brune }
11224989b8cSPeter Brune 
11324989b8cSPeter Brune /* This could restrict auxiliary information to the coarse level.
11424989b8cSPeter Brune  */
115942e3340SBarry Smith static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
11624989b8cSPeter Brune {
11724989b8cSPeter Brune   PetscFunctionBegin;
11824989b8cSPeter Brune   PetscFunctionReturn(0);
11924989b8cSPeter Brune }
12024989b8cSPeter Brune 
121258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
122258e1594SPeter Brune {
123258e1594SPeter Brune   PetscFunctionBegin;
1249566063dSJacob Faibussowitsch   PetscCall(DMCopyDMTS(dm,subdm));
125258e1594SPeter Brune   PetscFunctionReturn(0);
126258e1594SPeter Brune }
127258e1594SPeter Brune 
128258e1594SPeter Brune /* This could restrict auxiliary information to the coarse level.
129258e1594SPeter Brune  */
130258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
131258e1594SPeter Brune {
132258e1594SPeter Brune   PetscFunctionBegin;
133258e1594SPeter Brune   PetscFunctionReturn(0);
134258e1594SPeter Brune }
135258e1594SPeter Brune 
136d74926cbSBarry Smith /*@C
137d74926cbSBarry Smith    DMTSCopy - copies the information in a DMTS to another DMTS
138d74926cbSBarry Smith 
139d74926cbSBarry Smith    Not Collective
140d74926cbSBarry Smith 
1414165533cSJose E. Roman    Input Parameters:
142d74926cbSBarry Smith +  kdm - Original DMTS
143d74926cbSBarry Smith -  nkdm - DMTS to receive the data, should have been created with DMTSCreate()
144d74926cbSBarry Smith 
145d74926cbSBarry Smith    Level: developer
146d74926cbSBarry Smith 
147db781477SPatrick Sanan .seealso: `DMTSCreate()`, `DMTSDestroy()`
148d74926cbSBarry Smith @*/
149d74926cbSBarry Smith PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
15024989b8cSPeter Brune {
15124989b8cSPeter Brune   PetscFunctionBegin;
152d74926cbSBarry Smith   PetscValidHeaderSpecific(kdm,DMTS_CLASSID,1);
153d74926cbSBarry Smith   PetscValidHeaderSpecific(nkdm,DMTS_CLASSID,2);
154d74926cbSBarry Smith   nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
155d74926cbSBarry Smith   nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
156d74926cbSBarry Smith   nkdm->ops->ifunction   = kdm->ops->ifunction;
157d74926cbSBarry Smith   nkdm->ops->ijacobian   = kdm->ops->ijacobian;
158efe9872eSLisandro Dalcin   nkdm->ops->i2function  = kdm->ops->i2function;
159efe9872eSLisandro Dalcin   nkdm->ops->i2jacobian  = kdm->ops->i2jacobian;
160d74926cbSBarry Smith   nkdm->ops->solution    = kdm->ops->solution;
161d74926cbSBarry Smith   nkdm->ops->destroy     = kdm->ops->destroy;
162d74926cbSBarry Smith   nkdm->ops->duplicate   = kdm->ops->duplicate;
163d74926cbSBarry Smith 
164d74926cbSBarry Smith   nkdm->rhsfunctionctx = kdm->rhsfunctionctx;
165d74926cbSBarry Smith   nkdm->rhsjacobianctx = kdm->rhsjacobianctx;
166d74926cbSBarry Smith   nkdm->ifunctionctx   = kdm->ifunctionctx;
167d74926cbSBarry Smith   nkdm->ijacobianctx   = kdm->ijacobianctx;
168efe9872eSLisandro Dalcin   nkdm->i2functionctx  = kdm->i2functionctx;
169efe9872eSLisandro Dalcin   nkdm->i2jacobianctx  = kdm->i2jacobianctx;
170d74926cbSBarry Smith   nkdm->solutionctx    = kdm->solutionctx;
171d74926cbSBarry Smith 
172d74926cbSBarry Smith   nkdm->data = kdm->data;
173d74926cbSBarry Smith 
174d74926cbSBarry Smith   /*
175d74926cbSBarry Smith   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
176d74926cbSBarry Smith   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
177d74926cbSBarry Smith   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
178d74926cbSBarry Smith   */
179d74926cbSBarry Smith 
180d74926cbSBarry Smith   /* implementation specific copy hooks */
1819566063dSJacob Faibussowitsch   if (kdm->ops->duplicate) PetscCall((*kdm->ops->duplicate)(kdm,nkdm));
18224989b8cSPeter Brune   PetscFunctionReturn(0);
18324989b8cSPeter Brune }
18424989b8cSPeter Brune 
18524989b8cSPeter Brune /*@C
186942e3340SBarry Smith    DMGetDMTS - get read-only private DMTS context from a DM
18724989b8cSPeter Brune 
18824989b8cSPeter Brune    Not Collective
18924989b8cSPeter Brune 
1904165533cSJose E. Roman    Input Parameter:
19124989b8cSPeter Brune .  dm - DM to be used with TS
19224989b8cSPeter Brune 
1934165533cSJose E. Roman    Output Parameter:
194942e3340SBarry Smith .  tsdm - private DMTS context
19524989b8cSPeter Brune 
19624989b8cSPeter Brune    Level: developer
19724989b8cSPeter Brune 
19824989b8cSPeter Brune    Notes:
199942e3340SBarry Smith    Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.
20024989b8cSPeter Brune 
201db781477SPatrick Sanan .seealso: `DMGetDMTSWrite()`
20224989b8cSPeter Brune @*/
203942e3340SBarry Smith PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
20424989b8cSPeter Brune {
20524989b8cSPeter Brune   PetscFunctionBegin;
20624989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2072a34c10cSBarry Smith   *tsdm = (DMTS) dm->dmts;
208d74926cbSBarry Smith   if (!*tsdm) {
2099566063dSJacob Faibussowitsch     PetscCall(PetscInfo(dm,"Creating new DMTS\n"));
2109566063dSJacob Faibussowitsch     PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm));
2112a34c10cSBarry Smith     dm->dmts = (PetscObject) *tsdm;
2125c87d4f4SJunchao Zhang     (*tsdm)->originaldm = dm;
2139566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL));
2149566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL));
21524989b8cSPeter Brune   }
21624989b8cSPeter Brune   PetscFunctionReturn(0);
21724989b8cSPeter Brune }
21824989b8cSPeter Brune 
21924989b8cSPeter Brune /*@C
220942e3340SBarry Smith    DMGetDMTSWrite - get write access to private DMTS context from a DM
22124989b8cSPeter Brune 
22224989b8cSPeter Brune    Not Collective
22324989b8cSPeter Brune 
2244165533cSJose E. Roman    Input Parameter:
22524989b8cSPeter Brune .  dm - DM to be used with TS
22624989b8cSPeter Brune 
2274165533cSJose E. Roman    Output Parameter:
228942e3340SBarry Smith .  tsdm - private DMTS context
22924989b8cSPeter Brune 
23024989b8cSPeter Brune    Level: developer
23124989b8cSPeter Brune 
232db781477SPatrick Sanan .seealso: `DMGetDMTS()`
23324989b8cSPeter Brune @*/
234942e3340SBarry Smith PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
23524989b8cSPeter Brune {
236942e3340SBarry Smith   DMTS           sdm;
23724989b8cSPeter Brune 
23824989b8cSPeter Brune   PetscFunctionBegin;
23924989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2409566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&sdm));
2413c633725SBarry Smith   PetscCheck(sdm->originaldm,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DMTS has a NULL originaldm");
24224989b8cSPeter Brune   if (sdm->originaldm != dm) {  /* Copy on write */
2432a34c10cSBarry Smith     DMTS oldsdm = sdm;
2449566063dSJacob Faibussowitsch     PetscCall(PetscInfo(dm,"Copying DMTS due to write\n"));
2459566063dSJacob Faibussowitsch     PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm));
2469566063dSJacob Faibussowitsch     PetscCall(DMTSCopy(oldsdm,sdm));
2479566063dSJacob Faibussowitsch     PetscCall(DMTSDestroy((DMTS*)&dm->dmts));
2482a34c10cSBarry Smith     dm->dmts = (PetscObject) sdm;
2495c87d4f4SJunchao Zhang     sdm->originaldm = dm;
25024989b8cSPeter Brune   }
25124989b8cSPeter Brune   *tsdm = sdm;
25224989b8cSPeter Brune   PetscFunctionReturn(0);
25324989b8cSPeter Brune }
25424989b8cSPeter Brune 
25524989b8cSPeter Brune /*@C
256942e3340SBarry Smith    DMCopyDMTS - copies a DM context to a new DM
25724989b8cSPeter Brune 
25824989b8cSPeter Brune    Logically Collective
25924989b8cSPeter Brune 
2604165533cSJose E. Roman    Input Parameters:
26124989b8cSPeter Brune +  dmsrc - DM to obtain context from
26224989b8cSPeter Brune -  dmdest - DM to add context to
26324989b8cSPeter Brune 
26424989b8cSPeter Brune    Level: developer
26524989b8cSPeter Brune 
26624989b8cSPeter Brune    Note:
26724989b8cSPeter Brune    The context is copied by reference. This function does not ensure that a context exists.
26824989b8cSPeter Brune 
269db781477SPatrick Sanan .seealso: `DMGetDMTS()`, `TSSetDM()`
27024989b8cSPeter Brune @*/
271942e3340SBarry Smith PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
27224989b8cSPeter Brune {
27324989b8cSPeter Brune   PetscFunctionBegin;
27424989b8cSPeter Brune   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
27524989b8cSPeter Brune   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
2769566063dSJacob Faibussowitsch   PetscCall(DMTSDestroy((DMTS*)&dmdest->dmts));
2772a34c10cSBarry Smith   dmdest->dmts = dmsrc->dmts;
2789566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference(dmdest->dmts));
2799566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL));
2809566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL));
28124989b8cSPeter Brune   PetscFunctionReturn(0);
28224989b8cSPeter Brune }
28324989b8cSPeter Brune 
28424989b8cSPeter Brune /*@C
28524989b8cSPeter Brune    DMTSSetIFunction - set TS implicit function evaluation function
28624989b8cSPeter Brune 
28724989b8cSPeter Brune    Not Collective
28824989b8cSPeter Brune 
2894165533cSJose E. Roman    Input Parameters:
29024989b8cSPeter Brune +  dm - DM to be used with TS
291a96d6ef6SBarry Smith .  func - function evaluating f(t,u,u_t)
29224989b8cSPeter Brune -  ctx - context for residual evaluation
29324989b8cSPeter Brune 
294a96d6ef6SBarry Smith    Calling sequence of func:
295a96d6ef6SBarry Smith $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
296a96d6ef6SBarry Smith 
297a96d6ef6SBarry Smith +  t   - time at step/stage being solved
298a96d6ef6SBarry Smith .  u   - state vector
299a96d6ef6SBarry Smith .  u_t - time derivative of state vector
300a96d6ef6SBarry Smith .  F   - function vector
301a96d6ef6SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
302a96d6ef6SBarry Smith 
30324989b8cSPeter Brune    Level: advanced
30424989b8cSPeter Brune 
30524989b8cSPeter Brune    Note:
30624989b8cSPeter Brune    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
30724989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
30824989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
30924989b8cSPeter Brune 
310db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetIFunction()`, `DMTSSetJacobian()`
31124989b8cSPeter Brune @*/
31224989b8cSPeter Brune PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
31324989b8cSPeter Brune {
314942e3340SBarry Smith   DMTS           tsdm;
31524989b8cSPeter Brune 
31624989b8cSPeter Brune   PetscFunctionBegin;
31724989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3189566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
319d74926cbSBarry Smith   if (func) tsdm->ops->ifunction = func;
32024989b8cSPeter Brune   if (ctx)  tsdm->ifunctionctx = ctx;
32124989b8cSPeter Brune   PetscFunctionReturn(0);
32224989b8cSPeter Brune }
32324989b8cSPeter Brune 
32424989b8cSPeter Brune /*@C
32524989b8cSPeter Brune    DMTSGetIFunction - get TS implicit residual evaluation function
32624989b8cSPeter Brune 
32724989b8cSPeter Brune    Not Collective
32824989b8cSPeter Brune 
3294165533cSJose E. Roman    Input Parameter:
33024989b8cSPeter Brune .  dm - DM to be used with TS
33124989b8cSPeter Brune 
3324165533cSJose E. Roman    Output Parameters:
33324989b8cSPeter Brune +  func - function evaluation function, see TSSetIFunction() for calling sequence
33424989b8cSPeter Brune -  ctx - context for residual evaluation
33524989b8cSPeter Brune 
33624989b8cSPeter Brune    Level: advanced
33724989b8cSPeter Brune 
33824989b8cSPeter Brune    Note:
33924989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
34024989b8cSPeter Brune    associated with the DM.
34124989b8cSPeter Brune 
342db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `DMTSSetFunction()`, `TSSetFunction()`
34324989b8cSPeter Brune @*/
34424989b8cSPeter Brune PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
34524989b8cSPeter Brune {
346942e3340SBarry Smith   DMTS           tsdm;
34724989b8cSPeter Brune 
34824989b8cSPeter Brune   PetscFunctionBegin;
34924989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3509566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
351d74926cbSBarry Smith   if (func) *func = tsdm->ops->ifunction;
35224989b8cSPeter Brune   if (ctx)  *ctx = tsdm->ifunctionctx;
35324989b8cSPeter Brune   PetscFunctionReturn(0);
35424989b8cSPeter Brune }
35524989b8cSPeter Brune 
356efe9872eSLisandro Dalcin /*@C
357efe9872eSLisandro Dalcin    DMTSSetI2Function - set TS implicit function evaluation function for 2nd order systems
358efe9872eSLisandro Dalcin 
359efe9872eSLisandro Dalcin    Not Collective
360efe9872eSLisandro Dalcin 
3614165533cSJose E. Roman    Input Parameters:
362efe9872eSLisandro Dalcin +  dm - DM to be used with TS
363a96d6ef6SBarry Smith .  fun - function evaluation routine
364efe9872eSLisandro Dalcin -  ctx - context for residual evaluation
365efe9872eSLisandro Dalcin 
366a96d6ef6SBarry Smith    Calling sequence of fun:
367a96d6ef6SBarry Smith $     PetscErrorCode fun(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,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 .  U_tt - second time derivative of state vector
373a96d6ef6SBarry Smith .  F    - function vector
374a96d6ef6SBarry Smith -  ctx  - [optional] user-defined context for matrix evaluation routine (may be NULL)
375a96d6ef6SBarry Smith 
376efe9872eSLisandro Dalcin    Level: advanced
377efe9872eSLisandro Dalcin 
378efe9872eSLisandro Dalcin    Note:
379efe9872eSLisandro Dalcin    TSSetI2Function() is normally used, but it calls this function internally because the user context is actually
380efe9872eSLisandro Dalcin    associated with the DM.
381efe9872eSLisandro Dalcin 
382db781477SPatrick Sanan .seealso: `TSSetI2Function()`
383efe9872eSLisandro Dalcin @*/
384efe9872eSLisandro Dalcin PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx)
385efe9872eSLisandro Dalcin {
386efe9872eSLisandro Dalcin   DMTS           tsdm;
387efe9872eSLisandro Dalcin 
388efe9872eSLisandro Dalcin   PetscFunctionBegin;
389efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3909566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
391efe9872eSLisandro Dalcin   if (fun) tsdm->ops->i2function = fun;
392efe9872eSLisandro Dalcin   if (ctx) tsdm->i2functionctx   = ctx;
393efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
394efe9872eSLisandro Dalcin }
395efe9872eSLisandro Dalcin 
396efe9872eSLisandro Dalcin /*@C
397efe9872eSLisandro Dalcin    DMTSGetI2Function - get TS implicit residual evaluation function for 2nd order systems
398efe9872eSLisandro Dalcin 
399efe9872eSLisandro Dalcin    Not Collective
400efe9872eSLisandro Dalcin 
4014165533cSJose E. Roman    Input Parameter:
402efe9872eSLisandro Dalcin .  dm - DM to be used with TS
403efe9872eSLisandro Dalcin 
4044165533cSJose E. Roman    Output Parameters:
405efe9872eSLisandro Dalcin +  fun - function evaluation function, see TSSetI2Function() for calling sequence
406efe9872eSLisandro Dalcin -  ctx - context for residual evaluation
407efe9872eSLisandro Dalcin 
408efe9872eSLisandro Dalcin    Level: advanced
409efe9872eSLisandro Dalcin 
410efe9872eSLisandro Dalcin    Note:
411efe9872eSLisandro Dalcin    TSGetI2Function() is normally used, but it calls this function internally because the user context is actually
412efe9872eSLisandro Dalcin    associated with the DM.
413efe9872eSLisandro Dalcin 
414*c2e3fba1SPatrick Sanan .seealso: `DMTSSetI2Function()`, `TSGetI2Function()`
415efe9872eSLisandro Dalcin @*/
416efe9872eSLisandro Dalcin PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx)
417efe9872eSLisandro Dalcin {
418efe9872eSLisandro Dalcin   DMTS           tsdm;
419efe9872eSLisandro Dalcin 
420efe9872eSLisandro Dalcin   PetscFunctionBegin;
421efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4229566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
423efe9872eSLisandro Dalcin   if (fun) *fun = tsdm->ops->i2function;
424efe9872eSLisandro Dalcin   if (ctx) *ctx = tsdm->i2functionctx;
425efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
426efe9872eSLisandro Dalcin }
427efe9872eSLisandro Dalcin 
428efe9872eSLisandro Dalcin /*@C
429efe9872eSLisandro Dalcin    DMTSSetI2Jacobian - set TS implicit Jacobian evaluation function for 2nd order systems
430efe9872eSLisandro Dalcin 
431efe9872eSLisandro Dalcin    Not Collective
432efe9872eSLisandro Dalcin 
4334165533cSJose E. Roman    Input Parameters:
434efe9872eSLisandro Dalcin +  dm - DM to be used with TS
435a96d6ef6SBarry Smith .  fun - Jacobian evaluation routine
436efe9872eSLisandro Dalcin -  ctx - context for Jacobian evaluation
437efe9872eSLisandro Dalcin 
438a96d6ef6SBarry Smith    Calling sequence of jac:
439a96d6ef6SBarry 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);
440a96d6ef6SBarry Smith 
441a96d6ef6SBarry Smith +  t    - time at step/stage being solved
442a96d6ef6SBarry Smith .  U    - state vector
443a96d6ef6SBarry Smith .  U_t  - time derivative of state vector
444a96d6ef6SBarry Smith .  U_tt - second time derivative of state vector
445a96d6ef6SBarry Smith .  v    - shift for U_t
446a96d6ef6SBarry Smith .  a    - shift for U_tt
447a96d6ef6SBarry 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
448a96d6ef6SBarry Smith .  P    - preconditioning matrix for J, may be same as J
449a96d6ef6SBarry Smith -  ctx  - [optional] user-defined context for matrix evaluation routine
450a96d6ef6SBarry Smith 
451efe9872eSLisandro Dalcin    Level: advanced
452efe9872eSLisandro Dalcin 
453efe9872eSLisandro Dalcin    Note:
454efe9872eSLisandro Dalcin    TSSetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
455efe9872eSLisandro Dalcin    associated with the DM.
456efe9872eSLisandro Dalcin 
457db781477SPatrick Sanan .seealso: `TSSetI2Jacobian()`
458efe9872eSLisandro Dalcin @*/
459efe9872eSLisandro Dalcin PetscErrorCode DMTSSetI2Jacobian(DM dm,TSI2Jacobian jac,void *ctx)
460efe9872eSLisandro Dalcin {
461efe9872eSLisandro Dalcin   DMTS           tsdm;
462efe9872eSLisandro Dalcin 
463efe9872eSLisandro Dalcin   PetscFunctionBegin;
464efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4659566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
466efe9872eSLisandro Dalcin   if (jac) tsdm->ops->i2jacobian = jac;
467efe9872eSLisandro Dalcin   if (ctx) tsdm->i2jacobianctx   = ctx;
468efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
469efe9872eSLisandro Dalcin }
470efe9872eSLisandro Dalcin 
471efe9872eSLisandro Dalcin /*@C
472efe9872eSLisandro Dalcin    DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems
473efe9872eSLisandro Dalcin 
474efe9872eSLisandro Dalcin    Not Collective
475efe9872eSLisandro Dalcin 
4764165533cSJose E. Roman    Input Parameter:
477efe9872eSLisandro Dalcin .  dm - DM to be used with TS
478efe9872eSLisandro Dalcin 
4794165533cSJose E. Roman    Output Parameters:
480efe9872eSLisandro Dalcin +  jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
481efe9872eSLisandro Dalcin -  ctx - context for Jacobian evaluation
482efe9872eSLisandro Dalcin 
483efe9872eSLisandro Dalcin    Level: advanced
484efe9872eSLisandro Dalcin 
485efe9872eSLisandro Dalcin    Note:
486efe9872eSLisandro Dalcin    TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
487efe9872eSLisandro Dalcin    associated with the DM.
488efe9872eSLisandro Dalcin 
489*c2e3fba1SPatrick Sanan .seealso: `DMTSSetI2Jacobian()`, `TSGetI2Jacobian()`
490efe9872eSLisandro Dalcin @*/
491efe9872eSLisandro Dalcin PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx)
492efe9872eSLisandro Dalcin {
493efe9872eSLisandro Dalcin   DMTS           tsdm;
494efe9872eSLisandro Dalcin 
495efe9872eSLisandro Dalcin   PetscFunctionBegin;
496efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4979566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
498efe9872eSLisandro Dalcin   if (jac) *jac = tsdm->ops->i2jacobian;
499efe9872eSLisandro Dalcin   if (ctx) *ctx = tsdm->i2jacobianctx;
500efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
501efe9872eSLisandro Dalcin }
50224989b8cSPeter Brune 
50324989b8cSPeter Brune /*@C
50424989b8cSPeter Brune    DMTSSetRHSFunction - set TS explicit residual evaluation function
50524989b8cSPeter Brune 
50624989b8cSPeter Brune    Not Collective
50724989b8cSPeter Brune 
5084165533cSJose E. Roman    Input Parameters:
50924989b8cSPeter Brune +  dm - DM to be used with TS
510a96d6ef6SBarry Smith .  func - RHS function evaluation routine
51124989b8cSPeter Brune -  ctx - context for residual evaluation
51224989b8cSPeter Brune 
513a96d6ef6SBarry Smith     Calling sequence of func:
514a96d6ef6SBarry Smith $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec F,void *ctx);
515a96d6ef6SBarry Smith 
516a96d6ef6SBarry Smith +   ts - timestep context
517a96d6ef6SBarry Smith .   t - current timestep
518a96d6ef6SBarry Smith .   u - input vector
519a96d6ef6SBarry Smith .   F - function vector
520a96d6ef6SBarry Smith -   ctx - [optional] user-defined function context
521a96d6ef6SBarry Smith 
52224989b8cSPeter Brune    Level: advanced
52324989b8cSPeter Brune 
52424989b8cSPeter Brune    Note:
525ef20d060SBarry Smith    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
52624989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
52724989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
52824989b8cSPeter Brune 
529db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetJacobian()`
53024989b8cSPeter Brune @*/
53124989b8cSPeter Brune PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
53224989b8cSPeter Brune {
533942e3340SBarry Smith   DMTS           tsdm;
53424989b8cSPeter Brune 
53524989b8cSPeter Brune   PetscFunctionBegin;
53624989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5379566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
538d74926cbSBarry Smith   if (func) tsdm->ops->rhsfunction = func;
53924989b8cSPeter Brune   if (ctx)  tsdm->rhsfunctionctx = ctx;
54024989b8cSPeter Brune   PetscFunctionReturn(0);
54124989b8cSPeter Brune }
54224989b8cSPeter Brune 
543ef20d060SBarry Smith /*@C
544e3c11fc1SJed Brown    DMTSSetTransientVariable - sets function to transform from state to transient variables
545e3c11fc1SJed Brown 
546e3c11fc1SJed Brown    Logically Collective
547e3c11fc1SJed Brown 
5484165533cSJose E. Roman    Input Parameters:
549e3c11fc1SJed Brown +  dm - DM to be used with TS
550a96d6ef6SBarry Smith .  tvar - a function that transforms to transient variables
551e3c11fc1SJed Brown -  ctx - a context for tvar
552e3c11fc1SJed Brown 
553a96d6ef6SBarry Smith     Calling sequence of tvar:
554a96d6ef6SBarry Smith $     PetscErrorCode tvar(TS ts,Vec p,Vec c,void *ctx);
555a96d6ef6SBarry Smith 
556a96d6ef6SBarry Smith +   ts - timestep context
557a96d6ef6SBarry Smith .   p - input vector (primative form)
558a96d6ef6SBarry Smith .   c - output vector, transient variables (conservative form)
559a96d6ef6SBarry Smith -   ctx - [optional] user-defined function context
560a96d6ef6SBarry Smith 
561e3c11fc1SJed Brown    Level: advanced
562e3c11fc1SJed Brown 
563e3c11fc1SJed Brown    Notes:
564e3c11fc1SJed Brown    This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., TSBDF)
565e3c11fc1SJed Brown    can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
566e3c11fc1SJed Brown    well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
567e3c11fc1SJed Brown    C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
568e3c11fc1SJed Brown    evaluated via the chain rule, as in
569e3c11fc1SJed Brown 
570e3c11fc1SJed Brown      dF/dP + shift * dF/dCdot dC/dP.
571e3c11fc1SJed Brown 
572db781477SPatrick Sanan .seealso: `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()`
573e3c11fc1SJed Brown @*/
574e3c11fc1SJed Brown PetscErrorCode DMTSSetTransientVariable(DM dm,TSTransientVariable tvar,void *ctx)
575e3c11fc1SJed Brown {
576e3c11fc1SJed Brown   DMTS           dmts;
577e3c11fc1SJed Brown 
578e3c11fc1SJed Brown   PetscFunctionBegin;
579e3c11fc1SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5809566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&dmts));
581e3c11fc1SJed Brown   dmts->ops->transientvar = tvar;
582e3c11fc1SJed Brown   dmts->transientvarctx = ctx;
583e3c11fc1SJed Brown   PetscFunctionReturn(0);
584e3c11fc1SJed Brown }
585e3c11fc1SJed Brown 
586e3c11fc1SJed Brown /*@C
587e3c11fc1SJed Brown    DMTSGetTransientVariable - gets function to transform from state to transient variables
588e3c11fc1SJed Brown 
589e3c11fc1SJed Brown    Logically Collective
590e3c11fc1SJed Brown 
5914165533cSJose E. Roman    Input Parameter:
592e3c11fc1SJed Brown .  dm - DM to be used with TS
593e3c11fc1SJed Brown 
5944165533cSJose E. Roman    Output Parameters:
595a96d6ef6SBarry Smith +  tvar - a function that transforms to transient variables
596e3c11fc1SJed Brown -  ctx - a context for tvar
597e3c11fc1SJed Brown 
598e3c11fc1SJed Brown    Level: advanced
599e3c11fc1SJed Brown 
600db781477SPatrick Sanan .seealso: `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()`
601e3c11fc1SJed Brown @*/
602e3c11fc1SJed Brown PetscErrorCode DMTSGetTransientVariable(DM dm,TSTransientVariable *tvar,void *ctx)
603e3c11fc1SJed Brown {
604e3c11fc1SJed Brown   DMTS           dmts;
605e3c11fc1SJed Brown 
606e3c11fc1SJed Brown   PetscFunctionBegin;
607e3c11fc1SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6089566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&dmts));
609e3c11fc1SJed Brown   if (tvar) *tvar = dmts->ops->transientvar;
610e3c11fc1SJed Brown   if (ctx)  *(void**)ctx = dmts->transientvarctx;
611e3c11fc1SJed Brown   PetscFunctionReturn(0);
612e3c11fc1SJed Brown }
613e3c11fc1SJed Brown 
614e3c11fc1SJed Brown /*@C
615ef20d060SBarry Smith    DMTSGetSolutionFunction - gets the TS solution evaluation function
616ef20d060SBarry Smith 
617ef20d060SBarry Smith    Not Collective
618ef20d060SBarry Smith 
6194165533cSJose E. Roman    Input Parameter:
620ef20d060SBarry Smith .  dm - DM to be used with TS
621ef20d060SBarry Smith 
622ef20d060SBarry Smith    Output Parameters:
623ef20d060SBarry Smith +  func - solution function evaluation function, see TSSetSolution() for calling sequence
624ef20d060SBarry Smith -  ctx - context for solution evaluation
625ef20d060SBarry Smith 
626ef20d060SBarry Smith    Level: advanced
627ef20d060SBarry Smith 
628db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSSetSolutionFunction()`
629ef20d060SBarry Smith @*/
630ef20d060SBarry Smith PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
631ef20d060SBarry Smith {
632942e3340SBarry Smith   DMTS           tsdm;
633ef20d060SBarry Smith 
634ef20d060SBarry Smith   PetscFunctionBegin;
635ef20d060SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6369566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
637d74926cbSBarry Smith   if (func) *func = tsdm->ops->solution;
638ef20d060SBarry Smith   if (ctx)  *ctx  = tsdm->solutionctx;
639ef20d060SBarry Smith   PetscFunctionReturn(0);
640ef20d060SBarry Smith }
641ef20d060SBarry Smith 
642ef20d060SBarry Smith /*@C
643ef20d060SBarry Smith    DMTSSetSolutionFunction - set TS solution evaluation function
644ef20d060SBarry Smith 
645ef20d060SBarry Smith    Not Collective
646ef20d060SBarry Smith 
6474165533cSJose E. Roman    Input Parameters:
648ef20d060SBarry Smith +  dm - DM to be used with TS
649a96d6ef6SBarry Smith .  func - solution function evaluation routine
650ef20d060SBarry Smith -  ctx - context for solution evaluation
651ef20d060SBarry Smith 
652a96d6ef6SBarry Smith     Calling sequence of f:
653a96d6ef6SBarry Smith $     PetscErrorCode f(TS ts,PetscReal t,Vec u,void *ctx);
654a96d6ef6SBarry Smith 
655a96d6ef6SBarry Smith +   ts - timestep context
656a96d6ef6SBarry Smith .   t - current timestep
657a96d6ef6SBarry Smith .   u - output vector
658a96d6ef6SBarry Smith -   ctx - [optional] user-defined function context
659a96d6ef6SBarry Smith 
660ef20d060SBarry Smith    Level: advanced
661ef20d060SBarry Smith 
662ef20d060SBarry Smith    Note:
663ef20d060SBarry Smith    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
664ef20d060SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
665ef20d060SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
666ef20d060SBarry Smith 
667db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSGetSolutionFunction()`
668ef20d060SBarry Smith @*/
669ef20d060SBarry Smith PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
670ef20d060SBarry Smith {
671942e3340SBarry Smith   DMTS           tsdm;
672ef20d060SBarry Smith 
673ef20d060SBarry Smith   PetscFunctionBegin;
674ef20d060SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6759566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
676d74926cbSBarry Smith   if (func) tsdm->ops->solution = func;
677ef20d060SBarry Smith   if (ctx)  tsdm->solutionctx   = ctx;
678ef20d060SBarry Smith   PetscFunctionReturn(0);
679ef20d060SBarry Smith }
680ef20d060SBarry Smith 
6819b7cd975SBarry Smith /*@C
6829b7cd975SBarry Smith    DMTSSetForcingFunction - set TS forcing function evaluation function
6839b7cd975SBarry Smith 
6849b7cd975SBarry Smith    Not Collective
6859b7cd975SBarry Smith 
6864165533cSJose E. Roman    Input Parameters:
6879b7cd975SBarry Smith +  dm - DM to be used with TS
688a96d6ef6SBarry Smith .  f - forcing function evaluation routine
6899b7cd975SBarry Smith -  ctx - context for solution evaluation
6909b7cd975SBarry Smith 
691a96d6ef6SBarry Smith     Calling sequence of func:
692a96d6ef6SBarry Smith $     PetscErrorCode func (TS ts,PetscReal t,Vec f,void *ctx);
693a96d6ef6SBarry Smith 
694a96d6ef6SBarry Smith +   ts - timestep context
695a96d6ef6SBarry Smith .   t - current timestep
696a96d6ef6SBarry Smith .   f - output vector
697a96d6ef6SBarry Smith -   ctx - [optional] user-defined function context
698a96d6ef6SBarry Smith 
6999b7cd975SBarry Smith    Level: advanced
7009b7cd975SBarry Smith 
7019b7cd975SBarry Smith    Note:
7029b7cd975SBarry Smith    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
7039b7cd975SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
7049b7cd975SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
7059b7cd975SBarry Smith 
706db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
7079b7cd975SBarry Smith @*/
708d56366bfSLisandro Dalcin PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx)
7099b7cd975SBarry Smith {
7109b7cd975SBarry Smith   DMTS           tsdm;
7119b7cd975SBarry Smith 
7129b7cd975SBarry Smith   PetscFunctionBegin;
7139b7cd975SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7149566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
715f8b49ee9SBarry Smith   if (f)    tsdm->ops->forcing = f;
7169b7cd975SBarry Smith   if (ctx)  tsdm->forcingctx   = ctx;
7179b7cd975SBarry Smith   PetscFunctionReturn(0);
7189b7cd975SBarry Smith }
7199b7cd975SBarry Smith 
7209b7cd975SBarry Smith /*@C
7219b7cd975SBarry Smith    DMTSGetForcingFunction - get TS forcing function evaluation function
7229b7cd975SBarry Smith 
7239b7cd975SBarry Smith    Not Collective
7249b7cd975SBarry Smith 
7254165533cSJose E. Roman    Input Parameter:
7269b7cd975SBarry Smith .   dm - DM to be used with TS
7279b7cd975SBarry Smith 
7284165533cSJose E. Roman    Output Parameters:
729f8b49ee9SBarry Smith +  f - forcing function evaluation function; see TSForcingFunction for details
7309b7cd975SBarry Smith -  ctx - context for solution evaluation
7319b7cd975SBarry Smith 
7329b7cd975SBarry Smith    Level: advanced
7339b7cd975SBarry Smith 
7349b7cd975SBarry Smith    Note:
7359b7cd975SBarry Smith    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
7369b7cd975SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
7379b7cd975SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
7389b7cd975SBarry Smith 
739db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
7409b7cd975SBarry Smith @*/
741d56366bfSLisandro Dalcin PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
7429b7cd975SBarry Smith {
7439b7cd975SBarry Smith   DMTS           tsdm;
7449b7cd975SBarry Smith 
7459b7cd975SBarry Smith   PetscFunctionBegin;
7469b7cd975SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7479566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
748f8b49ee9SBarry Smith   if (f)   *f   = tsdm->ops->forcing;
7499b7cd975SBarry Smith   if (ctx) *ctx = tsdm->forcingctx;
7509b7cd975SBarry Smith   PetscFunctionReturn(0);
7519b7cd975SBarry Smith }
7529b7cd975SBarry Smith 
75324989b8cSPeter Brune /*@C
75424989b8cSPeter Brune    DMTSGetRHSFunction - get TS explicit residual evaluation function
75524989b8cSPeter Brune 
75624989b8cSPeter Brune    Not Collective
75724989b8cSPeter Brune 
7584165533cSJose E. Roman    Input Parameter:
75924989b8cSPeter Brune .  dm - DM to be used with TS
76024989b8cSPeter Brune 
7614165533cSJose E. Roman    Output Parameters:
76224989b8cSPeter Brune +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
76324989b8cSPeter Brune -  ctx - context for residual evaluation
76424989b8cSPeter Brune 
76524989b8cSPeter Brune    Level: advanced
76624989b8cSPeter Brune 
76724989b8cSPeter Brune    Note:
76824989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
76924989b8cSPeter Brune    associated with the DM.
77024989b8cSPeter Brune 
771db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
77224989b8cSPeter Brune @*/
77324989b8cSPeter Brune PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
77424989b8cSPeter Brune {
775942e3340SBarry Smith   DMTS           tsdm;
77624989b8cSPeter Brune 
77724989b8cSPeter Brune   PetscFunctionBegin;
77824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7799566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
780d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsfunction;
78124989b8cSPeter Brune   if (ctx)  *ctx = tsdm->rhsfunctionctx;
78224989b8cSPeter Brune   PetscFunctionReturn(0);
78324989b8cSPeter Brune }
78424989b8cSPeter Brune 
78524989b8cSPeter Brune /*@C
78624989b8cSPeter Brune    DMTSSetIJacobian - set TS Jacobian evaluation function
78724989b8cSPeter Brune 
78824989b8cSPeter Brune    Not Collective
78924989b8cSPeter Brune 
7904165533cSJose E. Roman    Input Parameters:
79124989b8cSPeter Brune +  dm - DM to be used with TS
792a96d6ef6SBarry Smith .  func - Jacobian evaluation routine
79324989b8cSPeter Brune -  ctx - context for residual evaluation
79424989b8cSPeter Brune 
795a96d6ef6SBarry Smith    Calling sequence of f:
796a96d6ef6SBarry Smith $    PetscErrorCode f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx);
797a96d6ef6SBarry Smith 
798a96d6ef6SBarry Smith +  t    - time at step/stage being solved
799a96d6ef6SBarry Smith .  U    - state vector
800a96d6ef6SBarry Smith .  U_t  - time derivative of state vector
801a96d6ef6SBarry Smith .  a    - shift
802a96d6ef6SBarry Smith .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
803a96d6ef6SBarry Smith .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
804a96d6ef6SBarry Smith -  ctx  - [optional] user-defined context for matrix evaluation routine
805a96d6ef6SBarry Smith 
80624989b8cSPeter Brune    Level: advanced
80724989b8cSPeter Brune 
80824989b8cSPeter Brune    Note:
80924989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
81024989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
81124989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
81224989b8cSPeter Brune 
813db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSGetJacobian()`, `TSSetIJacobian()`, `TSSetIFunction()`
81424989b8cSPeter Brune @*/
81524989b8cSPeter Brune PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
81624989b8cSPeter Brune {
817942e3340SBarry Smith   DMTS           sdm;
81824989b8cSPeter Brune 
81924989b8cSPeter Brune   PetscFunctionBegin;
82024989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
8219566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&sdm));
822d74926cbSBarry Smith   if (func) sdm->ops->ijacobian = func;
82324989b8cSPeter Brune   if (ctx)  sdm->ijacobianctx   = ctx;
82424989b8cSPeter Brune   PetscFunctionReturn(0);
82524989b8cSPeter Brune }
82624989b8cSPeter Brune 
82724989b8cSPeter Brune /*@C
82824989b8cSPeter Brune    DMTSGetIJacobian - get TS Jacobian evaluation function
82924989b8cSPeter Brune 
83024989b8cSPeter Brune    Not Collective
83124989b8cSPeter Brune 
8324165533cSJose E. Roman    Input Parameter:
83324989b8cSPeter Brune .  dm - DM to be used with TS
83424989b8cSPeter Brune 
8354165533cSJose E. Roman    Output Parameters:
83624989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
83724989b8cSPeter Brune -  ctx - context for residual evaluation
83824989b8cSPeter Brune 
83924989b8cSPeter Brune    Level: advanced
84024989b8cSPeter Brune 
84124989b8cSPeter Brune    Note:
84224989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
84324989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
84424989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
84524989b8cSPeter Brune 
846db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
84724989b8cSPeter Brune @*/
84824989b8cSPeter Brune PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
84924989b8cSPeter Brune {
850942e3340SBarry Smith   DMTS           tsdm;
85124989b8cSPeter Brune 
85224989b8cSPeter Brune   PetscFunctionBegin;
85324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
8549566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
855d74926cbSBarry Smith   if (func) *func = tsdm->ops->ijacobian;
85624989b8cSPeter Brune   if (ctx)  *ctx = tsdm->ijacobianctx;
85724989b8cSPeter Brune   PetscFunctionReturn(0);
85824989b8cSPeter Brune }
85924989b8cSPeter Brune 
86024989b8cSPeter Brune /*@C
86124989b8cSPeter Brune    DMTSSetRHSJacobian - set TS Jacobian evaluation function
86224989b8cSPeter Brune 
86324989b8cSPeter Brune    Not Collective
86424989b8cSPeter Brune 
8654165533cSJose E. Roman    Input Parameters:
86624989b8cSPeter Brune +  dm - DM to be used with TS
867a96d6ef6SBarry Smith .  func - Jacobian evaluation routine
86824989b8cSPeter Brune -  ctx - context for residual evaluation
86924989b8cSPeter Brune 
870a96d6ef6SBarry Smith    Calling sequence of func:
871a96d6ef6SBarry Smith $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);
872a96d6ef6SBarry Smith 
873a96d6ef6SBarry Smith +  t - current timestep
874a96d6ef6SBarry Smith .  u - input vector
875a96d6ef6SBarry Smith .  Amat - (approximate) Jacobian matrix
876a96d6ef6SBarry Smith .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
877a96d6ef6SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
878a96d6ef6SBarry Smith 
87924989b8cSPeter Brune    Level: advanced
88024989b8cSPeter Brune 
88124989b8cSPeter Brune    Note:
88224989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
88324989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
88424989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
88524989b8cSPeter Brune 
886db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSGetJacobian()`, `TSSetRHSJacobian()`
88724989b8cSPeter Brune @*/
88824989b8cSPeter Brune PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
88924989b8cSPeter Brune {
890942e3340SBarry Smith   DMTS           tsdm;
89124989b8cSPeter Brune 
89224989b8cSPeter Brune   PetscFunctionBegin;
89324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
8949566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
895d74926cbSBarry Smith   if (func) tsdm->ops->rhsjacobian = func;
89624989b8cSPeter Brune   if (ctx)  tsdm->rhsjacobianctx = ctx;
89724989b8cSPeter Brune   PetscFunctionReturn(0);
89824989b8cSPeter Brune }
89924989b8cSPeter Brune 
90024989b8cSPeter Brune /*@C
90124989b8cSPeter Brune    DMTSGetRHSJacobian - get TS Jacobian evaluation function
90224989b8cSPeter Brune 
90324989b8cSPeter Brune    Not Collective
90424989b8cSPeter Brune 
9054165533cSJose E. Roman    Input Parameter:
90624989b8cSPeter Brune .  dm - DM to be used with TS
90724989b8cSPeter Brune 
9084165533cSJose E. Roman    Output Parameters:
90924989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
91024989b8cSPeter Brune -  ctx - context for residual evaluation
91124989b8cSPeter Brune 
91224989b8cSPeter Brune    Level: advanced
91324989b8cSPeter Brune 
91424989b8cSPeter Brune    Note:
91524989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
91624989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
91724989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
91824989b8cSPeter Brune 
919db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
92024989b8cSPeter Brune @*/
92124989b8cSPeter Brune PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
92224989b8cSPeter Brune {
923942e3340SBarry Smith   DMTS           tsdm;
92424989b8cSPeter Brune 
92524989b8cSPeter Brune   PetscFunctionBegin;
92624989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
9279566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm,&tsdm));
928d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsjacobian;
92924989b8cSPeter Brune   if (ctx)  *ctx = tsdm->rhsjacobianctx;
93024989b8cSPeter Brune   PetscFunctionReturn(0);
93124989b8cSPeter Brune }
932ad6bc421SBarry Smith 
933ad6bc421SBarry Smith /*@C
934ad6bc421SBarry Smith    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
935ad6bc421SBarry Smith 
936ad6bc421SBarry Smith    Not Collective
937ad6bc421SBarry Smith 
9384165533cSJose E. Roman    Input Parameters:
939ad6bc421SBarry Smith +  dm - DM to be used with TS
940ad6bc421SBarry Smith .  view - viewer function
941ad6bc421SBarry Smith -  load - loading function
942ad6bc421SBarry Smith 
943ad6bc421SBarry Smith    Level: advanced
944ad6bc421SBarry Smith 
945db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
946ad6bc421SBarry Smith @*/
947ad6bc421SBarry Smith PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
948ad6bc421SBarry Smith {
949ad6bc421SBarry Smith   DMTS           tsdm;
950ad6bc421SBarry Smith 
951ad6bc421SBarry Smith   PetscFunctionBegin;
952ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
9539566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
954ad6bc421SBarry Smith   tsdm->ops->ifunctionview = view;
955ad6bc421SBarry Smith   tsdm->ops->ifunctionload = load;
956ad6bc421SBarry Smith   PetscFunctionReturn(0);
957ad6bc421SBarry Smith }
958ad6bc421SBarry Smith 
959ad6bc421SBarry Smith /*@C
960ad6bc421SBarry Smith    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
961ad6bc421SBarry Smith 
962ad6bc421SBarry Smith    Not Collective
963ad6bc421SBarry Smith 
9644165533cSJose E. Roman    Input Parameters:
965ad6bc421SBarry Smith +  dm - DM to be used with TS
966ad6bc421SBarry Smith .  view - viewer function
967ad6bc421SBarry Smith -  load - loading function
968ad6bc421SBarry Smith 
969ad6bc421SBarry Smith    Level: advanced
970ad6bc421SBarry Smith 
971db781477SPatrick Sanan .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`
972ad6bc421SBarry Smith @*/
973ad6bc421SBarry Smith PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
974ad6bc421SBarry Smith {
975ad6bc421SBarry Smith   DMTS           tsdm;
976ad6bc421SBarry Smith 
977ad6bc421SBarry Smith   PetscFunctionBegin;
978ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
9799566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm,&tsdm));
980ad6bc421SBarry Smith   tsdm->ops->ijacobianview = view;
981ad6bc421SBarry Smith   tsdm->ops->ijacobianload = load;
982ad6bc421SBarry Smith   PetscFunctionReturn(0);
983ad6bc421SBarry Smith }
984