xref: /petsc/src/ts/utils/dmts.c (revision efe9872e0e12ef39052464c4f0551d2124a04b33)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
2af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
324989b8cSPeter Brune 
4d74926cbSBarry Smith #undef __FUNCT__
5d74926cbSBarry Smith #define __FUNCT__ "DMTSDestroy"
6d74926cbSBarry Smith static PetscErrorCode DMTSDestroy(DMTS *kdm)
7d74926cbSBarry Smith {
8d74926cbSBarry Smith   PetscErrorCode ierr;
9d74926cbSBarry Smith 
10d74926cbSBarry Smith   PetscFunctionBegin;
11d74926cbSBarry Smith   if (!*kdm) PetscFunctionReturn(0);
12d74926cbSBarry Smith   PetscValidHeaderSpecific((*kdm),DMTS_CLASSID,1);
13d74926cbSBarry Smith   if (--((PetscObject)(*kdm))->refct > 0) {*kdm = 0; PetscFunctionReturn(0);}
14d74926cbSBarry Smith   if ((*kdm)->ops->destroy) {ierr = ((*kdm)->ops->destroy)(*kdm);CHKERRQ(ierr);}
15d74926cbSBarry Smith   ierr = PetscHeaderDestroy(kdm);CHKERRQ(ierr);
16d74926cbSBarry Smith   PetscFunctionReturn(0);
17d74926cbSBarry Smith }
18d74926cbSBarry Smith 
19d74926cbSBarry Smith #undef __FUNCT__
202d53ad75SBarry Smith #define __FUNCT__ "DMTSLoad"
212d53ad75SBarry Smith PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer)
222d53ad75SBarry Smith {
232d53ad75SBarry Smith   PetscErrorCode ierr;
242d53ad75SBarry Smith 
252d53ad75SBarry Smith   PetscFunctionBegin;
26060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
27060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
28060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
29ad6bc421SBarry Smith   if (kdm->ops->ifunctionload) {
30ad6bc421SBarry Smith     ierr = (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);CHKERRQ(ierr);
31ad6bc421SBarry Smith   }
32060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
33060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
34060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
35ad6bc421SBarry Smith   if (kdm->ops->ijacobianload) {
36ad6bc421SBarry Smith     ierr = (*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer);CHKERRQ(ierr);
37ad6bc421SBarry Smith   }
382d53ad75SBarry Smith   PetscFunctionReturn(0);
392d53ad75SBarry Smith }
402d53ad75SBarry Smith 
412d53ad75SBarry Smith #undef __FUNCT__
422d53ad75SBarry Smith #define __FUNCT__ "DMTSView"
432d53ad75SBarry Smith PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer)
442d53ad75SBarry Smith {
452d53ad75SBarry Smith   PetscErrorCode ierr;
462d53ad75SBarry Smith   PetscBool      isascii,isbinary;
472d53ad75SBarry Smith 
482d53ad75SBarry Smith   PetscFunctionBegin;
492d53ad75SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
502d53ad75SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
512d53ad75SBarry Smith   if (isascii) {
52c7a10e08SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS)
532d53ad75SBarry Smith     const char *fname;
542d53ad75SBarry Smith 
552d53ad75SBarry Smith     ierr = PetscFPTFind(kdm->ops->ifunction,&fname);CHKERRQ(ierr);
562d53ad75SBarry Smith     if (fname) {
572d53ad75SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  IFunction used by TS: %s\n",fname);CHKERRQ(ierr);
582d53ad75SBarry Smith     }
592d53ad75SBarry Smith     ierr = PetscFPTFind(kdm->ops->ijacobian,&fname);CHKERRQ(ierr);
602d53ad75SBarry Smith     if (fname) {
612d53ad75SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  IJacobian function used by TS: %s\n",fname);CHKERRQ(ierr);
622d53ad75SBarry Smith     }
63c7a10e08SBarry Smith #endif
642d53ad75SBarry Smith   } else if (isbinary) {
653964eb88SJed Brown     struct {
663964eb88SJed Brown       TSIFunction ifunction;
679200755eSBarry Smith     } funcstruct;
689200755eSBarry Smith     struct {
693964eb88SJed Brown       PetscErrorCode (*ifunctionview)(void*,PetscViewer);
709200755eSBarry Smith     } funcviewstruct;
719200755eSBarry Smith     struct {
723964eb88SJed Brown       PetscErrorCode (*ifunctionload)(void**,PetscViewer);
739200755eSBarry Smith     } funcloadstruct;
743964eb88SJed Brown     struct {
753964eb88SJed Brown       TSIJacobian ijacobian;
769200755eSBarry Smith     } jacstruct;
779200755eSBarry Smith     struct {
783964eb88SJed Brown       PetscErrorCode (*ijacobianview)(void*,PetscViewer);
799200755eSBarry Smith     } jacviewstruct;
809200755eSBarry Smith     struct {
813964eb88SJed Brown       PetscErrorCode (*ijacobianload)(void**,PetscViewer);
829200755eSBarry Smith     } jacloadstruct;
833964eb88SJed Brown 
849200755eSBarry Smith     funcstruct.ifunction         = kdm->ops->ifunction;
859200755eSBarry Smith     funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
869200755eSBarry Smith     funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
879200755eSBarry Smith     ierr = PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
889200755eSBarry Smith     ierr = PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
899200755eSBarry Smith     ierr = PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
90ad6bc421SBarry Smith     if (kdm->ops->ifunctionview) {
91ad6bc421SBarry Smith       ierr = (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);CHKERRQ(ierr);
92ad6bc421SBarry Smith     }
939200755eSBarry Smith     jacstruct.ijacobian = kdm->ops->ijacobian;
949200755eSBarry Smith     jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
959200755eSBarry Smith     jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
969200755eSBarry Smith     ierr = PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
979200755eSBarry Smith     ierr = PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
989200755eSBarry Smith     ierr = PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
99ad6bc421SBarry Smith     if (kdm->ops->ijacobianview) {
100ad6bc421SBarry Smith       ierr = (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);CHKERRQ(ierr);
101ad6bc421SBarry Smith     }
1022d53ad75SBarry Smith   }
1032d53ad75SBarry Smith   PetscFunctionReturn(0);
1042d53ad75SBarry Smith }
1052d53ad75SBarry Smith 
1062d53ad75SBarry Smith #undef __FUNCT__
107d74926cbSBarry Smith #define __FUNCT__ "DMTSCreate"
108d74926cbSBarry Smith static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
109d74926cbSBarry Smith {
110d74926cbSBarry Smith   PetscErrorCode ierr;
111d74926cbSBarry Smith 
112d74926cbSBarry Smith   PetscFunctionBegin;
113607a6623SBarry Smith   ierr = TSInitializePackage();CHKERRQ(ierr);
11473107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);CHKERRQ(ierr);
115d74926cbSBarry Smith   PetscFunctionReturn(0);
116d74926cbSBarry Smith }
11724989b8cSPeter Brune 
11824989b8cSPeter Brune #undef __FUNCT__
119942e3340SBarry Smith #define __FUNCT__ "DMCoarsenHook_DMTS"
1202a34c10cSBarry Smith /* Attaches the DMTS to the coarse level.
12124989b8cSPeter Brune  * Under what conditions should we copy versus duplicate?
12224989b8cSPeter Brune  */
123942e3340SBarry Smith static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
12424989b8cSPeter Brune {
12524989b8cSPeter Brune   PetscErrorCode ierr;
12624989b8cSPeter Brune 
12724989b8cSPeter Brune   PetscFunctionBegin;
128942e3340SBarry Smith   ierr = DMCopyDMTS(dm,dmc);CHKERRQ(ierr);
12924989b8cSPeter Brune   PetscFunctionReturn(0);
13024989b8cSPeter Brune }
13124989b8cSPeter Brune 
13224989b8cSPeter Brune #undef __FUNCT__
133942e3340SBarry Smith #define __FUNCT__ "DMRestrictHook_DMTS"
13424989b8cSPeter Brune /* This could restrict auxiliary information to the coarse level.
13524989b8cSPeter Brune  */
136942e3340SBarry Smith static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
13724989b8cSPeter Brune {
13824989b8cSPeter Brune 
13924989b8cSPeter Brune   PetscFunctionBegin;
14024989b8cSPeter Brune   PetscFunctionReturn(0);
14124989b8cSPeter Brune }
14224989b8cSPeter Brune 
14324989b8cSPeter Brune #undef __FUNCT__
144258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_DMTS"
145258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
146258e1594SPeter Brune {
147258e1594SPeter Brune   PetscErrorCode ierr;
148258e1594SPeter Brune 
149258e1594SPeter Brune   PetscFunctionBegin;
150258e1594SPeter Brune   ierr = DMCopyDMTS(dm,subdm);CHKERRQ(ierr);
151258e1594SPeter Brune   PetscFunctionReturn(0);
152258e1594SPeter Brune }
153258e1594SPeter Brune 
154258e1594SPeter Brune #undef __FUNCT__
155258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_DMTS"
156258e1594SPeter Brune /* This could restrict auxiliary information to the coarse level.
157258e1594SPeter Brune  */
158258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
159258e1594SPeter Brune {
160258e1594SPeter Brune   PetscFunctionBegin;
161258e1594SPeter Brune   PetscFunctionReturn(0);
162258e1594SPeter Brune }
163258e1594SPeter Brune 
164258e1594SPeter Brune #undef __FUNCT__
165d74926cbSBarry Smith #define __FUNCT__ "DMTSCopy"
166d74926cbSBarry Smith /*@C
167d74926cbSBarry Smith    DMTSCopy - copies the information in a DMTS to another DMTS
168d74926cbSBarry Smith 
169d74926cbSBarry Smith    Not Collective
170d74926cbSBarry Smith 
171d74926cbSBarry Smith    Input Argument:
172d74926cbSBarry Smith +  kdm - Original DMTS
173d74926cbSBarry Smith -  nkdm - DMTS to receive the data, should have been created with DMTSCreate()
174d74926cbSBarry Smith 
175d74926cbSBarry Smith    Level: developer
176d74926cbSBarry Smith 
177d74926cbSBarry Smith .seealso: DMTSCreate(), DMTSDestroy()
178d74926cbSBarry Smith @*/
179d74926cbSBarry Smith PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
18024989b8cSPeter Brune {
18124989b8cSPeter Brune   PetscErrorCode ierr;
18224989b8cSPeter Brune 
18324989b8cSPeter Brune   PetscFunctionBegin;
184d74926cbSBarry Smith   PetscValidHeaderSpecific(kdm,DMTS_CLASSID,1);
185d74926cbSBarry Smith   PetscValidHeaderSpecific(nkdm,DMTS_CLASSID,2);
186d74926cbSBarry Smith   nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
187d74926cbSBarry Smith   nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
188d74926cbSBarry Smith   nkdm->ops->ifunction   = kdm->ops->ifunction;
189d74926cbSBarry Smith   nkdm->ops->ijacobian   = kdm->ops->ijacobian;
190*efe9872eSLisandro Dalcin   nkdm->ops->i2function  = kdm->ops->i2function;
191*efe9872eSLisandro Dalcin   nkdm->ops->i2jacobian  = kdm->ops->i2jacobian;
192d74926cbSBarry Smith   nkdm->ops->solution    = kdm->ops->solution;
193d74926cbSBarry Smith   nkdm->ops->destroy     = kdm->ops->destroy;
194d74926cbSBarry Smith   nkdm->ops->duplicate   = kdm->ops->duplicate;
195d74926cbSBarry Smith 
196d74926cbSBarry Smith   nkdm->rhsfunctionctx = kdm->rhsfunctionctx;
197d74926cbSBarry Smith   nkdm->rhsjacobianctx = kdm->rhsjacobianctx;
198d74926cbSBarry Smith   nkdm->ifunctionctx   = kdm->ifunctionctx;
199d74926cbSBarry Smith   nkdm->ijacobianctx   = kdm->ijacobianctx;
200*efe9872eSLisandro Dalcin   nkdm->i2functionctx  = kdm->i2functionctx;
201*efe9872eSLisandro Dalcin   nkdm->i2jacobianctx  = kdm->i2jacobianctx;
202d74926cbSBarry Smith   nkdm->solutionctx    = kdm->solutionctx;
203d74926cbSBarry Smith 
204d74926cbSBarry Smith   nkdm->data = kdm->data;
205d74926cbSBarry Smith 
206d74926cbSBarry Smith   /*
207d74926cbSBarry Smith   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
208d74926cbSBarry Smith   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
209d74926cbSBarry Smith   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
210d74926cbSBarry Smith   */
211d74926cbSBarry Smith 
212d74926cbSBarry Smith   /* implementation specific copy hooks */
213d74926cbSBarry Smith   if (kdm->ops->duplicate) {ierr = (*kdm->ops->duplicate)(kdm,nkdm);CHKERRQ(ierr);}
21424989b8cSPeter Brune   PetscFunctionReturn(0);
21524989b8cSPeter Brune }
21624989b8cSPeter Brune 
21724989b8cSPeter Brune #undef __FUNCT__
218942e3340SBarry Smith #define __FUNCT__ "DMGetDMTS"
21924989b8cSPeter Brune /*@C
220942e3340SBarry Smith    DMGetDMTS - get read-only private DMTS context from a DM
22124989b8cSPeter Brune 
22224989b8cSPeter Brune    Not Collective
22324989b8cSPeter Brune 
22424989b8cSPeter Brune    Input Argument:
22524989b8cSPeter Brune .  dm - DM to be used with TS
22624989b8cSPeter Brune 
22724989b8cSPeter Brune    Output Argument:
228942e3340SBarry Smith .  tsdm - private DMTS context
22924989b8cSPeter Brune 
23024989b8cSPeter Brune    Level: developer
23124989b8cSPeter Brune 
23224989b8cSPeter Brune    Notes:
233942e3340SBarry Smith    Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.
23424989b8cSPeter Brune 
235942e3340SBarry Smith .seealso: DMGetDMTSWrite()
23624989b8cSPeter Brune @*/
237942e3340SBarry Smith PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
23824989b8cSPeter Brune {
23924989b8cSPeter Brune   PetscErrorCode ierr;
24024989b8cSPeter Brune 
24124989b8cSPeter Brune   PetscFunctionBegin;
24224989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2432a34c10cSBarry Smith   *tsdm = (DMTS) dm->dmts;
244d74926cbSBarry Smith   if (!*tsdm) {
245942e3340SBarry Smith     ierr = PetscInfo(dm,"Creating new DMTS\n");CHKERRQ(ierr);
246ce94432eSBarry Smith     ierr = DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm);CHKERRQ(ierr);
2472a34c10cSBarry Smith     dm->dmts = (PetscObject) *tsdm;
2480298fd71SBarry Smith     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr);
2490298fd71SBarry Smith     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr);
25024989b8cSPeter Brune   }
25124989b8cSPeter Brune   PetscFunctionReturn(0);
25224989b8cSPeter Brune }
25324989b8cSPeter Brune 
25424989b8cSPeter Brune #undef __FUNCT__
255942e3340SBarry Smith #define __FUNCT__ "DMGetDMTSWrite"
25624989b8cSPeter Brune /*@C
257942e3340SBarry Smith    DMGetDMTSWrite - get write access to private DMTS context from a DM
25824989b8cSPeter Brune 
25924989b8cSPeter Brune    Not Collective
26024989b8cSPeter Brune 
26124989b8cSPeter Brune    Input Argument:
26224989b8cSPeter Brune .  dm - DM to be used with TS
26324989b8cSPeter Brune 
26424989b8cSPeter Brune    Output Argument:
265942e3340SBarry Smith .  tsdm - private DMTS context
26624989b8cSPeter Brune 
26724989b8cSPeter Brune    Level: developer
26824989b8cSPeter Brune 
269942e3340SBarry Smith .seealso: DMGetDMTS()
27024989b8cSPeter Brune @*/
271942e3340SBarry Smith PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
27224989b8cSPeter Brune {
27324989b8cSPeter Brune   PetscErrorCode ierr;
274942e3340SBarry Smith   DMTS           sdm;
27524989b8cSPeter Brune 
27624989b8cSPeter Brune   PetscFunctionBegin;
27724989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
278942e3340SBarry Smith   ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr);
27924989b8cSPeter Brune   if (!sdm->originaldm) sdm->originaldm = dm;
28024989b8cSPeter Brune   if (sdm->originaldm != dm) {  /* Copy on write */
2812a34c10cSBarry Smith     DMTS oldsdm = sdm;
282942e3340SBarry Smith     ierr     = PetscInfo(dm,"Copying DMTS due to write\n");CHKERRQ(ierr);
283ce94432eSBarry Smith     ierr     = DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);CHKERRQ(ierr);
284d74926cbSBarry Smith     ierr     = DMTSCopy(oldsdm,sdm);CHKERRQ(ierr);
2852a34c10cSBarry Smith     ierr     = DMTSDestroy((DMTS*)&dm->dmts);CHKERRQ(ierr);
2862a34c10cSBarry Smith     dm->dmts = (PetscObject) sdm;
28724989b8cSPeter Brune   }
28824989b8cSPeter Brune   *tsdm = sdm;
28924989b8cSPeter Brune   PetscFunctionReturn(0);
29024989b8cSPeter Brune }
29124989b8cSPeter Brune 
29224989b8cSPeter Brune #undef __FUNCT__
293942e3340SBarry Smith #define __FUNCT__ "DMCopyDMTS"
29424989b8cSPeter Brune /*@C
295942e3340SBarry Smith    DMCopyDMTS - copies a DM context to a new DM
29624989b8cSPeter Brune 
29724989b8cSPeter Brune    Logically Collective
29824989b8cSPeter Brune 
29924989b8cSPeter Brune    Input Arguments:
30024989b8cSPeter Brune +  dmsrc - DM to obtain context from
30124989b8cSPeter Brune -  dmdest - DM to add context to
30224989b8cSPeter Brune 
30324989b8cSPeter Brune    Level: developer
30424989b8cSPeter Brune 
30524989b8cSPeter Brune    Note:
30624989b8cSPeter Brune    The context is copied by reference. This function does not ensure that a context exists.
30724989b8cSPeter Brune 
308942e3340SBarry Smith .seealso: DMGetDMTS(), TSSetDM()
30924989b8cSPeter Brune @*/
310942e3340SBarry Smith PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
31124989b8cSPeter Brune {
31224989b8cSPeter Brune   PetscErrorCode ierr;
31324989b8cSPeter Brune 
31424989b8cSPeter Brune   PetscFunctionBegin;
31524989b8cSPeter Brune   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
31624989b8cSPeter Brune   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
3172a34c10cSBarry Smith   ierr         = DMTSDestroy((DMTS*)&dmdest->dmts);CHKERRQ(ierr);
3182a34c10cSBarry Smith   dmdest->dmts = dmsrc->dmts;
3192a34c10cSBarry Smith   ierr         = PetscObjectReference(dmdest->dmts);CHKERRQ(ierr);
3200298fd71SBarry Smith   ierr         = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr);
3210298fd71SBarry Smith   ierr         = DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr);
32224989b8cSPeter Brune   PetscFunctionReturn(0);
32324989b8cSPeter Brune }
32424989b8cSPeter Brune 
32524989b8cSPeter Brune #undef __FUNCT__
32624989b8cSPeter Brune #define __FUNCT__ "DMTSSetIFunction"
32724989b8cSPeter Brune /*@C
32824989b8cSPeter Brune    DMTSSetIFunction - set TS implicit function evaluation function
32924989b8cSPeter Brune 
33024989b8cSPeter Brune    Not Collective
33124989b8cSPeter Brune 
33224989b8cSPeter Brune    Input Arguments:
33324989b8cSPeter Brune +  dm - DM to be used with TS
33424989b8cSPeter Brune .  func - function evaluation function, see TSSetIFunction() for calling sequence
33524989b8cSPeter Brune -  ctx - context for residual evaluation
33624989b8cSPeter Brune 
33724989b8cSPeter Brune    Level: advanced
33824989b8cSPeter Brune 
33924989b8cSPeter Brune    Note:
34024989b8cSPeter Brune    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
34124989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
34224989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
34324989b8cSPeter Brune 
34424989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
34524989b8cSPeter Brune @*/
34624989b8cSPeter Brune PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
34724989b8cSPeter Brune {
34824989b8cSPeter Brune   PetscErrorCode ierr;
349942e3340SBarry Smith   DMTS           tsdm;
35024989b8cSPeter Brune 
35124989b8cSPeter Brune   PetscFunctionBegin;
35224989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
353942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
354d74926cbSBarry Smith   if (func) tsdm->ops->ifunction = func;
35524989b8cSPeter Brune   if (ctx)  tsdm->ifunctionctx = ctx;
35624989b8cSPeter Brune   PetscFunctionReturn(0);
35724989b8cSPeter Brune }
35824989b8cSPeter Brune 
35924989b8cSPeter Brune #undef __FUNCT__
36024989b8cSPeter Brune #define __FUNCT__ "DMTSGetIFunction"
36124989b8cSPeter Brune /*@C
36224989b8cSPeter Brune    DMTSGetIFunction - get TS implicit residual evaluation function
36324989b8cSPeter Brune 
36424989b8cSPeter Brune    Not Collective
36524989b8cSPeter Brune 
36624989b8cSPeter Brune    Input Argument:
36724989b8cSPeter Brune .  dm - DM to be used with TS
36824989b8cSPeter Brune 
36924989b8cSPeter Brune    Output Arguments:
37024989b8cSPeter Brune +  func - function evaluation function, see TSSetIFunction() for calling sequence
37124989b8cSPeter Brune -  ctx - context for residual evaluation
37224989b8cSPeter Brune 
37324989b8cSPeter Brune    Level: advanced
37424989b8cSPeter Brune 
37524989b8cSPeter Brune    Note:
37624989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
37724989b8cSPeter Brune    associated with the DM.
37824989b8cSPeter Brune 
37924989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
38024989b8cSPeter Brune @*/
38124989b8cSPeter Brune PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
38224989b8cSPeter Brune {
38324989b8cSPeter Brune   PetscErrorCode ierr;
384942e3340SBarry Smith   DMTS           tsdm;
38524989b8cSPeter Brune 
38624989b8cSPeter Brune   PetscFunctionBegin;
38724989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
388942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
389d74926cbSBarry Smith   if (func) *func = tsdm->ops->ifunction;
39024989b8cSPeter Brune   if (ctx)  *ctx = tsdm->ifunctionctx;
39124989b8cSPeter Brune   PetscFunctionReturn(0);
39224989b8cSPeter Brune }
39324989b8cSPeter Brune 
394*efe9872eSLisandro Dalcin #undef __FUNCT__
395*efe9872eSLisandro Dalcin #define __FUNCT__ "DMTSSetI2Function"
396*efe9872eSLisandro Dalcin /*@C
397*efe9872eSLisandro Dalcin    DMTSSetI2Function - set TS implicit function evaluation function for 2nd order systems
398*efe9872eSLisandro Dalcin 
399*efe9872eSLisandro Dalcin    Not Collective
400*efe9872eSLisandro Dalcin 
401*efe9872eSLisandro Dalcin    Input Arguments:
402*efe9872eSLisandro Dalcin +  dm - DM to be used with TS
403*efe9872eSLisandro Dalcin .  fun - function evaluation function, see TSSetI2Function() for calling sequence
404*efe9872eSLisandro Dalcin -  ctx - context for residual evaluation
405*efe9872eSLisandro Dalcin 
406*efe9872eSLisandro Dalcin    Level: advanced
407*efe9872eSLisandro Dalcin 
408*efe9872eSLisandro Dalcin    Note:
409*efe9872eSLisandro Dalcin    TSSetI2Function() is normally used, but it calls this function internally because the user context is actually
410*efe9872eSLisandro Dalcin    associated with the DM.
411*efe9872eSLisandro Dalcin 
412*efe9872eSLisandro Dalcin .seealso: TSSetI2Function()
413*efe9872eSLisandro Dalcin @*/
414*efe9872eSLisandro Dalcin PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx)
415*efe9872eSLisandro Dalcin {
416*efe9872eSLisandro Dalcin   DMTS           tsdm;
417*efe9872eSLisandro Dalcin   PetscErrorCode ierr;
418*efe9872eSLisandro Dalcin 
419*efe9872eSLisandro Dalcin   PetscFunctionBegin;
420*efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
421*efe9872eSLisandro Dalcin   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
422*efe9872eSLisandro Dalcin   if (fun) tsdm->ops->i2function = fun;
423*efe9872eSLisandro Dalcin   if (ctx) tsdm->i2functionctx   = ctx;
424*efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
425*efe9872eSLisandro Dalcin }
426*efe9872eSLisandro Dalcin 
427*efe9872eSLisandro Dalcin #undef __FUNCT__
428*efe9872eSLisandro Dalcin #define __FUNCT__ "DMTSGetI2Function"
429*efe9872eSLisandro Dalcin /*@C
430*efe9872eSLisandro Dalcin    DMTSGetI2Function - get TS implicit residual evaluation function for 2nd order systems
431*efe9872eSLisandro Dalcin 
432*efe9872eSLisandro Dalcin    Not Collective
433*efe9872eSLisandro Dalcin 
434*efe9872eSLisandro Dalcin    Input Argument:
435*efe9872eSLisandro Dalcin .  dm - DM to be used with TS
436*efe9872eSLisandro Dalcin 
437*efe9872eSLisandro Dalcin    Output Arguments:
438*efe9872eSLisandro Dalcin +  fun - function evaluation function, see TSSetI2Function() for calling sequence
439*efe9872eSLisandro Dalcin -  ctx - context for residual evaluation
440*efe9872eSLisandro Dalcin 
441*efe9872eSLisandro Dalcin    Level: advanced
442*efe9872eSLisandro Dalcin 
443*efe9872eSLisandro Dalcin    Note:
444*efe9872eSLisandro Dalcin    TSGetI2Function() is normally used, but it calls this function internally because the user context is actually
445*efe9872eSLisandro Dalcin    associated with the DM.
446*efe9872eSLisandro Dalcin 
447*efe9872eSLisandro Dalcin .seealso: DMTSSetI2Function(),TSGetI2Function()
448*efe9872eSLisandro Dalcin @*/
449*efe9872eSLisandro Dalcin PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx)
450*efe9872eSLisandro Dalcin {
451*efe9872eSLisandro Dalcin   DMTS           tsdm;
452*efe9872eSLisandro Dalcin   PetscErrorCode ierr;
453*efe9872eSLisandro Dalcin 
454*efe9872eSLisandro Dalcin   PetscFunctionBegin;
455*efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
456*efe9872eSLisandro Dalcin   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
457*efe9872eSLisandro Dalcin   if (fun) *fun = tsdm->ops->i2function;
458*efe9872eSLisandro Dalcin   if (ctx) *ctx = tsdm->i2functionctx;
459*efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
460*efe9872eSLisandro Dalcin }
461*efe9872eSLisandro Dalcin 
462*efe9872eSLisandro Dalcin #undef __FUNCT__
463*efe9872eSLisandro Dalcin #define __FUNCT__ "DMTSSetI2Jacobian"
464*efe9872eSLisandro Dalcin /*@C
465*efe9872eSLisandro Dalcin    DMTSSetI2Jacobian - set TS implicit Jacobian evaluation function for 2nd order systems
466*efe9872eSLisandro Dalcin 
467*efe9872eSLisandro Dalcin    Not Collective
468*efe9872eSLisandro Dalcin 
469*efe9872eSLisandro Dalcin    Input Arguments:
470*efe9872eSLisandro Dalcin +  dm - DM to be used with TS
471*efe9872eSLisandro Dalcin .  fun - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
472*efe9872eSLisandro Dalcin -  ctx - context for Jacobian evaluation
473*efe9872eSLisandro Dalcin 
474*efe9872eSLisandro Dalcin    Level: advanced
475*efe9872eSLisandro Dalcin 
476*efe9872eSLisandro Dalcin    Note:
477*efe9872eSLisandro Dalcin    TSSetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
478*efe9872eSLisandro Dalcin    associated with the DM.
479*efe9872eSLisandro Dalcin 
480*efe9872eSLisandro Dalcin .seealso: TSSetI2Jacobian()
481*efe9872eSLisandro Dalcin @*/
482*efe9872eSLisandro Dalcin PetscErrorCode DMTSSetI2Jacobian(DM dm,TSI2Jacobian jac,void *ctx)
483*efe9872eSLisandro Dalcin {
484*efe9872eSLisandro Dalcin   DMTS           tsdm;
485*efe9872eSLisandro Dalcin   PetscErrorCode ierr;
486*efe9872eSLisandro Dalcin 
487*efe9872eSLisandro Dalcin   PetscFunctionBegin;
488*efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
489*efe9872eSLisandro Dalcin   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
490*efe9872eSLisandro Dalcin   if (jac) tsdm->ops->i2jacobian = jac;
491*efe9872eSLisandro Dalcin   if (ctx) tsdm->i2jacobianctx   = ctx;
492*efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
493*efe9872eSLisandro Dalcin }
494*efe9872eSLisandro Dalcin 
495*efe9872eSLisandro Dalcin #undef __FUNCT__
496*efe9872eSLisandro Dalcin #define __FUNCT__ "DMTSGetI2Jacobian"
497*efe9872eSLisandro Dalcin /*@C
498*efe9872eSLisandro Dalcin    DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems
499*efe9872eSLisandro Dalcin 
500*efe9872eSLisandro Dalcin    Not Collective
501*efe9872eSLisandro Dalcin 
502*efe9872eSLisandro Dalcin    Input Argument:
503*efe9872eSLisandro Dalcin .  dm - DM to be used with TS
504*efe9872eSLisandro Dalcin 
505*efe9872eSLisandro Dalcin    Output Arguments:
506*efe9872eSLisandro Dalcin +  jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
507*efe9872eSLisandro Dalcin -  ctx - context for Jacobian evaluation
508*efe9872eSLisandro Dalcin 
509*efe9872eSLisandro Dalcin    Level: advanced
510*efe9872eSLisandro Dalcin 
511*efe9872eSLisandro Dalcin    Note:
512*efe9872eSLisandro Dalcin    TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
513*efe9872eSLisandro Dalcin    associated with the DM.
514*efe9872eSLisandro Dalcin 
515*efe9872eSLisandro Dalcin .seealso: DMTSSetI2Jacobian(),TSGetI2Jacobian()
516*efe9872eSLisandro Dalcin @*/
517*efe9872eSLisandro Dalcin PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx)
518*efe9872eSLisandro Dalcin {
519*efe9872eSLisandro Dalcin   DMTS           tsdm;
520*efe9872eSLisandro Dalcin   PetscErrorCode ierr;
521*efe9872eSLisandro Dalcin 
522*efe9872eSLisandro Dalcin   PetscFunctionBegin;
523*efe9872eSLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
524*efe9872eSLisandro Dalcin   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
525*efe9872eSLisandro Dalcin   if (jac) *jac = tsdm->ops->i2jacobian;
526*efe9872eSLisandro Dalcin   if (ctx) *ctx = tsdm->i2jacobianctx;
527*efe9872eSLisandro Dalcin   PetscFunctionReturn(0);
528*efe9872eSLisandro Dalcin }
52924989b8cSPeter Brune 
53024989b8cSPeter Brune #undef __FUNCT__
53124989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSFunction"
53224989b8cSPeter Brune /*@C
53324989b8cSPeter Brune    DMTSSetRHSFunction - set TS explicit residual evaluation function
53424989b8cSPeter Brune 
53524989b8cSPeter Brune    Not Collective
53624989b8cSPeter Brune 
53724989b8cSPeter Brune    Input Arguments:
53824989b8cSPeter Brune +  dm - DM to be used with TS
53924989b8cSPeter Brune .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
54024989b8cSPeter Brune -  ctx - context for residual evaluation
54124989b8cSPeter Brune 
54224989b8cSPeter Brune    Level: advanced
54324989b8cSPeter Brune 
54424989b8cSPeter Brune    Note:
545ef20d060SBarry Smith    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
54624989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
54724989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
54824989b8cSPeter Brune 
54924989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
55024989b8cSPeter Brune @*/
55124989b8cSPeter Brune PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
55224989b8cSPeter Brune {
55324989b8cSPeter Brune   PetscErrorCode ierr;
554942e3340SBarry Smith   DMTS           tsdm;
55524989b8cSPeter Brune 
55624989b8cSPeter Brune   PetscFunctionBegin;
55724989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
558942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
559d74926cbSBarry Smith   if (func) tsdm->ops->rhsfunction = func;
56024989b8cSPeter Brune   if (ctx)  tsdm->rhsfunctionctx = ctx;
56124989b8cSPeter Brune   PetscFunctionReturn(0);
56224989b8cSPeter Brune }
56324989b8cSPeter Brune 
56424989b8cSPeter Brune #undef __FUNCT__
565ef20d060SBarry Smith #define __FUNCT__ "DMTSGetSolutionFunction"
566ef20d060SBarry Smith /*@C
567ef20d060SBarry Smith    DMTSGetSolutionFunction - gets the TS solution evaluation function
568ef20d060SBarry Smith 
569ef20d060SBarry Smith    Not Collective
570ef20d060SBarry Smith 
571ef20d060SBarry Smith    Input Arguments:
572ef20d060SBarry Smith .  dm - DM to be used with TS
573ef20d060SBarry Smith 
574ef20d060SBarry Smith    Output Parameters:
575ef20d060SBarry Smith +  func - solution function evaluation function, see TSSetSolution() for calling sequence
576ef20d060SBarry Smith -  ctx - context for solution evaluation
577ef20d060SBarry Smith 
578ef20d060SBarry Smith    Level: advanced
579ef20d060SBarry Smith 
580ef20d060SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
581ef20d060SBarry Smith @*/
582ef20d060SBarry Smith PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
583ef20d060SBarry Smith {
584ef20d060SBarry Smith   PetscErrorCode ierr;
585942e3340SBarry Smith   DMTS           tsdm;
586ef20d060SBarry Smith 
587ef20d060SBarry Smith   PetscFunctionBegin;
588ef20d060SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
589942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
590d74926cbSBarry Smith   if (func) *func = tsdm->ops->solution;
591ef20d060SBarry Smith   if (ctx)  *ctx  = tsdm->solutionctx;
592ef20d060SBarry Smith   PetscFunctionReturn(0);
593ef20d060SBarry Smith }
594ef20d060SBarry Smith 
595ef20d060SBarry Smith #undef __FUNCT__
596ef20d060SBarry Smith #define __FUNCT__ "DMTSSetSolutionFunction"
597ef20d060SBarry Smith /*@C
598ef20d060SBarry Smith    DMTSSetSolutionFunction - set TS solution evaluation function
599ef20d060SBarry Smith 
600ef20d060SBarry Smith    Not Collective
601ef20d060SBarry Smith 
602ef20d060SBarry Smith    Input Arguments:
603ef20d060SBarry Smith +  dm - DM to be used with TS
604ef20d060SBarry Smith .  func - solution function evaluation function, see TSSetSolution() for calling sequence
605ef20d060SBarry Smith -  ctx - context for solution evaluation
606ef20d060SBarry Smith 
607ef20d060SBarry Smith    Level: advanced
608ef20d060SBarry Smith 
609ef20d060SBarry Smith    Note:
610ef20d060SBarry Smith    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
611ef20d060SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
612ef20d060SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
613ef20d060SBarry Smith 
614ef20d060SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
615ef20d060SBarry Smith @*/
616ef20d060SBarry Smith PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
617ef20d060SBarry Smith {
618ef20d060SBarry Smith   PetscErrorCode ierr;
619942e3340SBarry Smith   DMTS           tsdm;
620ef20d060SBarry Smith 
621ef20d060SBarry Smith   PetscFunctionBegin;
622ef20d060SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
623942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
624d74926cbSBarry Smith   if (func) tsdm->ops->solution = func;
625ef20d060SBarry Smith   if (ctx)  tsdm->solutionctx   = ctx;
626ef20d060SBarry Smith   PetscFunctionReturn(0);
627ef20d060SBarry Smith }
628ef20d060SBarry Smith 
629ef20d060SBarry Smith #undef __FUNCT__
6309b7cd975SBarry Smith #define __FUNCT__ "DMTSSetForcingFunction"
6319b7cd975SBarry Smith /*@C
6329b7cd975SBarry Smith    DMTSSetForcingFunction - set TS forcing function evaluation function
6339b7cd975SBarry Smith 
6349b7cd975SBarry Smith    Not Collective
6359b7cd975SBarry Smith 
6369b7cd975SBarry Smith    Input Arguments:
6379b7cd975SBarry Smith +  dm - DM to be used with TS
638f8b49ee9SBarry Smith .  f - forcing function evaluation function; see TSForcingFunction
6399b7cd975SBarry Smith -  ctx - context for solution evaluation
6409b7cd975SBarry Smith 
6419b7cd975SBarry Smith    Level: advanced
6429b7cd975SBarry Smith 
6439b7cd975SBarry Smith    Note:
6449b7cd975SBarry Smith    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
6459b7cd975SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
6469b7cd975SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
6479b7cd975SBarry Smith 
6489b7cd975SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
6499b7cd975SBarry Smith @*/
650d56366bfSLisandro Dalcin PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx)
6519b7cd975SBarry Smith {
6529b7cd975SBarry Smith   PetscErrorCode ierr;
6539b7cd975SBarry Smith   DMTS           tsdm;
6549b7cd975SBarry Smith 
6559b7cd975SBarry Smith   PetscFunctionBegin;
6569b7cd975SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6579b7cd975SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
658f8b49ee9SBarry Smith   if (f)    tsdm->ops->forcing = f;
6599b7cd975SBarry Smith   if (ctx)  tsdm->forcingctx   = ctx;
6609b7cd975SBarry Smith   PetscFunctionReturn(0);
6619b7cd975SBarry Smith }
6629b7cd975SBarry Smith 
6639b7cd975SBarry Smith 
6649b7cd975SBarry Smith #undef __FUNCT__
6659b7cd975SBarry Smith #define __FUNCT__ "DMTSGetForcingFunction"
6669b7cd975SBarry Smith /*@C
6679b7cd975SBarry Smith    DMTSGetForcingFunction - get TS forcing function evaluation function
6689b7cd975SBarry Smith 
6699b7cd975SBarry Smith    Not Collective
6709b7cd975SBarry Smith 
6719b7cd975SBarry Smith    Input Argument:
6729b7cd975SBarry Smith .   dm - DM to be used with TS
6739b7cd975SBarry Smith 
6749b7cd975SBarry Smith    Output Arguments:
675f8b49ee9SBarry Smith +  f - forcing function evaluation function; see TSForcingFunction for details
6769b7cd975SBarry Smith -  ctx - context for solution evaluation
6779b7cd975SBarry Smith 
6789b7cd975SBarry Smith    Level: advanced
6799b7cd975SBarry Smith 
6809b7cd975SBarry Smith    Note:
6819b7cd975SBarry Smith    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
6829b7cd975SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
6839b7cd975SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
6849b7cd975SBarry Smith 
6859b7cd975SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
6869b7cd975SBarry Smith @*/
687d56366bfSLisandro Dalcin PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
6889b7cd975SBarry Smith {
6899b7cd975SBarry Smith   PetscErrorCode ierr;
6909b7cd975SBarry Smith   DMTS           tsdm;
6919b7cd975SBarry Smith 
6929b7cd975SBarry Smith   PetscFunctionBegin;
6939b7cd975SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6949b7cd975SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
695f8b49ee9SBarry Smith   if (f)   *f   = tsdm->ops->forcing;
6969b7cd975SBarry Smith   if (ctx) *ctx = tsdm->forcingctx;
6979b7cd975SBarry Smith   PetscFunctionReturn(0);
6989b7cd975SBarry Smith }
6999b7cd975SBarry Smith 
7009b7cd975SBarry Smith #undef __FUNCT__
70124989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSFunction"
70224989b8cSPeter Brune /*@C
70324989b8cSPeter Brune    DMTSGetRHSFunction - get TS explicit residual evaluation function
70424989b8cSPeter Brune 
70524989b8cSPeter Brune    Not Collective
70624989b8cSPeter Brune 
70724989b8cSPeter Brune    Input Argument:
70824989b8cSPeter Brune .  dm - DM to be used with TS
70924989b8cSPeter Brune 
71024989b8cSPeter Brune    Output Arguments:
71124989b8cSPeter Brune +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
71224989b8cSPeter Brune -  ctx - context for residual evaluation
71324989b8cSPeter Brune 
71424989b8cSPeter Brune    Level: advanced
71524989b8cSPeter Brune 
71624989b8cSPeter Brune    Note:
71724989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
71824989b8cSPeter Brune    associated with the DM.
71924989b8cSPeter Brune 
72024989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
72124989b8cSPeter Brune @*/
72224989b8cSPeter Brune PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
72324989b8cSPeter Brune {
72424989b8cSPeter Brune   PetscErrorCode ierr;
725942e3340SBarry Smith   DMTS           tsdm;
72624989b8cSPeter Brune 
72724989b8cSPeter Brune   PetscFunctionBegin;
72824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
729942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
730d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsfunction;
73124989b8cSPeter Brune   if (ctx)  *ctx = tsdm->rhsfunctionctx;
73224989b8cSPeter Brune   PetscFunctionReturn(0);
73324989b8cSPeter Brune }
73424989b8cSPeter Brune 
73524989b8cSPeter Brune #undef __FUNCT__
73624989b8cSPeter Brune #define __FUNCT__ "DMTSSetIJacobian"
73724989b8cSPeter Brune /*@C
73824989b8cSPeter Brune    DMTSSetIJacobian - set TS Jacobian evaluation function
73924989b8cSPeter Brune 
74024989b8cSPeter Brune    Not Collective
74124989b8cSPeter Brune 
74224989b8cSPeter Brune    Input Argument:
74324989b8cSPeter Brune +  dm - DM to be used with TS
74424989b8cSPeter Brune .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
74524989b8cSPeter Brune -  ctx - context for residual evaluation
74624989b8cSPeter Brune 
74724989b8cSPeter Brune    Level: advanced
74824989b8cSPeter Brune 
74924989b8cSPeter Brune    Note:
75024989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
75124989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
75224989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
75324989b8cSPeter Brune 
75424989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
75524989b8cSPeter Brune @*/
75624989b8cSPeter Brune PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
75724989b8cSPeter Brune {
75824989b8cSPeter Brune   PetscErrorCode ierr;
759942e3340SBarry Smith   DMTS           sdm;
76024989b8cSPeter Brune 
76124989b8cSPeter Brune   PetscFunctionBegin;
76224989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
763942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
764d74926cbSBarry Smith   if (func) sdm->ops->ijacobian = func;
76524989b8cSPeter Brune   if (ctx)  sdm->ijacobianctx   = ctx;
76624989b8cSPeter Brune   PetscFunctionReturn(0);
76724989b8cSPeter Brune }
76824989b8cSPeter Brune 
76924989b8cSPeter Brune #undef __FUNCT__
77024989b8cSPeter Brune #define __FUNCT__ "DMTSGetIJacobian"
77124989b8cSPeter Brune /*@C
77224989b8cSPeter Brune    DMTSGetIJacobian - get TS Jacobian evaluation function
77324989b8cSPeter Brune 
77424989b8cSPeter Brune    Not Collective
77524989b8cSPeter Brune 
77624989b8cSPeter Brune    Input Argument:
77724989b8cSPeter Brune .  dm - DM to be used with TS
77824989b8cSPeter Brune 
77924989b8cSPeter Brune    Output Arguments:
78024989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
78124989b8cSPeter Brune -  ctx - context for residual evaluation
78224989b8cSPeter Brune 
78324989b8cSPeter Brune    Level: advanced
78424989b8cSPeter Brune 
78524989b8cSPeter Brune    Note:
78624989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
78724989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
78824989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
78924989b8cSPeter Brune 
79024989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
79124989b8cSPeter Brune @*/
79224989b8cSPeter Brune PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
79324989b8cSPeter Brune {
79424989b8cSPeter Brune   PetscErrorCode ierr;
795942e3340SBarry Smith   DMTS           tsdm;
79624989b8cSPeter Brune 
79724989b8cSPeter Brune   PetscFunctionBegin;
79824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
799942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
800d74926cbSBarry Smith   if (func) *func = tsdm->ops->ijacobian;
80124989b8cSPeter Brune   if (ctx)  *ctx = tsdm->ijacobianctx;
80224989b8cSPeter Brune   PetscFunctionReturn(0);
80324989b8cSPeter Brune }
80424989b8cSPeter Brune 
80524989b8cSPeter Brune 
80624989b8cSPeter Brune #undef __FUNCT__
80724989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSJacobian"
80824989b8cSPeter Brune /*@C
80924989b8cSPeter Brune    DMTSSetRHSJacobian - set TS Jacobian evaluation function
81024989b8cSPeter Brune 
81124989b8cSPeter Brune    Not Collective
81224989b8cSPeter Brune 
81324989b8cSPeter Brune    Input Argument:
81424989b8cSPeter Brune +  dm - DM to be used with TS
81524989b8cSPeter Brune .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
81624989b8cSPeter Brune -  ctx - context for residual evaluation
81724989b8cSPeter Brune 
81824989b8cSPeter Brune    Level: advanced
81924989b8cSPeter Brune 
82024989b8cSPeter Brune    Note:
82124989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
82224989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
82324989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
82424989b8cSPeter Brune 
82524989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
82624989b8cSPeter Brune @*/
82724989b8cSPeter Brune PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
82824989b8cSPeter Brune {
82924989b8cSPeter Brune   PetscErrorCode ierr;
830942e3340SBarry Smith   DMTS           tsdm;
83124989b8cSPeter Brune 
83224989b8cSPeter Brune   PetscFunctionBegin;
83324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
834942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
835d74926cbSBarry Smith   if (func) tsdm->ops->rhsjacobian = func;
83624989b8cSPeter Brune   if (ctx)  tsdm->rhsjacobianctx = ctx;
83724989b8cSPeter Brune   PetscFunctionReturn(0);
83824989b8cSPeter Brune }
83924989b8cSPeter Brune 
84024989b8cSPeter Brune #undef __FUNCT__
84124989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSJacobian"
84224989b8cSPeter Brune /*@C
84324989b8cSPeter Brune    DMTSGetRHSJacobian - get TS Jacobian evaluation function
84424989b8cSPeter Brune 
84524989b8cSPeter Brune    Not Collective
84624989b8cSPeter Brune 
84724989b8cSPeter Brune    Input Argument:
84824989b8cSPeter Brune .  dm - DM to be used with TS
84924989b8cSPeter Brune 
85024989b8cSPeter Brune    Output Arguments:
85124989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
85224989b8cSPeter Brune -  ctx - context for residual evaluation
85324989b8cSPeter Brune 
85424989b8cSPeter Brune    Level: advanced
85524989b8cSPeter Brune 
85624989b8cSPeter Brune    Note:
85724989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
85824989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
85924989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
86024989b8cSPeter Brune 
86124989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
86224989b8cSPeter Brune @*/
86324989b8cSPeter Brune PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
86424989b8cSPeter Brune {
86524989b8cSPeter Brune   PetscErrorCode ierr;
866942e3340SBarry Smith   DMTS           tsdm;
86724989b8cSPeter Brune 
86824989b8cSPeter Brune   PetscFunctionBegin;
86924989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
870942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
871d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsjacobian;
87224989b8cSPeter Brune   if (ctx)  *ctx = tsdm->rhsjacobianctx;
87324989b8cSPeter Brune   PetscFunctionReturn(0);
87424989b8cSPeter Brune }
875ad6bc421SBarry Smith 
876ad6bc421SBarry Smith #undef __FUNCT__
877ad6bc421SBarry Smith #define __FUNCT__ "DMTSSetIFunctionSerialize"
878ad6bc421SBarry Smith /*@C
879ad6bc421SBarry Smith    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
880ad6bc421SBarry Smith 
881ad6bc421SBarry Smith    Not Collective
882ad6bc421SBarry Smith 
883ad6bc421SBarry Smith    Input Arguments:
884ad6bc421SBarry Smith +  dm - DM to be used with TS
885ad6bc421SBarry Smith .  view - viewer function
886ad6bc421SBarry Smith -  load - loading function
887ad6bc421SBarry Smith 
888ad6bc421SBarry Smith    Level: advanced
889ad6bc421SBarry Smith 
890ad6bc421SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
891ad6bc421SBarry Smith @*/
892ad6bc421SBarry Smith PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
893ad6bc421SBarry Smith {
894ad6bc421SBarry Smith   PetscErrorCode ierr;
895ad6bc421SBarry Smith   DMTS           tsdm;
896ad6bc421SBarry Smith 
897ad6bc421SBarry Smith   PetscFunctionBegin;
898ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
899ad6bc421SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
900ad6bc421SBarry Smith   tsdm->ops->ifunctionview = view;
901ad6bc421SBarry Smith   tsdm->ops->ifunctionload = load;
902ad6bc421SBarry Smith   PetscFunctionReturn(0);
903ad6bc421SBarry Smith }
904ad6bc421SBarry Smith 
905ad6bc421SBarry Smith #undef __FUNCT__
906ad6bc421SBarry Smith #define __FUNCT__ "DMTSSetIJacobianSerialize"
907ad6bc421SBarry Smith /*@C
908ad6bc421SBarry Smith    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
909ad6bc421SBarry Smith 
910ad6bc421SBarry Smith    Not Collective
911ad6bc421SBarry Smith 
912ad6bc421SBarry Smith    Input Arguments:
913ad6bc421SBarry Smith +  dm - DM to be used with TS
914ad6bc421SBarry Smith .  view - viewer function
915ad6bc421SBarry Smith -  load - loading function
916ad6bc421SBarry Smith 
917ad6bc421SBarry Smith    Level: advanced
918ad6bc421SBarry Smith 
919ad6bc421SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
920ad6bc421SBarry Smith @*/
921ad6bc421SBarry Smith PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
922ad6bc421SBarry Smith {
923ad6bc421SBarry Smith   PetscErrorCode ierr;
924ad6bc421SBarry Smith   DMTS           tsdm;
925ad6bc421SBarry Smith 
926ad6bc421SBarry Smith   PetscFunctionBegin;
927ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
928ad6bc421SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
929ad6bc421SBarry Smith   tsdm->ops->ijacobianview = view;
930ad6bc421SBarry Smith   tsdm->ops->ijacobianload = load;
931ad6bc421SBarry Smith   PetscFunctionReturn(0);
932ad6bc421SBarry Smith }
933