xref: /petsc/src/ts/utils/dmts.c (revision 31ad16ddc0dbd1b33ade343d3ea8579aceeea9d4)
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__
624989b8cSPeter Brune #define __FUNCT__ "DMCoarsenHook_TSDM"
724989b8cSPeter Brune /* Attaches the SNESDM to the coarse level.
824989b8cSPeter Brune  * Under what conditions should we copy versus duplicate?
924989b8cSPeter Brune  */
1024989b8cSPeter Brune static PetscErrorCode DMCoarsenHook_TSDM(DM dm,DM dmc,void *ctx)
1124989b8cSPeter Brune {
1224989b8cSPeter Brune   PetscErrorCode ierr;
1324989b8cSPeter Brune 
1424989b8cSPeter Brune   PetscFunctionBegin;
1524989b8cSPeter Brune   ierr = DMTSCopyContext(dm,dmc);CHKERRQ(ierr);
1624989b8cSPeter Brune   PetscFunctionReturn(0);
1724989b8cSPeter Brune }
1824989b8cSPeter Brune 
1924989b8cSPeter Brune #undef __FUNCT__
20*31ad16ddSPeter Brune #define __FUNCT__ "DMRestrictHook_TSDM"
2124989b8cSPeter Brune /* This could restrict auxiliary information to the coarse level.
2224989b8cSPeter Brune  */
2324989b8cSPeter Brune static PetscErrorCode DMRestrictHook_TSDM(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__
3224989b8cSPeter Brune #define __FUNCT__ "PetscContainerDestroy_TSDM"
3324989b8cSPeter Brune static PetscErrorCode PetscContainerDestroy_TSDM(void *ctx)
3424989b8cSPeter Brune {
3524989b8cSPeter Brune   PetscErrorCode ierr;
3624989b8cSPeter Brune   TSDM tsdm = (TSDM)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__
4524989b8cSPeter Brune #define __FUNCT__ "DMTSGetContext"
4624989b8cSPeter Brune /*@C
4724989b8cSPeter Brune    DMTSGetContext - get read-only private TSDM 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:
5524989b8cSPeter Brune .  tsdm - private TSDM context
5624989b8cSPeter Brune 
5724989b8cSPeter Brune    Level: developer
5824989b8cSPeter Brune 
5924989b8cSPeter Brune    Notes:
6024989b8cSPeter Brune    Use DMTSGetContextWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.
6124989b8cSPeter Brune 
6224989b8cSPeter Brune .seealso: DMTSGetContextWrite()
6324989b8cSPeter Brune @*/
6424989b8cSPeter Brune PetscErrorCode DMTSGetContext(DM dm,TSDM *tsdm)
6524989b8cSPeter Brune {
6624989b8cSPeter Brune   PetscErrorCode ierr;
6724989b8cSPeter Brune   PetscContainer container;
6824989b8cSPeter Brune   TSDM           tsdmnew;
6924989b8cSPeter Brune 
7024989b8cSPeter Brune 
7124989b8cSPeter Brune   PetscFunctionBegin;
7224989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7324989b8cSPeter Brune   ierr = PetscObjectQuery((PetscObject)dm,"TSDM",(PetscObject*)&container);CHKERRQ(ierr);
7424989b8cSPeter Brune   if (container) {
7524989b8cSPeter Brune     ierr = PetscContainerGetPointer(container,(void**)tsdm);CHKERRQ(ierr);
7624989b8cSPeter Brune   } else {
7724989b8cSPeter Brune     ierr = PetscInfo(dm,"Creating new TSDM\n");CHKERRQ(ierr);
7824989b8cSPeter Brune     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
7924989b8cSPeter Brune     ierr = PetscNewLog(dm,struct _n_TSDM,&tsdmnew);CHKERRQ(ierr);
8024989b8cSPeter Brune     ierr = PetscContainerSetPointer(container,tsdmnew);CHKERRQ(ierr);
8124989b8cSPeter Brune     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_TSDM);CHKERRQ(ierr);
8224989b8cSPeter Brune     ierr = PetscObjectCompose((PetscObject)dm,"TSDM",(PetscObject)container);CHKERRQ(ierr);
83*31ad16ddSPeter Brune     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSDM,DMRestrictHook_TSDM,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__
9124989b8cSPeter Brune #define __FUNCT__ "DMTSGetContextWrite"
9224989b8cSPeter Brune /*@C
9324989b8cSPeter Brune    DMTSGetContextWrite - get write access to private TSDM 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:
10124989b8cSPeter Brune .  tsdm - private TSDM context
10224989b8cSPeter Brune 
10324989b8cSPeter Brune    Level: developer
10424989b8cSPeter Brune 
10524989b8cSPeter Brune .seealso: DMTSGetContext()
10624989b8cSPeter Brune @*/
10724989b8cSPeter Brune PetscErrorCode DMTSGetContextWrite(DM dm,TSDM *tsdm)
10824989b8cSPeter Brune {
10924989b8cSPeter Brune   PetscErrorCode ierr;
11024989b8cSPeter Brune   TSDM           sdm;
11124989b8cSPeter Brune 
11224989b8cSPeter Brune   PetscFunctionBegin;
11324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
11424989b8cSPeter Brune   ierr = DMTSGetContext(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;
11824989b8cSPeter Brune     TSDM         oldsdm = sdm;
11924989b8cSPeter Brune     ierr = PetscInfo(dm,"Copying TSDM due to write\n");CHKERRQ(ierr);
12024989b8cSPeter Brune     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
12124989b8cSPeter Brune     ierr = PetscNewLog(dm,struct _n_TSDM,&sdm);CHKERRQ(ierr);
12224989b8cSPeter Brune     ierr = PetscMemcpy(sdm,oldsdm,sizeof *sdm);CHKERRQ(ierr);
12324989b8cSPeter Brune     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
12424989b8cSPeter Brune     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_TSDM);CHKERRQ(ierr);
12524989b8cSPeter Brune     ierr = PetscObjectCompose((PetscObject)dm,"TSDM",(PetscObject)container);CHKERRQ(ierr);
12624989b8cSPeter Brune     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
12724989b8cSPeter Brune   }
12824989b8cSPeter Brune   *tsdm = sdm;
12924989b8cSPeter Brune   PetscFunctionReturn(0);
13024989b8cSPeter Brune }
13124989b8cSPeter Brune 
13224989b8cSPeter Brune #undef __FUNCT__
13324989b8cSPeter Brune #define __FUNCT__ "DMTSCopyContext"
13424989b8cSPeter Brune /*@C
13524989b8cSPeter Brune    DMTSCopyContext - copies a DM context to a new DM
13624989b8cSPeter Brune 
13724989b8cSPeter Brune    Logically Collective
13824989b8cSPeter Brune 
13924989b8cSPeter Brune    Input Arguments:
14024989b8cSPeter Brune +  dmsrc - DM to obtain context from
14124989b8cSPeter Brune -  dmdest - DM to add context to
14224989b8cSPeter Brune 
14324989b8cSPeter Brune    Level: developer
14424989b8cSPeter Brune 
14524989b8cSPeter Brune    Note:
14624989b8cSPeter Brune    The context is copied by reference. This function does not ensure that a context exists.
14724989b8cSPeter Brune 
14824989b8cSPeter Brune .seealso: DMTSGetContext(), TSSetDM()
14924989b8cSPeter Brune @*/
15024989b8cSPeter Brune PetscErrorCode DMTSCopyContext(DM dmsrc,DM dmdest)
15124989b8cSPeter Brune {
15224989b8cSPeter Brune   PetscErrorCode ierr;
15324989b8cSPeter Brune   PetscContainer container;
15424989b8cSPeter Brune 
15524989b8cSPeter Brune   PetscFunctionBegin;
15624989b8cSPeter Brune   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
15724989b8cSPeter Brune   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
15824989b8cSPeter Brune   ierr = PetscObjectQuery((PetscObject)dmsrc,"TSDM",(PetscObject*)&container);CHKERRQ(ierr);
15924989b8cSPeter Brune   if (container) {
16024989b8cSPeter Brune     ierr = PetscObjectCompose((PetscObject)dmdest,"TSDM",(PetscObject)container);CHKERRQ(ierr);
16124989b8cSPeter Brune     ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_TSDM,DMRestrictHook_TSDM,PETSC_NULL);CHKERRQ(ierr);
16224989b8cSPeter Brune   }
16324989b8cSPeter Brune   PetscFunctionReturn(0);
16424989b8cSPeter Brune }
16524989b8cSPeter Brune 
16624989b8cSPeter Brune #undef __FUNCT__
16724989b8cSPeter Brune #define __FUNCT__ "DMTSSetIFunction"
16824989b8cSPeter Brune /*@C
16924989b8cSPeter Brune    DMTSSetIFunction - set TS implicit function evaluation function
17024989b8cSPeter Brune 
17124989b8cSPeter Brune    Not Collective
17224989b8cSPeter Brune 
17324989b8cSPeter Brune    Input Arguments:
17424989b8cSPeter Brune +  dm - DM to be used with TS
17524989b8cSPeter Brune .  func - function evaluation function, see TSSetIFunction() for calling sequence
17624989b8cSPeter Brune -  ctx - context for residual evaluation
17724989b8cSPeter Brune 
17824989b8cSPeter Brune    Level: advanced
17924989b8cSPeter Brune 
18024989b8cSPeter Brune    Note:
18124989b8cSPeter Brune    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
18224989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
18324989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
18424989b8cSPeter Brune 
18524989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
18624989b8cSPeter Brune @*/
18724989b8cSPeter Brune PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
18824989b8cSPeter Brune {
18924989b8cSPeter Brune   PetscErrorCode ierr;
19024989b8cSPeter Brune   TSDM tsdm;
19124989b8cSPeter Brune 
19224989b8cSPeter Brune   PetscFunctionBegin;
19324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
19424989b8cSPeter Brune   ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr);
19524989b8cSPeter Brune   if (func) tsdm->ifunction = func;
19624989b8cSPeter Brune   if (ctx)  tsdm->ifunctionctx = ctx;
19724989b8cSPeter Brune   PetscFunctionReturn(0);
19824989b8cSPeter Brune }
19924989b8cSPeter Brune 
20024989b8cSPeter Brune #undef __FUNCT__
20124989b8cSPeter Brune #define __FUNCT__ "DMTSGetIFunction"
20224989b8cSPeter Brune /*@C
20324989b8cSPeter Brune    DMTSGetIFunction - get TS implicit residual evaluation function
20424989b8cSPeter Brune 
20524989b8cSPeter Brune    Not Collective
20624989b8cSPeter Brune 
20724989b8cSPeter Brune    Input Argument:
20824989b8cSPeter Brune .  dm - DM to be used with TS
20924989b8cSPeter Brune 
21024989b8cSPeter Brune    Output Arguments:
21124989b8cSPeter Brune +  func - function evaluation function, see TSSetIFunction() for calling sequence
21224989b8cSPeter Brune -  ctx - context for residual evaluation
21324989b8cSPeter Brune 
21424989b8cSPeter Brune    Level: advanced
21524989b8cSPeter Brune 
21624989b8cSPeter Brune    Note:
21724989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
21824989b8cSPeter Brune    associated with the DM.
21924989b8cSPeter Brune 
22024989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
22124989b8cSPeter Brune @*/
22224989b8cSPeter Brune PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
22324989b8cSPeter Brune {
22424989b8cSPeter Brune   PetscErrorCode ierr;
22524989b8cSPeter Brune   TSDM tsdm;
22624989b8cSPeter Brune 
22724989b8cSPeter Brune   PetscFunctionBegin;
22824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
22924989b8cSPeter Brune   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
23024989b8cSPeter Brune   if (func) *func = tsdm->ifunction;
23124989b8cSPeter Brune   if (ctx)  *ctx = tsdm->ifunctionctx;
23224989b8cSPeter Brune   PetscFunctionReturn(0);
23324989b8cSPeter Brune }
23424989b8cSPeter Brune 
23524989b8cSPeter Brune 
23624989b8cSPeter Brune #undef __FUNCT__
23724989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSFunction"
23824989b8cSPeter Brune /*@C
23924989b8cSPeter Brune    DMTSSetRHSFunction - set TS explicit residual evaluation function
24024989b8cSPeter Brune 
24124989b8cSPeter Brune    Not Collective
24224989b8cSPeter Brune 
24324989b8cSPeter Brune    Input Arguments:
24424989b8cSPeter Brune +  dm - DM to be used with TS
24524989b8cSPeter Brune .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
24624989b8cSPeter Brune -  ctx - context for residual evaluation
24724989b8cSPeter Brune 
24824989b8cSPeter Brune    Level: advanced
24924989b8cSPeter Brune 
25024989b8cSPeter Brune    Note:
25124989b8cSPeter Brune    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
25224989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
25324989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
25424989b8cSPeter Brune 
25524989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
25624989b8cSPeter Brune @*/
25724989b8cSPeter Brune PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
25824989b8cSPeter Brune {
25924989b8cSPeter Brune   PetscErrorCode ierr;
26024989b8cSPeter Brune   TSDM tsdm;
26124989b8cSPeter Brune 
26224989b8cSPeter Brune   PetscFunctionBegin;
26324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
26424989b8cSPeter Brune   ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr);
26524989b8cSPeter Brune   if (func) tsdm->rhsfunction = func;
26624989b8cSPeter Brune   if (ctx)  tsdm->rhsfunctionctx = ctx;
26724989b8cSPeter Brune   PetscFunctionReturn(0);
26824989b8cSPeter Brune }
26924989b8cSPeter Brune 
27024989b8cSPeter Brune #undef __FUNCT__
27124989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSFunction"
27224989b8cSPeter Brune /*@C
27324989b8cSPeter Brune    DMTSGetRHSFunction - get TS explicit residual evaluation function
27424989b8cSPeter Brune 
27524989b8cSPeter Brune    Not Collective
27624989b8cSPeter Brune 
27724989b8cSPeter Brune    Input Argument:
27824989b8cSPeter Brune .  dm - DM to be used with TS
27924989b8cSPeter Brune 
28024989b8cSPeter Brune    Output Arguments:
28124989b8cSPeter Brune +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
28224989b8cSPeter Brune -  ctx - context for residual evaluation
28324989b8cSPeter Brune 
28424989b8cSPeter Brune    Level: advanced
28524989b8cSPeter Brune 
28624989b8cSPeter Brune    Note:
28724989b8cSPeter Brune    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
28824989b8cSPeter Brune    associated with the DM.
28924989b8cSPeter Brune 
29024989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
29124989b8cSPeter Brune @*/
29224989b8cSPeter Brune PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
29324989b8cSPeter Brune {
29424989b8cSPeter Brune   PetscErrorCode ierr;
29524989b8cSPeter Brune   TSDM tsdm;
29624989b8cSPeter Brune 
29724989b8cSPeter Brune   PetscFunctionBegin;
29824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
29924989b8cSPeter Brune   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
30024989b8cSPeter Brune   if (func) *func = tsdm->rhsfunction;
30124989b8cSPeter Brune   if (ctx)  *ctx = tsdm->rhsfunctionctx;
30224989b8cSPeter Brune   PetscFunctionReturn(0);
30324989b8cSPeter Brune }
30424989b8cSPeter Brune 
30524989b8cSPeter Brune #undef __FUNCT__
30624989b8cSPeter Brune #define __FUNCT__ "DMTSSetIJacobian"
30724989b8cSPeter Brune /*@C
30824989b8cSPeter Brune    DMTSSetIJacobian - set TS Jacobian evaluation function
30924989b8cSPeter Brune 
31024989b8cSPeter Brune    Not Collective
31124989b8cSPeter Brune 
31224989b8cSPeter Brune    Input Argument:
31324989b8cSPeter Brune +  dm - DM to be used with TS
31424989b8cSPeter Brune .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
31524989b8cSPeter Brune -  ctx - context for residual evaluation
31624989b8cSPeter Brune 
31724989b8cSPeter Brune    Level: advanced
31824989b8cSPeter Brune 
31924989b8cSPeter Brune    Note:
32024989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
32124989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
32224989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
32324989b8cSPeter Brune 
32424989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
32524989b8cSPeter Brune @*/
32624989b8cSPeter Brune PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
32724989b8cSPeter Brune {
32824989b8cSPeter Brune   PetscErrorCode ierr;
32924989b8cSPeter Brune   TSDM sdm;
33024989b8cSPeter Brune 
33124989b8cSPeter Brune   PetscFunctionBegin;
33224989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
33324989b8cSPeter Brune   ierr = DMTSGetContextWrite(dm,&sdm);CHKERRQ(ierr);
33424989b8cSPeter Brune   if (func) sdm->ijacobian = func;
33524989b8cSPeter Brune   if (ctx)  sdm->ijacobianctx = ctx;
33624989b8cSPeter Brune   PetscFunctionReturn(0);
33724989b8cSPeter Brune }
33824989b8cSPeter Brune 
33924989b8cSPeter Brune #undef __FUNCT__
34024989b8cSPeter Brune #define __FUNCT__ "DMTSGetIJacobian"
34124989b8cSPeter Brune /*@C
34224989b8cSPeter Brune    DMTSGetIJacobian - get TS Jacobian evaluation function
34324989b8cSPeter Brune 
34424989b8cSPeter Brune    Not Collective
34524989b8cSPeter Brune 
34624989b8cSPeter Brune    Input Argument:
34724989b8cSPeter Brune .  dm - DM to be used with TS
34824989b8cSPeter Brune 
34924989b8cSPeter Brune    Output Arguments:
35024989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
35124989b8cSPeter Brune -  ctx - context for residual evaluation
35224989b8cSPeter Brune 
35324989b8cSPeter Brune    Level: advanced
35424989b8cSPeter Brune 
35524989b8cSPeter Brune    Note:
35624989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
35724989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
35824989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
35924989b8cSPeter Brune 
36024989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
36124989b8cSPeter Brune @*/
36224989b8cSPeter Brune PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
36324989b8cSPeter Brune {
36424989b8cSPeter Brune   PetscErrorCode ierr;
36524989b8cSPeter Brune   TSDM tsdm;
36624989b8cSPeter Brune 
36724989b8cSPeter Brune   PetscFunctionBegin;
36824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36924989b8cSPeter Brune   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
37024989b8cSPeter Brune   if (func) *func = tsdm->ijacobian;
37124989b8cSPeter Brune   if (ctx)  *ctx = tsdm->ijacobianctx;
37224989b8cSPeter Brune   PetscFunctionReturn(0);
37324989b8cSPeter Brune }
37424989b8cSPeter Brune 
37524989b8cSPeter Brune 
37624989b8cSPeter Brune #undef __FUNCT__
37724989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSJacobian"
37824989b8cSPeter Brune /*@C
37924989b8cSPeter Brune    DMTSSetRHSJacobian - set TS Jacobian evaluation function
38024989b8cSPeter Brune 
38124989b8cSPeter Brune    Not Collective
38224989b8cSPeter Brune 
38324989b8cSPeter Brune    Input Argument:
38424989b8cSPeter Brune +  dm - DM to be used with TS
38524989b8cSPeter Brune .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
38624989b8cSPeter Brune -  ctx - context for residual evaluation
38724989b8cSPeter Brune 
38824989b8cSPeter Brune    Level: advanced
38924989b8cSPeter Brune 
39024989b8cSPeter Brune    Note:
39124989b8cSPeter Brune    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
39224989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
39324989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
39424989b8cSPeter Brune 
39524989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
39624989b8cSPeter Brune @*/
39724989b8cSPeter Brune PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
39824989b8cSPeter Brune {
39924989b8cSPeter Brune   PetscErrorCode ierr;
40024989b8cSPeter Brune   TSDM tsdm;
40124989b8cSPeter Brune 
40224989b8cSPeter Brune   PetscFunctionBegin;
40324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
40424989b8cSPeter Brune   ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr);
40524989b8cSPeter Brune   if (func) tsdm->rhsjacobian = func;
40624989b8cSPeter Brune   if (ctx)  tsdm->rhsjacobianctx = ctx;
40724989b8cSPeter Brune   PetscFunctionReturn(0);
40824989b8cSPeter Brune }
40924989b8cSPeter Brune 
41024989b8cSPeter Brune #undef __FUNCT__
41124989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSJacobian"
41224989b8cSPeter Brune /*@C
41324989b8cSPeter Brune    DMTSGetRHSJacobian - get TS Jacobian evaluation function
41424989b8cSPeter Brune 
41524989b8cSPeter Brune    Not Collective
41624989b8cSPeter Brune 
41724989b8cSPeter Brune    Input Argument:
41824989b8cSPeter Brune .  dm - DM to be used with TS
41924989b8cSPeter Brune 
42024989b8cSPeter Brune    Output Arguments:
42124989b8cSPeter Brune +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
42224989b8cSPeter Brune -  ctx - context for residual evaluation
42324989b8cSPeter Brune 
42424989b8cSPeter Brune    Level: advanced
42524989b8cSPeter Brune 
42624989b8cSPeter Brune    Note:
42724989b8cSPeter Brune    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
42824989b8cSPeter Brune    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
42924989b8cSPeter Brune    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
43024989b8cSPeter Brune 
43124989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
43224989b8cSPeter Brune @*/
43324989b8cSPeter Brune PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
43424989b8cSPeter Brune {
43524989b8cSPeter Brune   PetscErrorCode ierr;
43624989b8cSPeter Brune   TSDM tsdm;
43724989b8cSPeter Brune 
43824989b8cSPeter Brune   PetscFunctionBegin;
43924989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
44024989b8cSPeter Brune   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
44124989b8cSPeter Brune   if (func) *func = tsdm->rhsjacobian;
44224989b8cSPeter Brune   if (ctx)  *ctx = tsdm->rhsjacobianctx;
44324989b8cSPeter Brune   PetscFunctionReturn(0);
44424989b8cSPeter Brune }
445