xref: /petsc/src/ts/utils/dmts.c (revision 942e334039a1a47ff909c032afdf0581bd880e9f)
124989b8cSPeter Brune #include <petsc-private/tsimpl.h>     /*I "petscts.h" I*/
224989b8cSPeter Brune #include <petsc-private/dmimpl.h>     /*I "petscdm.h" I*/
324989b8cSPeter Brune 
424989b8cSPeter Brune 
524989b8cSPeter Brune #undef __FUNCT__
6*942e3340SBarry Smith #define __FUNCT__ "DMCoarsenHook_DMTS"
7*942e3340SBarry Smith /* Attaches the DMSNES to the coarse level.
824989b8cSPeter Brune  * Under what conditions should we copy versus duplicate?
924989b8cSPeter Brune  */
10*942e3340SBarry Smith static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
1124989b8cSPeter Brune {
1224989b8cSPeter Brune   PetscErrorCode ierr;
1324989b8cSPeter Brune 
1424989b8cSPeter Brune   PetscFunctionBegin;
15*942e3340SBarry Smith   ierr = DMCopyDMTS(dm,dmc);CHKERRQ(ierr);
1624989b8cSPeter Brune   PetscFunctionReturn(0);
1724989b8cSPeter Brune }
1824989b8cSPeter Brune 
1924989b8cSPeter Brune #undef __FUNCT__
20*942e3340SBarry Smith #define __FUNCT__ "DMRestrictHook_DMTS"
2124989b8cSPeter Brune /* This could restrict auxiliary information to the coarse level.
2224989b8cSPeter Brune  */
23*942e3340SBarry Smith static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
2424989b8cSPeter Brune {
2524989b8cSPeter Brune 
2624989b8cSPeter Brune   PetscFunctionBegin;
2724989b8cSPeter Brune   PetscFunctionReturn(0);
2824989b8cSPeter Brune }
2924989b8cSPeter Brune 
3024989b8cSPeter Brune 
3124989b8cSPeter Brune #undef __FUNCT__
32*942e3340SBarry Smith #define __FUNCT__ "PetscContainerDestroy_DMTS"
33*942e3340SBarry Smith static PetscErrorCode PetscContainerDestroy_DMTS(void *ctx)
3424989b8cSPeter Brune {
3524989b8cSPeter Brune   PetscErrorCode ierr;
36*942e3340SBarry Smith   DMTS           tsdm = (DMTS)ctx;
3724989b8cSPeter Brune 
3824989b8cSPeter Brune   PetscFunctionBegin;
3924989b8cSPeter Brune   if (tsdm->destroy) {ierr = (*tsdm->destroy)(tsdm);CHKERRQ(ierr);}
4024989b8cSPeter Brune   ierr = PetscFree(tsdm);CHKERRQ(ierr);
4124989b8cSPeter Brune   PetscFunctionReturn(0);
4224989b8cSPeter Brune }
4324989b8cSPeter Brune 
4424989b8cSPeter Brune #undef __FUNCT__
45*942e3340SBarry Smith #define __FUNCT__ "DMGetDMTS"
4624989b8cSPeter Brune /*@C
47*942e3340SBarry Smith    DMGetDMTS - get read-only private DMTS context from a DM
4824989b8cSPeter Brune 
4924989b8cSPeter Brune    Not Collective
5024989b8cSPeter Brune 
5124989b8cSPeter Brune    Input Argument:
5224989b8cSPeter Brune .  dm - DM to be used with TS
5324989b8cSPeter Brune 
5424989b8cSPeter Brune    Output Argument:
55*942e3340SBarry Smith .  tsdm - private DMTS context
5624989b8cSPeter Brune 
5724989b8cSPeter Brune    Level: developer
5824989b8cSPeter Brune 
5924989b8cSPeter Brune    Notes:
60*942e3340SBarry Smith    Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.
6124989b8cSPeter Brune 
62*942e3340SBarry Smith .seealso: DMGetDMTSWrite()
6324989b8cSPeter Brune @*/
64*942e3340SBarry Smith PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
6524989b8cSPeter Brune {
6624989b8cSPeter Brune   PetscErrorCode ierr;
6724989b8cSPeter Brune   PetscContainer container;
68*942e3340SBarry Smith   DMTS           tsdmnew;
6924989b8cSPeter Brune 
7024989b8cSPeter Brune 
7124989b8cSPeter Brune   PetscFunctionBegin;
7224989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
73*942e3340SBarry Smith   ierr = PetscObjectQuery((PetscObject)dm,"DMTS",(PetscObject*)&container);CHKERRQ(ierr);
7424989b8cSPeter Brune   if (container) {
7524989b8cSPeter Brune     ierr = PetscContainerGetPointer(container,(void**)tsdm);CHKERRQ(ierr);
7624989b8cSPeter Brune   } else {
77*942e3340SBarry Smith     ierr = PetscInfo(dm,"Creating new DMTS\n");CHKERRQ(ierr);
7824989b8cSPeter Brune     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
79*942e3340SBarry Smith     ierr = PetscNewLog(dm,struct _n_DMTS,&tsdmnew);CHKERRQ(ierr);
8024989b8cSPeter Brune     ierr = PetscContainerSetPointer(container,tsdmnew);CHKERRQ(ierr);
81*942e3340SBarry Smith     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_DMTS);CHKERRQ(ierr);
82*942e3340SBarry Smith     ierr = PetscObjectCompose((PetscObject)dm,"DMTS",(PetscObject)container);CHKERRQ(ierr);
83*942e3340SBarry Smith     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr);
8424989b8cSPeter Brune     ierr = PetscContainerGetPointer(container,(void**)tsdm);CHKERRQ(ierr);
8524989b8cSPeter Brune     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
8624989b8cSPeter Brune   }
8724989b8cSPeter Brune   PetscFunctionReturn(0);
8824989b8cSPeter Brune }
8924989b8cSPeter Brune 
9024989b8cSPeter Brune #undef __FUNCT__
91*942e3340SBarry Smith #define __FUNCT__ "DMGetDMTSWrite"
9224989b8cSPeter Brune /*@C
93*942e3340SBarry Smith    DMGetDMTSWrite - get write access to private DMTS context from a DM
9424989b8cSPeter Brune 
9524989b8cSPeter Brune    Not Collective
9624989b8cSPeter Brune 
9724989b8cSPeter Brune    Input Argument:
9824989b8cSPeter Brune .  dm - DM to be used with TS
9924989b8cSPeter Brune 
10024989b8cSPeter Brune    Output Argument:
101*942e3340SBarry Smith .  tsdm - private DMTS context
10224989b8cSPeter Brune 
10324989b8cSPeter Brune    Level: developer
10424989b8cSPeter Brune 
105*942e3340SBarry Smith .seealso: DMGetDMTS()
10624989b8cSPeter Brune @*/
107*942e3340SBarry Smith PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
10824989b8cSPeter Brune {
10924989b8cSPeter Brune   PetscErrorCode ierr;
110*942e3340SBarry Smith   DMTS           sdm;
11124989b8cSPeter Brune 
11224989b8cSPeter Brune   PetscFunctionBegin;
11324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
114*942e3340SBarry Smith   ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr);
11524989b8cSPeter Brune   if (!sdm->originaldm) sdm->originaldm = dm;
11624989b8cSPeter Brune   if (sdm->originaldm != dm) {  /* Copy on write */
11724989b8cSPeter Brune     PetscContainer container;
118*942e3340SBarry Smith     DMTS           oldsdm = sdm;
119*942e3340SBarry Smith     ierr = PetscInfo(dm,"Copying DMTS due to write\n");CHKERRQ(ierr);
12024989b8cSPeter Brune     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
121*942e3340SBarry Smith     ierr = PetscNewLog(dm,struct _n_DMTS,&sdm);CHKERRQ(ierr);
1228caf3d72SBarry Smith     ierr = PetscMemcpy(sdm,oldsdm,sizeof(*sdm));CHKERRQ(ierr);
12324989b8cSPeter Brune     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
124*942e3340SBarry Smith     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_DMTS);CHKERRQ(ierr);
125*942e3340SBarry Smith     ierr = PetscObjectCompose((PetscObject)dm,"DMTS",(PetscObject)container);CHKERRQ(ierr);
12624989b8cSPeter Brune     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
1270313d0e0SJed Brown     /* implementation specific copy hooks */
128c9058dfbSJed Brown     if (sdm->duplicate) {ierr = (*sdm->duplicate)(oldsdm,dm);CHKERRQ(ierr);}
12924989b8cSPeter Brune   }
13024989b8cSPeter Brune   *tsdm = sdm;
13124989b8cSPeter Brune   PetscFunctionReturn(0);
13224989b8cSPeter Brune }
13324989b8cSPeter Brune 
13424989b8cSPeter Brune #undef __FUNCT__
135*942e3340SBarry Smith #define __FUNCT__ "DMCopyDMTS"
13624989b8cSPeter Brune /*@C
137*942e3340SBarry Smith    DMCopyDMTS - copies a DM context to a new DM
13824989b8cSPeter Brune 
13924989b8cSPeter Brune    Logically Collective
14024989b8cSPeter Brune 
14124989b8cSPeter Brune    Input Arguments:
14224989b8cSPeter Brune +  dmsrc - DM to obtain context from
14324989b8cSPeter Brune -  dmdest - DM to add context to
14424989b8cSPeter Brune 
14524989b8cSPeter Brune    Level: developer
14624989b8cSPeter Brune 
14724989b8cSPeter Brune    Note:
14824989b8cSPeter Brune    The context is copied by reference. This function does not ensure that a context exists.
14924989b8cSPeter Brune 
150*942e3340SBarry Smith .seealso: DMGetDMTS(), TSSetDM()
15124989b8cSPeter Brune @*/
152*942e3340SBarry Smith PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
15324989b8cSPeter Brune {
15424989b8cSPeter Brune   PetscErrorCode ierr;
15524989b8cSPeter Brune   PetscContainer container;
15624989b8cSPeter Brune 
15724989b8cSPeter Brune   PetscFunctionBegin;
15824989b8cSPeter Brune   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
15924989b8cSPeter Brune   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
160*942e3340SBarry Smith   ierr = PetscObjectQuery((PetscObject)dmsrc,"DMTS",(PetscObject*)&container);CHKERRQ(ierr);
16124989b8cSPeter Brune   if (container) {
162*942e3340SBarry Smith     ierr = PetscObjectCompose((PetscObject)dmdest,"DMTS",(PetscObject)container);CHKERRQ(ierr);
163*942e3340SBarry Smith     ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr);
16424989b8cSPeter Brune   }
16524989b8cSPeter Brune   PetscFunctionReturn(0);
16624989b8cSPeter Brune }
16724989b8cSPeter Brune 
16824989b8cSPeter Brune #undef __FUNCT__
16924989b8cSPeter Brune #define __FUNCT__ "DMTSSetIFunction"
17024989b8cSPeter Brune /*@C
17124989b8cSPeter Brune    DMTSSetIFunction - set TS implicit function evaluation function
17224989b8cSPeter Brune 
17324989b8cSPeter Brune    Not Collective
17424989b8cSPeter Brune 
17524989b8cSPeter Brune    Input Arguments:
17624989b8cSPeter Brune +  dm - DM to be used with TS
17724989b8cSPeter Brune .  func - function evaluation function, see TSSetIFunction() for calling sequence
17824989b8cSPeter Brune -  ctx - context for residual evaluation
17924989b8cSPeter Brune 
18024989b8cSPeter Brune    Level: advanced
18124989b8cSPeter Brune 
18224989b8cSPeter Brune    Note:
18324989b8cSPeter Brune    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
18424989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
18524989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
18624989b8cSPeter Brune 
18724989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
18824989b8cSPeter Brune @*/
18924989b8cSPeter Brune PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
19024989b8cSPeter Brune {
19124989b8cSPeter Brune   PetscErrorCode ierr;
192*942e3340SBarry Smith   DMTS           tsdm;
19324989b8cSPeter Brune 
19424989b8cSPeter Brune   PetscFunctionBegin;
19524989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
196*942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
19724989b8cSPeter Brune   if (func) tsdm->ifunction = func;
19824989b8cSPeter Brune   if (ctx)  tsdm->ifunctionctx = ctx;
19924989b8cSPeter Brune   PetscFunctionReturn(0);
20024989b8cSPeter Brune }
20124989b8cSPeter Brune 
20224989b8cSPeter Brune #undef __FUNCT__
20324989b8cSPeter Brune #define __FUNCT__ "DMTSGetIFunction"
20424989b8cSPeter Brune /*@C
20524989b8cSPeter Brune    DMTSGetIFunction - get TS implicit residual evaluation function
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 Arguments:
21324989b8cSPeter Brune +  func - function evaluation function, see TSSetIFunction() for calling sequence
21424989b8cSPeter Brune -  ctx - context for residual evaluation
21524989b8cSPeter Brune 
21624989b8cSPeter Brune    Level: advanced
21724989b8cSPeter Brune 
21824989b8cSPeter Brune    Note:
21924989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
22024989b8cSPeter Brune    associated with the DM.
22124989b8cSPeter Brune 
22224989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
22324989b8cSPeter Brune @*/
22424989b8cSPeter Brune PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
22524989b8cSPeter Brune {
22624989b8cSPeter Brune   PetscErrorCode ierr;
227*942e3340SBarry Smith   DMTS           tsdm;
22824989b8cSPeter Brune 
22924989b8cSPeter Brune   PetscFunctionBegin;
23024989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
231*942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
23224989b8cSPeter Brune   if (func) *func = tsdm->ifunction;
23324989b8cSPeter Brune   if (ctx)  *ctx = tsdm->ifunctionctx;
23424989b8cSPeter Brune   PetscFunctionReturn(0);
23524989b8cSPeter Brune }
23624989b8cSPeter Brune 
23724989b8cSPeter Brune 
23824989b8cSPeter Brune #undef __FUNCT__
23924989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSFunction"
24024989b8cSPeter Brune /*@C
24124989b8cSPeter Brune    DMTSSetRHSFunction - set TS explicit residual evaluation function
24224989b8cSPeter Brune 
24324989b8cSPeter Brune    Not Collective
24424989b8cSPeter Brune 
24524989b8cSPeter Brune    Input Arguments:
24624989b8cSPeter Brune +  dm - DM to be used with TS
24724989b8cSPeter Brune .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
24824989b8cSPeter Brune -  ctx - context for residual evaluation
24924989b8cSPeter Brune 
25024989b8cSPeter Brune    Level: advanced
25124989b8cSPeter Brune 
25224989b8cSPeter Brune    Note:
253ef20d060SBarry Smith    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
25424989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
25524989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
25624989b8cSPeter Brune 
25724989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
25824989b8cSPeter Brune @*/
25924989b8cSPeter Brune PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
26024989b8cSPeter Brune {
26124989b8cSPeter Brune   PetscErrorCode ierr;
262*942e3340SBarry Smith   DMTS           tsdm;
26324989b8cSPeter Brune 
26424989b8cSPeter Brune   PetscFunctionBegin;
26524989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
266*942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
26724989b8cSPeter Brune   if (func) tsdm->rhsfunction = func;
26824989b8cSPeter Brune   if (ctx)  tsdm->rhsfunctionctx = ctx;
26924989b8cSPeter Brune   PetscFunctionReturn(0);
27024989b8cSPeter Brune }
27124989b8cSPeter Brune 
27224989b8cSPeter Brune #undef __FUNCT__
273ef20d060SBarry Smith #define __FUNCT__ "DMTSGetSolutionFunction"
274ef20d060SBarry Smith /*@C
275ef20d060SBarry Smith    DMTSGetSolutionFunction - gets the TS solution evaluation function
276ef20d060SBarry Smith 
277ef20d060SBarry Smith    Not Collective
278ef20d060SBarry Smith 
279ef20d060SBarry Smith    Input Arguments:
280ef20d060SBarry Smith .  dm - DM to be used with TS
281ef20d060SBarry Smith 
282ef20d060SBarry Smith    Output Parameters:
283ef20d060SBarry Smith +  func - solution function evaluation function, see TSSetSolution() for calling sequence
284ef20d060SBarry Smith -  ctx - context for solution evaluation
285ef20d060SBarry Smith 
286ef20d060SBarry Smith    Level: advanced
287ef20d060SBarry Smith 
288ef20d060SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
289ef20d060SBarry Smith @*/
290ef20d060SBarry Smith PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
291ef20d060SBarry Smith {
292ef20d060SBarry Smith   PetscErrorCode ierr;
293*942e3340SBarry Smith   DMTS           tsdm;
294ef20d060SBarry Smith 
295ef20d060SBarry Smith   PetscFunctionBegin;
296ef20d060SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
297*942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
298ef20d060SBarry Smith   if (func) *func = tsdm->solution;
299ef20d060SBarry Smith   if (ctx)  *ctx  = tsdm->solutionctx;
300ef20d060SBarry Smith   PetscFunctionReturn(0);
301ef20d060SBarry Smith }
302ef20d060SBarry Smith 
303ef20d060SBarry Smith #undef __FUNCT__
304ef20d060SBarry Smith #define __FUNCT__ "DMTSSetSolutionFunction"
305ef20d060SBarry Smith /*@C
306ef20d060SBarry Smith    DMTSSetSolutionFunction - set TS solution evaluation function
307ef20d060SBarry Smith 
308ef20d060SBarry Smith    Not Collective
309ef20d060SBarry Smith 
310ef20d060SBarry Smith    Input Arguments:
311ef20d060SBarry Smith +  dm - DM to be used with TS
312ef20d060SBarry Smith .  func - solution function evaluation function, see TSSetSolution() for calling sequence
313ef20d060SBarry Smith -  ctx - context for solution evaluation
314ef20d060SBarry Smith 
315ef20d060SBarry Smith    Level: advanced
316ef20d060SBarry Smith 
317ef20d060SBarry Smith    Note:
318ef20d060SBarry Smith    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
319ef20d060SBarry Smith    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
320ef20d060SBarry Smith    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
321ef20d060SBarry Smith 
322ef20d060SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
323ef20d060SBarry Smith @*/
324ef20d060SBarry Smith PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
325ef20d060SBarry Smith {
326ef20d060SBarry Smith   PetscErrorCode ierr;
327*942e3340SBarry Smith   DMTS           tsdm;
328ef20d060SBarry Smith 
329ef20d060SBarry Smith   PetscFunctionBegin;
330ef20d060SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
331*942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
332ef20d060SBarry Smith   if (func) tsdm->solution    = func;
333ef20d060SBarry Smith   if (ctx)  tsdm->solutionctx = ctx;
334ef20d060SBarry Smith   PetscFunctionReturn(0);
335ef20d060SBarry Smith }
336ef20d060SBarry Smith 
337ef20d060SBarry Smith #undef __FUNCT__
33824989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSFunction"
33924989b8cSPeter Brune /*@C
34024989b8cSPeter Brune    DMTSGetRHSFunction - get TS explicit residual evaluation function
34124989b8cSPeter Brune 
34224989b8cSPeter Brune    Not Collective
34324989b8cSPeter Brune 
34424989b8cSPeter Brune    Input Argument:
34524989b8cSPeter Brune .  dm - DM to be used with TS
34624989b8cSPeter Brune 
34724989b8cSPeter Brune    Output Arguments:
34824989b8cSPeter Brune +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
34924989b8cSPeter Brune -  ctx - context for residual evaluation
35024989b8cSPeter Brune 
35124989b8cSPeter Brune    Level: advanced
35224989b8cSPeter Brune 
35324989b8cSPeter Brune    Note:
35424989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
35524989b8cSPeter Brune    associated with the DM.
35624989b8cSPeter Brune 
35724989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
35824989b8cSPeter Brune @*/
35924989b8cSPeter Brune PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
36024989b8cSPeter Brune {
36124989b8cSPeter Brune   PetscErrorCode ierr;
362*942e3340SBarry Smith   DMTS           tsdm;
36324989b8cSPeter Brune 
36424989b8cSPeter Brune   PetscFunctionBegin;
36524989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
366*942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
36724989b8cSPeter Brune   if (func) *func = tsdm->rhsfunction;
36824989b8cSPeter Brune   if (ctx)  *ctx = tsdm->rhsfunctionctx;
36924989b8cSPeter Brune   PetscFunctionReturn(0);
37024989b8cSPeter Brune }
37124989b8cSPeter Brune 
37224989b8cSPeter Brune #undef __FUNCT__
37324989b8cSPeter Brune #define __FUNCT__ "DMTSSetIJacobian"
37424989b8cSPeter Brune /*@C
37524989b8cSPeter Brune    DMTSSetIJacobian - set TS Jacobian evaluation function
37624989b8cSPeter Brune 
37724989b8cSPeter Brune    Not Collective
37824989b8cSPeter Brune 
37924989b8cSPeter Brune    Input Argument:
38024989b8cSPeter Brune +  dm - DM to be used with TS
38124989b8cSPeter Brune .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
38224989b8cSPeter Brune -  ctx - context for residual evaluation
38324989b8cSPeter Brune 
38424989b8cSPeter Brune    Level: advanced
38524989b8cSPeter Brune 
38624989b8cSPeter Brune    Note:
38724989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
38824989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
38924989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
39024989b8cSPeter Brune 
39124989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
39224989b8cSPeter Brune @*/
39324989b8cSPeter Brune PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
39424989b8cSPeter Brune {
39524989b8cSPeter Brune   PetscErrorCode ierr;
396*942e3340SBarry Smith   DMTS           sdm;
39724989b8cSPeter Brune 
39824989b8cSPeter Brune   PetscFunctionBegin;
39924989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
400*942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
40124989b8cSPeter Brune   if (func) sdm->ijacobian = func;
40224989b8cSPeter Brune   if (ctx)  sdm->ijacobianctx = ctx;
40324989b8cSPeter Brune   PetscFunctionReturn(0);
40424989b8cSPeter Brune }
40524989b8cSPeter Brune 
40624989b8cSPeter Brune #undef __FUNCT__
40724989b8cSPeter Brune #define __FUNCT__ "DMTSGetIJacobian"
40824989b8cSPeter Brune /*@C
40924989b8cSPeter Brune    DMTSGetIJacobian - get TS Jacobian evaluation function
41024989b8cSPeter Brune 
41124989b8cSPeter Brune    Not Collective
41224989b8cSPeter Brune 
41324989b8cSPeter Brune    Input Argument:
41424989b8cSPeter Brune .  dm - DM to be used with TS
41524989b8cSPeter Brune 
41624989b8cSPeter Brune    Output Arguments:
41724989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
41824989b8cSPeter Brune -  ctx - context for residual evaluation
41924989b8cSPeter Brune 
42024989b8cSPeter Brune    Level: advanced
42124989b8cSPeter Brune 
42224989b8cSPeter Brune    Note:
42324989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
42424989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
42524989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
42624989b8cSPeter Brune 
42724989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
42824989b8cSPeter Brune @*/
42924989b8cSPeter Brune PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
43024989b8cSPeter Brune {
43124989b8cSPeter Brune   PetscErrorCode ierr;
432*942e3340SBarry Smith   DMTS           tsdm;
43324989b8cSPeter Brune 
43424989b8cSPeter Brune   PetscFunctionBegin;
43524989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
436*942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
43724989b8cSPeter Brune   if (func) *func = tsdm->ijacobian;
43824989b8cSPeter Brune   if (ctx)  *ctx = tsdm->ijacobianctx;
43924989b8cSPeter Brune   PetscFunctionReturn(0);
44024989b8cSPeter Brune }
44124989b8cSPeter Brune 
44224989b8cSPeter Brune 
44324989b8cSPeter Brune #undef __FUNCT__
44424989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSJacobian"
44524989b8cSPeter Brune /*@C
44624989b8cSPeter Brune    DMTSSetRHSJacobian - set TS Jacobian evaluation function
44724989b8cSPeter Brune 
44824989b8cSPeter Brune    Not Collective
44924989b8cSPeter Brune 
45024989b8cSPeter Brune    Input Argument:
45124989b8cSPeter Brune +  dm - DM to be used with TS
45224989b8cSPeter Brune .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
45324989b8cSPeter Brune -  ctx - context for residual evaluation
45424989b8cSPeter Brune 
45524989b8cSPeter Brune    Level: advanced
45624989b8cSPeter Brune 
45724989b8cSPeter Brune    Note:
45824989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
45924989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
46024989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
46124989b8cSPeter Brune 
46224989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
46324989b8cSPeter Brune @*/
46424989b8cSPeter Brune PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
46524989b8cSPeter Brune {
46624989b8cSPeter Brune   PetscErrorCode ierr;
467*942e3340SBarry Smith   DMTS           tsdm;
46824989b8cSPeter Brune 
46924989b8cSPeter Brune   PetscFunctionBegin;
47024989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
471*942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
47224989b8cSPeter Brune   if (func) tsdm->rhsjacobian = func;
47324989b8cSPeter Brune   if (ctx)  tsdm->rhsjacobianctx = ctx;
47424989b8cSPeter Brune   PetscFunctionReturn(0);
47524989b8cSPeter Brune }
47624989b8cSPeter Brune 
47724989b8cSPeter Brune #undef __FUNCT__
47824989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSJacobian"
47924989b8cSPeter Brune /*@C
48024989b8cSPeter Brune    DMTSGetRHSJacobian - get TS Jacobian evaluation function
48124989b8cSPeter Brune 
48224989b8cSPeter Brune    Not Collective
48324989b8cSPeter Brune 
48424989b8cSPeter Brune    Input Argument:
48524989b8cSPeter Brune .  dm - DM to be used with TS
48624989b8cSPeter Brune 
48724989b8cSPeter Brune    Output Arguments:
48824989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
48924989b8cSPeter Brune -  ctx - context for residual evaluation
49024989b8cSPeter Brune 
49124989b8cSPeter Brune    Level: advanced
49224989b8cSPeter Brune 
49324989b8cSPeter Brune    Note:
49424989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
49524989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
49624989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
49724989b8cSPeter Brune 
49824989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
49924989b8cSPeter Brune @*/
50024989b8cSPeter Brune PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
50124989b8cSPeter Brune {
50224989b8cSPeter Brune   PetscErrorCode ierr;
503*942e3340SBarry Smith   DMTS           tsdm;
50424989b8cSPeter Brune 
50524989b8cSPeter Brune   PetscFunctionBegin;
50624989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
507*942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
50824989b8cSPeter Brune   if (func) *func = tsdm->rhsjacobian;
50924989b8cSPeter Brune   if (ctx)  *ctx = tsdm->rhsjacobianctx;
51024989b8cSPeter Brune   PetscFunctionReturn(0);
51124989b8cSPeter Brune }
512