xref: /petsc/src/ts/utils/dmts.c (revision 3964eb8815a373d98cd54ebd45324cfa27b228df)
124989b8cSPeter Brune #include <petsc-private/tsimpl.h>     /*I "petscts.h" I*/
207475bc1SBarry 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;
262d53ad75SBarry Smith   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,PETSC_FUNCTION);CHKERRQ(ierr);
27ad6bc421SBarry Smith   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,PETSC_FUNCTION);CHKERRQ(ierr);
28ad6bc421SBarry Smith   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,PETSC_FUNCTION);CHKERRQ(ierr);
29ad6bc421SBarry Smith   if (kdm->ops->ifunctionload) {
30ad6bc421SBarry Smith     ierr = (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);CHKERRQ(ierr);
31ad6bc421SBarry Smith   }
322d53ad75SBarry Smith   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,PETSC_FUNCTION);CHKERRQ(ierr);
33ad6bc421SBarry Smith   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,PETSC_FUNCTION);CHKERRQ(ierr);
34ad6bc421SBarry Smith   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,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) {
65*3964eb88SJed Brown     struct {
66*3964eb88SJed Brown       TSIFunction ifunction;
67*3964eb88SJed Brown       PetscErrorCode (*ifunctionview)(void*,PetscViewer);
68*3964eb88SJed Brown       PetscErrorCode (*ifunctionload)(void**,PetscViewer);
69*3964eb88SJed Brown     } funcstruct = {kdm->ops->ifunction,
70*3964eb88SJed Brown                     kdm->ops->ifunctionview,
71*3964eb88SJed Brown                     kdm->ops->ifunctionload};
72*3964eb88SJed Brown     struct {
73*3964eb88SJed Brown       TSIJacobian ijacobian;
74*3964eb88SJed Brown       PetscErrorCode (*ijacobianview)(void*,PetscViewer);
75*3964eb88SJed Brown       PetscErrorCode (*ijacobianload)(void**,PetscViewer);
76*3964eb88SJed Brown     } jacstruct = {kdm->ops->ijacobian,
77*3964eb88SJed Brown                    kdm->ops->ijacobianview,
78*3964eb88SJed Brown                    kdm->ops->ijacobianload};
79*3964eb88SJed Brown 
80*3964eb88SJed Brown     ierr = PetscViewerBinaryWrite(viewer,&funcstruct,3,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
81ad6bc421SBarry Smith     if (kdm->ops->ifunctionview) {
82ad6bc421SBarry Smith       ierr = (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);CHKERRQ(ierr);
83ad6bc421SBarry Smith     }
84*3964eb88SJed Brown     ierr = PetscViewerBinaryWrite(viewer,&jacstruct,3,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
85ad6bc421SBarry Smith     if (kdm->ops->ijacobianview) {
86ad6bc421SBarry Smith       ierr = (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);CHKERRQ(ierr);
87ad6bc421SBarry Smith     }
882d53ad75SBarry Smith   }
892d53ad75SBarry Smith   PetscFunctionReturn(0);
902d53ad75SBarry Smith }
912d53ad75SBarry Smith 
922d53ad75SBarry Smith #undef __FUNCT__
93d74926cbSBarry Smith #define __FUNCT__ "DMTSCreate"
94d74926cbSBarry Smith static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
95d74926cbSBarry Smith {
96d74926cbSBarry Smith   PetscErrorCode ierr;
97d74926cbSBarry Smith 
98d74926cbSBarry Smith   PetscFunctionBegin;
99519f805aSKarl Rupp #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1000298fd71SBarry Smith   ierr = TSInitializePackage(NULL);CHKERRQ(ierr);
101d74926cbSBarry Smith #endif
10267c2884eSBarry Smith   ierr = PetscHeaderCreate(*kdm, _p_DMTS, struct _DMTSOps, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);CHKERRQ(ierr);
103d74926cbSBarry Smith   ierr = PetscMemzero((*kdm)->ops, sizeof(struct _DMTSOps));CHKERRQ(ierr);
104d74926cbSBarry Smith   PetscFunctionReturn(0);
105d74926cbSBarry Smith }
10624989b8cSPeter Brune 
10724989b8cSPeter Brune #undef __FUNCT__
108942e3340SBarry Smith #define __FUNCT__ "DMCoarsenHook_DMTS"
1092a34c10cSBarry Smith /* Attaches the DMTS to the coarse level.
11024989b8cSPeter Brune  * Under what conditions should we copy versus duplicate?
11124989b8cSPeter Brune  */
112942e3340SBarry Smith static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
11324989b8cSPeter Brune {
11424989b8cSPeter Brune   PetscErrorCode ierr;
11524989b8cSPeter Brune 
11624989b8cSPeter Brune   PetscFunctionBegin;
117942e3340SBarry Smith   ierr = DMCopyDMTS(dm,dmc);CHKERRQ(ierr);
11824989b8cSPeter Brune   PetscFunctionReturn(0);
11924989b8cSPeter Brune }
12024989b8cSPeter Brune 
12124989b8cSPeter Brune #undef __FUNCT__
122942e3340SBarry Smith #define __FUNCT__ "DMRestrictHook_DMTS"
12324989b8cSPeter Brune /* This could restrict auxiliary information to the coarse level.
12424989b8cSPeter Brune  */
125942e3340SBarry Smith static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
12624989b8cSPeter Brune {
12724989b8cSPeter Brune 
12824989b8cSPeter Brune   PetscFunctionBegin;
12924989b8cSPeter Brune   PetscFunctionReturn(0);
13024989b8cSPeter Brune }
13124989b8cSPeter Brune 
13224989b8cSPeter Brune #undef __FUNCT__
133258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_DMTS"
134258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
135258e1594SPeter Brune {
136258e1594SPeter Brune   PetscErrorCode ierr;
137258e1594SPeter Brune 
138258e1594SPeter Brune   PetscFunctionBegin;
139258e1594SPeter Brune   ierr = DMCopyDMTS(dm,subdm);CHKERRQ(ierr);
140258e1594SPeter Brune   PetscFunctionReturn(0);
141258e1594SPeter Brune }
142258e1594SPeter Brune 
143258e1594SPeter Brune #undef __FUNCT__
144258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_DMTS"
145258e1594SPeter Brune /* This could restrict auxiliary information to the coarse level.
146258e1594SPeter Brune  */
147258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
148258e1594SPeter Brune {
149258e1594SPeter Brune   PetscFunctionBegin;
150258e1594SPeter Brune   PetscFunctionReturn(0);
151258e1594SPeter Brune }
152258e1594SPeter Brune 
153258e1594SPeter Brune #undef __FUNCT__
154d74926cbSBarry Smith #define __FUNCT__ "DMTSCopy"
155d74926cbSBarry Smith /*@C
156d74926cbSBarry Smith    DMTSCopy - copies the information in a DMTS to another DMTS
157d74926cbSBarry Smith 
158d74926cbSBarry Smith    Not Collective
159d74926cbSBarry Smith 
160d74926cbSBarry Smith    Input Argument:
161d74926cbSBarry Smith +  kdm - Original DMTS
162d74926cbSBarry Smith -  nkdm - DMTS to receive the data, should have been created with DMTSCreate()
163d74926cbSBarry Smith 
164d74926cbSBarry Smith    Level: developer
165d74926cbSBarry Smith 
166d74926cbSBarry Smith .seealso: DMTSCreate(), DMTSDestroy()
167d74926cbSBarry Smith @*/
168d74926cbSBarry Smith PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
16924989b8cSPeter Brune {
17024989b8cSPeter Brune   PetscErrorCode ierr;
17124989b8cSPeter Brune 
17224989b8cSPeter Brune   PetscFunctionBegin;
173d74926cbSBarry Smith   PetscValidHeaderSpecific(kdm,DMTS_CLASSID,1);
174d74926cbSBarry Smith   PetscValidHeaderSpecific(nkdm,DMTS_CLASSID,2);
175d74926cbSBarry Smith   nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
176d74926cbSBarry Smith   nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
177d74926cbSBarry Smith   nkdm->ops->ifunction   = kdm->ops->ifunction;
178d74926cbSBarry Smith   nkdm->ops->ijacobian   = kdm->ops->ijacobian;
179d74926cbSBarry Smith   nkdm->ops->solution    = kdm->ops->solution;
180d74926cbSBarry Smith   nkdm->ops->destroy     = kdm->ops->destroy;
181d74926cbSBarry Smith   nkdm->ops->duplicate   = kdm->ops->duplicate;
182d74926cbSBarry Smith 
183d74926cbSBarry Smith   nkdm->rhsfunctionctx = kdm->rhsfunctionctx;
184d74926cbSBarry Smith   nkdm->rhsjacobianctx = kdm->rhsjacobianctx;
185d74926cbSBarry Smith   nkdm->ifunctionctx   = kdm->ifunctionctx;
186d74926cbSBarry Smith   nkdm->ijacobianctx   = kdm->ijacobianctx;
187d74926cbSBarry Smith   nkdm->solutionctx    = kdm->solutionctx;
188d74926cbSBarry Smith 
189d74926cbSBarry Smith   nkdm->data = kdm->data;
190d74926cbSBarry Smith 
191d74926cbSBarry Smith   /*
192d74926cbSBarry Smith   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
193d74926cbSBarry Smith   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
194d74926cbSBarry Smith   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
195d74926cbSBarry Smith   */
196d74926cbSBarry Smith 
197d74926cbSBarry Smith   /* implementation specific copy hooks */
198d74926cbSBarry Smith   if (kdm->ops->duplicate) {ierr = (*kdm->ops->duplicate)(kdm,nkdm);CHKERRQ(ierr);}
19924989b8cSPeter Brune   PetscFunctionReturn(0);
20024989b8cSPeter Brune }
20124989b8cSPeter Brune 
20224989b8cSPeter Brune #undef __FUNCT__
203942e3340SBarry Smith #define __FUNCT__ "DMGetDMTS"
20424989b8cSPeter Brune /*@C
205942e3340SBarry Smith    DMGetDMTS - get read-only private DMTS context from a DM
20624989b8cSPeter Brune 
20724989b8cSPeter Brune    Not Collective
20824989b8cSPeter Brune 
20924989b8cSPeter Brune    Input Argument:
21024989b8cSPeter Brune .  dm - DM to be used with TS
21124989b8cSPeter Brune 
21224989b8cSPeter Brune    Output Argument:
213942e3340SBarry Smith .  tsdm - private DMTS context
21424989b8cSPeter Brune 
21524989b8cSPeter Brune    Level: developer
21624989b8cSPeter Brune 
21724989b8cSPeter Brune    Notes:
218942e3340SBarry Smith    Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.
21924989b8cSPeter Brune 
220942e3340SBarry Smith .seealso: DMGetDMTSWrite()
22124989b8cSPeter Brune @*/
222942e3340SBarry Smith PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
22324989b8cSPeter Brune {
22424989b8cSPeter Brune   PetscErrorCode ierr;
22524989b8cSPeter Brune 
22624989b8cSPeter Brune   PetscFunctionBegin;
22724989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2282a34c10cSBarry Smith   *tsdm = (DMTS) dm->dmts;
229d74926cbSBarry Smith   if (!*tsdm) {
230942e3340SBarry Smith     ierr = PetscInfo(dm,"Creating new DMTS\n");CHKERRQ(ierr);
231ce94432eSBarry Smith     ierr = DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm);CHKERRQ(ierr);
2322a34c10cSBarry Smith     dm->dmts = (PetscObject) *tsdm;
2330298fd71SBarry Smith     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr);
2340298fd71SBarry Smith     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr);
23524989b8cSPeter Brune   }
23624989b8cSPeter Brune   PetscFunctionReturn(0);
23724989b8cSPeter Brune }
23824989b8cSPeter Brune 
23924989b8cSPeter Brune #undef __FUNCT__
240942e3340SBarry Smith #define __FUNCT__ "DMGetDMTSWrite"
24124989b8cSPeter Brune /*@C
242942e3340SBarry Smith    DMGetDMTSWrite - get write access to private DMTS context from a DM
24324989b8cSPeter Brune 
24424989b8cSPeter Brune    Not Collective
24524989b8cSPeter Brune 
24624989b8cSPeter Brune    Input Argument:
24724989b8cSPeter Brune .  dm - DM to be used with TS
24824989b8cSPeter Brune 
24924989b8cSPeter Brune    Output Argument:
250942e3340SBarry Smith .  tsdm - private DMTS context
25124989b8cSPeter Brune 
25224989b8cSPeter Brune    Level: developer
25324989b8cSPeter Brune 
254942e3340SBarry Smith .seealso: DMGetDMTS()
25524989b8cSPeter Brune @*/
256942e3340SBarry Smith PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
25724989b8cSPeter Brune {
25824989b8cSPeter Brune   PetscErrorCode ierr;
259942e3340SBarry Smith   DMTS           sdm;
26024989b8cSPeter Brune 
26124989b8cSPeter Brune   PetscFunctionBegin;
26224989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
263942e3340SBarry Smith   ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr);
26424989b8cSPeter Brune   if (!sdm->originaldm) sdm->originaldm = dm;
26524989b8cSPeter Brune   if (sdm->originaldm != dm) {  /* Copy on write */
2662a34c10cSBarry Smith     DMTS oldsdm = sdm;
267942e3340SBarry Smith     ierr     = PetscInfo(dm,"Copying DMTS due to write\n");CHKERRQ(ierr);
268ce94432eSBarry Smith     ierr     = DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);CHKERRQ(ierr);
269d74926cbSBarry Smith     ierr     = DMTSCopy(oldsdm,sdm);CHKERRQ(ierr);
2702a34c10cSBarry Smith     ierr     = DMTSDestroy((DMTS*)&dm->dmts);CHKERRQ(ierr);
2712a34c10cSBarry Smith     dm->dmts = (PetscObject) sdm;
27224989b8cSPeter Brune   }
27324989b8cSPeter Brune   *tsdm = sdm;
27424989b8cSPeter Brune   PetscFunctionReturn(0);
27524989b8cSPeter Brune }
27624989b8cSPeter Brune 
27724989b8cSPeter Brune #undef __FUNCT__
278942e3340SBarry Smith #define __FUNCT__ "DMCopyDMTS"
27924989b8cSPeter Brune /*@C
280942e3340SBarry Smith    DMCopyDMTS - copies a DM context to a new DM
28124989b8cSPeter Brune 
28224989b8cSPeter Brune    Logically Collective
28324989b8cSPeter Brune 
28424989b8cSPeter Brune    Input Arguments:
28524989b8cSPeter Brune +  dmsrc - DM to obtain context from
28624989b8cSPeter Brune -  dmdest - DM to add context to
28724989b8cSPeter Brune 
28824989b8cSPeter Brune    Level: developer
28924989b8cSPeter Brune 
29024989b8cSPeter Brune    Note:
29124989b8cSPeter Brune    The context is copied by reference. This function does not ensure that a context exists.
29224989b8cSPeter Brune 
293942e3340SBarry Smith .seealso: DMGetDMTS(), TSSetDM()
29424989b8cSPeter Brune @*/
295942e3340SBarry Smith PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
29624989b8cSPeter Brune {
29724989b8cSPeter Brune   PetscErrorCode ierr;
29824989b8cSPeter Brune 
29924989b8cSPeter Brune   PetscFunctionBegin;
30024989b8cSPeter Brune   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
30124989b8cSPeter Brune   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
3022a34c10cSBarry Smith   ierr         = DMTSDestroy((DMTS*)&dmdest->dmts);CHKERRQ(ierr);
3032a34c10cSBarry Smith   dmdest->dmts = dmsrc->dmts;
3042a34c10cSBarry Smith   ierr         = PetscObjectReference(dmdest->dmts);CHKERRQ(ierr);
3050298fd71SBarry Smith   ierr         = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr);
3060298fd71SBarry Smith   ierr         = DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr);
30724989b8cSPeter Brune   PetscFunctionReturn(0);
30824989b8cSPeter Brune }
30924989b8cSPeter Brune 
31024989b8cSPeter Brune #undef __FUNCT__
31124989b8cSPeter Brune #define __FUNCT__ "DMTSSetIFunction"
31224989b8cSPeter Brune /*@C
31324989b8cSPeter Brune    DMTSSetIFunction - set TS implicit function evaluation function
31424989b8cSPeter Brune 
31524989b8cSPeter Brune    Not Collective
31624989b8cSPeter Brune 
31724989b8cSPeter Brune    Input Arguments:
31824989b8cSPeter Brune +  dm - DM to be used with TS
31924989b8cSPeter Brune .  func - function evaluation function, see TSSetIFunction() for calling sequence
32024989b8cSPeter Brune -  ctx - context for residual evaluation
32124989b8cSPeter Brune 
32224989b8cSPeter Brune    Level: advanced
32324989b8cSPeter Brune 
32424989b8cSPeter Brune    Note:
32524989b8cSPeter Brune    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
32624989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
32724989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
32824989b8cSPeter Brune 
32924989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
33024989b8cSPeter Brune @*/
33124989b8cSPeter Brune PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
33224989b8cSPeter Brune {
33324989b8cSPeter Brune   PetscErrorCode ierr;
334942e3340SBarry Smith   DMTS           tsdm;
33524989b8cSPeter Brune 
33624989b8cSPeter Brune   PetscFunctionBegin;
33724989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
338942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
339d74926cbSBarry Smith   if (func) tsdm->ops->ifunction = func;
34024989b8cSPeter Brune   if (ctx)  tsdm->ifunctionctx = ctx;
34124989b8cSPeter Brune   PetscFunctionReturn(0);
34224989b8cSPeter Brune }
34324989b8cSPeter Brune 
34424989b8cSPeter Brune #undef __FUNCT__
34524989b8cSPeter Brune #define __FUNCT__ "DMTSGetIFunction"
34624989b8cSPeter Brune /*@C
34724989b8cSPeter Brune    DMTSGetIFunction - get TS implicit residual evaluation function
34824989b8cSPeter Brune 
34924989b8cSPeter Brune    Not Collective
35024989b8cSPeter Brune 
35124989b8cSPeter Brune    Input Argument:
35224989b8cSPeter Brune .  dm - DM to be used with TS
35324989b8cSPeter Brune 
35424989b8cSPeter Brune    Output Arguments:
35524989b8cSPeter Brune +  func - function evaluation function, see TSSetIFunction() for calling sequence
35624989b8cSPeter Brune -  ctx - context for residual evaluation
35724989b8cSPeter Brune 
35824989b8cSPeter Brune    Level: advanced
35924989b8cSPeter Brune 
36024989b8cSPeter Brune    Note:
36124989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
36224989b8cSPeter Brune    associated with the DM.
36324989b8cSPeter Brune 
36424989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
36524989b8cSPeter Brune @*/
36624989b8cSPeter Brune PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
36724989b8cSPeter Brune {
36824989b8cSPeter Brune   PetscErrorCode ierr;
369942e3340SBarry Smith   DMTS           tsdm;
37024989b8cSPeter Brune 
37124989b8cSPeter Brune   PetscFunctionBegin;
37224989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
373942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
374d74926cbSBarry Smith   if (func) *func = tsdm->ops->ifunction;
37524989b8cSPeter Brune   if (ctx)  *ctx = tsdm->ifunctionctx;
37624989b8cSPeter Brune   PetscFunctionReturn(0);
37724989b8cSPeter Brune }
37824989b8cSPeter Brune 
37924989b8cSPeter Brune 
38024989b8cSPeter Brune #undef __FUNCT__
38124989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSFunction"
38224989b8cSPeter Brune /*@C
38324989b8cSPeter Brune    DMTSSetRHSFunction - set TS explicit residual evaluation function
38424989b8cSPeter Brune 
38524989b8cSPeter Brune    Not Collective
38624989b8cSPeter Brune 
38724989b8cSPeter Brune    Input Arguments:
38824989b8cSPeter Brune +  dm - DM to be used with TS
38924989b8cSPeter Brune .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
39024989b8cSPeter Brune -  ctx - context for residual evaluation
39124989b8cSPeter Brune 
39224989b8cSPeter Brune    Level: advanced
39324989b8cSPeter Brune 
39424989b8cSPeter Brune    Note:
395ef20d060SBarry Smith    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
39624989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
39724989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
39824989b8cSPeter Brune 
39924989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
40024989b8cSPeter Brune @*/
40124989b8cSPeter Brune PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
40224989b8cSPeter Brune {
40324989b8cSPeter Brune   PetscErrorCode ierr;
404942e3340SBarry Smith   DMTS           tsdm;
40524989b8cSPeter Brune 
40624989b8cSPeter Brune   PetscFunctionBegin;
40724989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
408942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
409d74926cbSBarry Smith   if (func) tsdm->ops->rhsfunction = func;
41024989b8cSPeter Brune   if (ctx)  tsdm->rhsfunctionctx = ctx;
41124989b8cSPeter Brune   PetscFunctionReturn(0);
41224989b8cSPeter Brune }
41324989b8cSPeter Brune 
41424989b8cSPeter Brune #undef __FUNCT__
415ef20d060SBarry Smith #define __FUNCT__ "DMTSGetSolutionFunction"
416ef20d060SBarry Smith /*@C
417ef20d060SBarry Smith    DMTSGetSolutionFunction - gets the TS solution evaluation function
418ef20d060SBarry Smith 
419ef20d060SBarry Smith    Not Collective
420ef20d060SBarry Smith 
421ef20d060SBarry Smith    Input Arguments:
422ef20d060SBarry Smith .  dm - DM to be used with TS
423ef20d060SBarry Smith 
424ef20d060SBarry Smith    Output Parameters:
425ef20d060SBarry Smith +  func - solution function evaluation function, see TSSetSolution() for calling sequence
426ef20d060SBarry Smith -  ctx - context for solution evaluation
427ef20d060SBarry Smith 
428ef20d060SBarry Smith    Level: advanced
429ef20d060SBarry Smith 
430ef20d060SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
431ef20d060SBarry Smith @*/
432ef20d060SBarry Smith PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
433ef20d060SBarry Smith {
434ef20d060SBarry Smith   PetscErrorCode ierr;
435942e3340SBarry Smith   DMTS           tsdm;
436ef20d060SBarry Smith 
437ef20d060SBarry Smith   PetscFunctionBegin;
438ef20d060SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
439942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
440d74926cbSBarry Smith   if (func) *func = tsdm->ops->solution;
441ef20d060SBarry Smith   if (ctx)  *ctx  = tsdm->solutionctx;
442ef20d060SBarry Smith   PetscFunctionReturn(0);
443ef20d060SBarry Smith }
444ef20d060SBarry Smith 
445ef20d060SBarry Smith #undef __FUNCT__
446ef20d060SBarry Smith #define __FUNCT__ "DMTSSetSolutionFunction"
447ef20d060SBarry Smith /*@C
448ef20d060SBarry Smith    DMTSSetSolutionFunction - set TS solution evaluation function
449ef20d060SBarry Smith 
450ef20d060SBarry Smith    Not Collective
451ef20d060SBarry Smith 
452ef20d060SBarry Smith    Input Arguments:
453ef20d060SBarry Smith +  dm - DM to be used with TS
454ef20d060SBarry Smith .  func - solution function evaluation function, see TSSetSolution() for calling sequence
455ef20d060SBarry Smith -  ctx - context for solution evaluation
456ef20d060SBarry Smith 
457ef20d060SBarry Smith    Level: advanced
458ef20d060SBarry Smith 
459ef20d060SBarry Smith    Note:
460ef20d060SBarry Smith    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
461ef20d060SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
462ef20d060SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
463ef20d060SBarry Smith 
464ef20d060SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
465ef20d060SBarry Smith @*/
466ef20d060SBarry Smith PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
467ef20d060SBarry Smith {
468ef20d060SBarry Smith   PetscErrorCode ierr;
469942e3340SBarry Smith   DMTS           tsdm;
470ef20d060SBarry Smith 
471ef20d060SBarry Smith   PetscFunctionBegin;
472ef20d060SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
473942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
474d74926cbSBarry Smith   if (func) tsdm->ops->solution = func;
475ef20d060SBarry Smith   if (ctx)  tsdm->solutionctx   = ctx;
476ef20d060SBarry Smith   PetscFunctionReturn(0);
477ef20d060SBarry Smith }
478ef20d060SBarry Smith 
479ef20d060SBarry Smith #undef __FUNCT__
4809b7cd975SBarry Smith #define __FUNCT__ "DMTSSetForcingFunction"
4819b7cd975SBarry Smith /*@C
4829b7cd975SBarry Smith    DMTSSetForcingFunction - set TS forcing function evaluation function
4839b7cd975SBarry Smith 
4849b7cd975SBarry Smith    Not Collective
4859b7cd975SBarry Smith 
4869b7cd975SBarry Smith    Input Arguments:
4879b7cd975SBarry Smith +  dm - DM to be used with TS
4889b7cd975SBarry Smith .  TSForcingFunction - forcing function evaluation function
4899b7cd975SBarry Smith -  ctx - context for solution evaluation
4909b7cd975SBarry Smith 
4919b7cd975SBarry Smith    Level: advanced
4929b7cd975SBarry Smith 
4939b7cd975SBarry Smith    Note:
4949b7cd975SBarry Smith    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
4959b7cd975SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
4969b7cd975SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
4979b7cd975SBarry Smith 
4989b7cd975SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
4999b7cd975SBarry Smith @*/
5009b7cd975SBarry Smith PetscErrorCode DMTSSetForcingFunction(DM dm,PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*),void *ctx)
5019b7cd975SBarry Smith {
5029b7cd975SBarry Smith   PetscErrorCode ierr;
5039b7cd975SBarry Smith   DMTS           tsdm;
5049b7cd975SBarry Smith 
5059b7cd975SBarry Smith   PetscFunctionBegin;
5069b7cd975SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5079b7cd975SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
5089b7cd975SBarry Smith   if (TSForcingFunction) tsdm->ops->forcing = TSForcingFunction;
5099b7cd975SBarry Smith   if (ctx)  tsdm->forcingctx   = ctx;
5109b7cd975SBarry Smith   PetscFunctionReturn(0);
5119b7cd975SBarry Smith }
5129b7cd975SBarry Smith 
5139b7cd975SBarry Smith 
5149b7cd975SBarry Smith #undef __FUNCT__
5159b7cd975SBarry Smith #define __FUNCT__ "DMTSGetForcingFunction"
5169b7cd975SBarry Smith /*@C
5179b7cd975SBarry Smith    DMTSGetForcingFunction - get TS forcing function evaluation function
5189b7cd975SBarry Smith 
5199b7cd975SBarry Smith    Not Collective
5209b7cd975SBarry Smith 
5219b7cd975SBarry Smith    Input Argument:
5229b7cd975SBarry Smith .   dm - DM to be used with TS
5239b7cd975SBarry Smith 
5249b7cd975SBarry Smith    Output Arguments:
5259b7cd975SBarry Smith +  TSForcingFunction - forcing function evaluation function
5269b7cd975SBarry Smith -  ctx - context for solution evaluation
5279b7cd975SBarry Smith 
5289b7cd975SBarry Smith    Level: advanced
5299b7cd975SBarry Smith 
5309b7cd975SBarry Smith    Note:
5319b7cd975SBarry Smith    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
5329b7cd975SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
5339b7cd975SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
5349b7cd975SBarry Smith 
5359b7cd975SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
5369b7cd975SBarry Smith @*/
5379b7cd975SBarry Smith PetscErrorCode DMTSGetForcingFunction(DM dm,PetscErrorCode (**TSForcingFunction)(TS,PetscReal,Vec,void*),void **ctx)
5389b7cd975SBarry Smith {
5399b7cd975SBarry Smith   PetscErrorCode ierr;
5409b7cd975SBarry Smith   DMTS           tsdm;
5419b7cd975SBarry Smith 
5429b7cd975SBarry Smith   PetscFunctionBegin;
5439b7cd975SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5449b7cd975SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
5459b7cd975SBarry Smith   if (TSForcingFunction) *TSForcingFunction = tsdm->ops->forcing;
5469b7cd975SBarry Smith   if (ctx) *ctx = tsdm->forcingctx;
5479b7cd975SBarry Smith   PetscFunctionReturn(0);
5489b7cd975SBarry Smith }
5499b7cd975SBarry Smith 
5509b7cd975SBarry Smith #undef __FUNCT__
55124989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSFunction"
55224989b8cSPeter Brune /*@C
55324989b8cSPeter Brune    DMTSGetRHSFunction - get TS explicit residual evaluation function
55424989b8cSPeter Brune 
55524989b8cSPeter Brune    Not Collective
55624989b8cSPeter Brune 
55724989b8cSPeter Brune    Input Argument:
55824989b8cSPeter Brune .  dm - DM to be used with TS
55924989b8cSPeter Brune 
56024989b8cSPeter Brune    Output Arguments:
56124989b8cSPeter Brune +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
56224989b8cSPeter Brune -  ctx - context for residual evaluation
56324989b8cSPeter Brune 
56424989b8cSPeter Brune    Level: advanced
56524989b8cSPeter Brune 
56624989b8cSPeter Brune    Note:
56724989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
56824989b8cSPeter Brune    associated with the DM.
56924989b8cSPeter Brune 
57024989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
57124989b8cSPeter Brune @*/
57224989b8cSPeter Brune PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
57324989b8cSPeter Brune {
57424989b8cSPeter Brune   PetscErrorCode ierr;
575942e3340SBarry Smith   DMTS           tsdm;
57624989b8cSPeter Brune 
57724989b8cSPeter Brune   PetscFunctionBegin;
57824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
579942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
580d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsfunction;
58124989b8cSPeter Brune   if (ctx)  *ctx = tsdm->rhsfunctionctx;
58224989b8cSPeter Brune   PetscFunctionReturn(0);
58324989b8cSPeter Brune }
58424989b8cSPeter Brune 
58524989b8cSPeter Brune #undef __FUNCT__
58624989b8cSPeter Brune #define __FUNCT__ "DMTSSetIJacobian"
58724989b8cSPeter Brune /*@C
58824989b8cSPeter Brune    DMTSSetIJacobian - set TS Jacobian evaluation function
58924989b8cSPeter Brune 
59024989b8cSPeter Brune    Not Collective
59124989b8cSPeter Brune 
59224989b8cSPeter Brune    Input Argument:
59324989b8cSPeter Brune +  dm - DM to be used with TS
59424989b8cSPeter Brune .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
59524989b8cSPeter Brune -  ctx - context for residual evaluation
59624989b8cSPeter Brune 
59724989b8cSPeter Brune    Level: advanced
59824989b8cSPeter Brune 
59924989b8cSPeter Brune    Note:
60024989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
60124989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
60224989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
60324989b8cSPeter Brune 
60424989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
60524989b8cSPeter Brune @*/
60624989b8cSPeter Brune PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
60724989b8cSPeter Brune {
60824989b8cSPeter Brune   PetscErrorCode ierr;
609942e3340SBarry Smith   DMTS           sdm;
61024989b8cSPeter Brune 
61124989b8cSPeter Brune   PetscFunctionBegin;
61224989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
613942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
614d74926cbSBarry Smith   if (func) sdm->ops->ijacobian = func;
61524989b8cSPeter Brune   if (ctx)  sdm->ijacobianctx   = ctx;
61624989b8cSPeter Brune   PetscFunctionReturn(0);
61724989b8cSPeter Brune }
61824989b8cSPeter Brune 
61924989b8cSPeter Brune #undef __FUNCT__
62024989b8cSPeter Brune #define __FUNCT__ "DMTSGetIJacobian"
62124989b8cSPeter Brune /*@C
62224989b8cSPeter Brune    DMTSGetIJacobian - get TS Jacobian evaluation function
62324989b8cSPeter Brune 
62424989b8cSPeter Brune    Not Collective
62524989b8cSPeter Brune 
62624989b8cSPeter Brune    Input Argument:
62724989b8cSPeter Brune .  dm - DM to be used with TS
62824989b8cSPeter Brune 
62924989b8cSPeter Brune    Output Arguments:
63024989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
63124989b8cSPeter Brune -  ctx - context for residual evaluation
63224989b8cSPeter Brune 
63324989b8cSPeter Brune    Level: advanced
63424989b8cSPeter Brune 
63524989b8cSPeter Brune    Note:
63624989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
63724989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
63824989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
63924989b8cSPeter Brune 
64024989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
64124989b8cSPeter Brune @*/
64224989b8cSPeter Brune PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
64324989b8cSPeter Brune {
64424989b8cSPeter Brune   PetscErrorCode ierr;
645942e3340SBarry Smith   DMTS           tsdm;
64624989b8cSPeter Brune 
64724989b8cSPeter Brune   PetscFunctionBegin;
64824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
649942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
650d74926cbSBarry Smith   if (func) *func = tsdm->ops->ijacobian;
65124989b8cSPeter Brune   if (ctx)  *ctx = tsdm->ijacobianctx;
65224989b8cSPeter Brune   PetscFunctionReturn(0);
65324989b8cSPeter Brune }
65424989b8cSPeter Brune 
65524989b8cSPeter Brune 
65624989b8cSPeter Brune #undef __FUNCT__
65724989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSJacobian"
65824989b8cSPeter Brune /*@C
65924989b8cSPeter Brune    DMTSSetRHSJacobian - set TS Jacobian evaluation function
66024989b8cSPeter Brune 
66124989b8cSPeter Brune    Not Collective
66224989b8cSPeter Brune 
66324989b8cSPeter Brune    Input Argument:
66424989b8cSPeter Brune +  dm - DM to be used with TS
66524989b8cSPeter Brune .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
66624989b8cSPeter Brune -  ctx - context for residual evaluation
66724989b8cSPeter Brune 
66824989b8cSPeter Brune    Level: advanced
66924989b8cSPeter Brune 
67024989b8cSPeter Brune    Note:
67124989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
67224989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
67324989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
67424989b8cSPeter Brune 
67524989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
67624989b8cSPeter Brune @*/
67724989b8cSPeter Brune PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
67824989b8cSPeter Brune {
67924989b8cSPeter Brune   PetscErrorCode ierr;
680942e3340SBarry Smith   DMTS           tsdm;
68124989b8cSPeter Brune 
68224989b8cSPeter Brune   PetscFunctionBegin;
68324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
684942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
685d74926cbSBarry Smith   if (func) tsdm->ops->rhsjacobian = func;
68624989b8cSPeter Brune   if (ctx)  tsdm->rhsjacobianctx = ctx;
68724989b8cSPeter Brune   PetscFunctionReturn(0);
68824989b8cSPeter Brune }
68924989b8cSPeter Brune 
69024989b8cSPeter Brune #undef __FUNCT__
69124989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSJacobian"
69224989b8cSPeter Brune /*@C
69324989b8cSPeter Brune    DMTSGetRHSJacobian - get TS Jacobian evaluation function
69424989b8cSPeter Brune 
69524989b8cSPeter Brune    Not Collective
69624989b8cSPeter Brune 
69724989b8cSPeter Brune    Input Argument:
69824989b8cSPeter Brune .  dm - DM to be used with TS
69924989b8cSPeter Brune 
70024989b8cSPeter Brune    Output Arguments:
70124989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
70224989b8cSPeter Brune -  ctx - context for residual evaluation
70324989b8cSPeter Brune 
70424989b8cSPeter Brune    Level: advanced
70524989b8cSPeter Brune 
70624989b8cSPeter Brune    Note:
70724989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
70824989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
70924989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
71024989b8cSPeter Brune 
71124989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
71224989b8cSPeter Brune @*/
71324989b8cSPeter Brune PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
71424989b8cSPeter Brune {
71524989b8cSPeter Brune   PetscErrorCode ierr;
716942e3340SBarry Smith   DMTS           tsdm;
71724989b8cSPeter Brune 
71824989b8cSPeter Brune   PetscFunctionBegin;
71924989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
720942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
721d74926cbSBarry Smith   if (func) *func = tsdm->ops->rhsjacobian;
72224989b8cSPeter Brune   if (ctx)  *ctx = tsdm->rhsjacobianctx;
72324989b8cSPeter Brune   PetscFunctionReturn(0);
72424989b8cSPeter Brune }
725ad6bc421SBarry Smith 
726ad6bc421SBarry Smith #undef __FUNCT__
727ad6bc421SBarry Smith #define __FUNCT__ "DMTSSetIFunctionSerialize"
728ad6bc421SBarry Smith /*@C
729ad6bc421SBarry Smith    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
730ad6bc421SBarry Smith 
731ad6bc421SBarry Smith    Not Collective
732ad6bc421SBarry Smith 
733ad6bc421SBarry Smith    Input Arguments:
734ad6bc421SBarry Smith +  dm - DM to be used with TS
735ad6bc421SBarry Smith .  view - viewer function
736ad6bc421SBarry Smith -  load - loading function
737ad6bc421SBarry Smith 
738ad6bc421SBarry Smith    Level: advanced
739ad6bc421SBarry Smith 
740ad6bc421SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
741ad6bc421SBarry Smith @*/
742ad6bc421SBarry Smith PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
743ad6bc421SBarry Smith {
744ad6bc421SBarry Smith   PetscErrorCode ierr;
745ad6bc421SBarry Smith   DMTS           tsdm;
746ad6bc421SBarry Smith 
747ad6bc421SBarry Smith   PetscFunctionBegin;
748ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
749ad6bc421SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
750ad6bc421SBarry Smith   tsdm->ops->ifunctionview = view;
751ad6bc421SBarry Smith   tsdm->ops->ifunctionload = load;
752ad6bc421SBarry Smith   PetscFunctionReturn(0);
753ad6bc421SBarry Smith }
754ad6bc421SBarry Smith 
755ad6bc421SBarry Smith #undef __FUNCT__
756ad6bc421SBarry Smith #define __FUNCT__ "DMTSSetIJacobianSerialize"
757ad6bc421SBarry Smith /*@C
758ad6bc421SBarry Smith    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
759ad6bc421SBarry Smith 
760ad6bc421SBarry Smith    Not Collective
761ad6bc421SBarry Smith 
762ad6bc421SBarry Smith    Input Arguments:
763ad6bc421SBarry Smith +  dm - DM to be used with TS
764ad6bc421SBarry Smith .  view - viewer function
765ad6bc421SBarry Smith -  load - loading function
766ad6bc421SBarry Smith 
767ad6bc421SBarry Smith    Level: advanced
768ad6bc421SBarry Smith 
769ad6bc421SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
770ad6bc421SBarry Smith @*/
771ad6bc421SBarry Smith PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
772ad6bc421SBarry Smith {
773ad6bc421SBarry Smith   PetscErrorCode ierr;
774ad6bc421SBarry Smith   DMTS           tsdm;
775ad6bc421SBarry Smith 
776ad6bc421SBarry Smith   PetscFunctionBegin;
777ad6bc421SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
778ad6bc421SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
779ad6bc421SBarry Smith   tsdm->ops->ijacobianview = view;
780ad6bc421SBarry Smith   tsdm->ops->ijacobianload = load;
781ad6bc421SBarry Smith   PetscFunctionReturn(0);
782ad6bc421SBarry Smith }
783