xref: /petsc/src/ts/utils/dmdats.c (revision 6c6b9e74c1fbf133f434703ab2311cf7e8baf98d)
124989b8cSPeter Brune #include <petscdmda.h>          /*I "petscdmda.h" I*/
224989b8cSPeter Brune #include <petsc-private/dmimpl.h>
324989b8cSPeter Brune #include <petsc-private/tsimpl.h>   /*I "petscts.h" I*/
424989b8cSPeter Brune 
524989b8cSPeter Brune /* This structure holds the user-provided DMDA callbacks */
624989b8cSPeter Brune typedef struct {
724989b8cSPeter Brune   PetscErrorCode (*ifunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
824989b8cSPeter Brune   PetscErrorCode (*rhsfunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
924989b8cSPeter Brune   PetscErrorCode (*ijacobianlocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*);
1024989b8cSPeter Brune   PetscErrorCode (*rhsjacobianlocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*);
1124989b8cSPeter Brune   void *ifunctionlocalctx;
1224989b8cSPeter Brune   void *ijacobianlocalctx;
1324989b8cSPeter Brune   void *rhsfunctionlocalctx;
1424989b8cSPeter Brune   void *rhsjacobianlocalctx;
1524989b8cSPeter Brune   InsertMode ifunctionlocalimode;
1624989b8cSPeter Brune   InsertMode rhsfunctionlocalimode;
1724989b8cSPeter Brune } DM_DA_TS;
1824989b8cSPeter Brune 
1924989b8cSPeter Brune #undef __FUNCT__
2024989b8cSPeter Brune #define __FUNCT__ "TSDMDestroy_DMDA"
2124989b8cSPeter Brune static PetscErrorCode TSDMDestroy_DMDA(TSDM sdm)
2224989b8cSPeter Brune {
2324989b8cSPeter Brune   PetscErrorCode ierr;
2424989b8cSPeter Brune 
2524989b8cSPeter Brune   PetscFunctionBegin;
2624989b8cSPeter Brune   ierr = PetscFree(sdm->data);CHKERRQ(ierr);
2724989b8cSPeter Brune   PetscFunctionReturn(0);
2824989b8cSPeter Brune }
2924989b8cSPeter Brune 
3024989b8cSPeter Brune #undef __FUNCT__
3124989b8cSPeter Brune #define __FUNCT__ "DMDATSGetContext"
3224989b8cSPeter Brune static PetscErrorCode DMDATSGetContext(DM dm,TSDM sdm,DM_DA_TS **dmdats)
3324989b8cSPeter Brune {
3424989b8cSPeter Brune   PetscErrorCode ierr;
3524989b8cSPeter Brune 
3624989b8cSPeter Brune   PetscFunctionBegin;
3724989b8cSPeter Brune   *dmdats = PETSC_NULL;
3824989b8cSPeter Brune   if (!sdm->data) {
3924989b8cSPeter Brune     ierr = PetscNewLog(dm,DM_DA_TS,&sdm->data);CHKERRQ(ierr);
4024989b8cSPeter Brune     sdm->destroy = TSDMDestroy_DMDA;
4124989b8cSPeter Brune   }
4224989b8cSPeter Brune   *dmdats = (DM_DA_TS*)sdm->data;
4324989b8cSPeter Brune   PetscFunctionReturn(0);
4424989b8cSPeter Brune }
4524989b8cSPeter Brune 
4624989b8cSPeter Brune #undef __FUNCT__
4724989b8cSPeter Brune #define __FUNCT__ "TSComputeIFunction_DMDA"
4824989b8cSPeter Brune /*
4924989b8cSPeter Brune   This function should eventually replace:
5024989b8cSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
5124989b8cSPeter Brune  */
5224989b8cSPeter Brune static PetscErrorCode TSComputeIFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *ctx)
5324989b8cSPeter Brune {
5424989b8cSPeter Brune   PetscErrorCode ierr;
5524989b8cSPeter Brune   DM             dm;
5624989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
5724989b8cSPeter Brune   DMDALocalInfo  info;
5824989b8cSPeter Brune   Vec            Xloc;
5924989b8cSPeter Brune   void           *x,*f,*xdot;
6024989b8cSPeter Brune 
6124989b8cSPeter Brune   PetscFunctionBegin;
6224989b8cSPeter Brune   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6324989b8cSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
6424989b8cSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
6524989b8cSPeter Brune   if (!dmdats->ifunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
6624989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
6724989b8cSPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
6824989b8cSPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
6924989b8cSPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
7024989b8cSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
7124989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
7224989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr);
7324989b8cSPeter Brune   switch (dmdats->ifunctionlocalimode) {
7424989b8cSPeter Brune   case INSERT_VALUES: {
7524989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
7624989b8cSPeter Brune     CHKMEMQ;
7724989b8cSPeter Brune     ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr);
7824989b8cSPeter Brune     CHKMEMQ;
7924989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
8024989b8cSPeter Brune   } break;
8124989b8cSPeter Brune   case ADD_VALUES: {
8224989b8cSPeter Brune     Vec Floc;
8324989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
8424989b8cSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
8524989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
8624989b8cSPeter Brune     CHKMEMQ;
8724989b8cSPeter Brune     ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr);
8824989b8cSPeter Brune     CHKMEMQ;
8924989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
9024989b8cSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
9124989b8cSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
9224989b8cSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
9324989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
9424989b8cSPeter Brune   } break;
9524989b8cSPeter Brune   default: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->ifunctionlocalimode);
9624989b8cSPeter Brune   }
9724989b8cSPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
9824989b8cSPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
9924989b8cSPeter Brune   PetscFunctionReturn(0);
10024989b8cSPeter Brune }
10124989b8cSPeter Brune 
10224989b8cSPeter Brune #undef __FUNCT__
10324989b8cSPeter Brune #define __FUNCT__ "TSComputeIJacobian_DMDA"
10424989b8cSPeter Brune /*
10524989b8cSPeter Brune   This function should eventually replace:
10624989b8cSPeter Brune     DMComputeJacobian() and DMDAComputeJacobian1()
10724989b8cSPeter Brune  */
10824989b8cSPeter Brune static PetscErrorCode TSComputeIJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
10924989b8cSPeter Brune {
11024989b8cSPeter Brune   PetscErrorCode ierr;
11124989b8cSPeter Brune   DM             dm;
11224989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
11324989b8cSPeter Brune   DMDALocalInfo  info;
11424989b8cSPeter Brune   Vec            Xloc;
11524989b8cSPeter Brune   void           *x,*xdot;
11624989b8cSPeter Brune 
11724989b8cSPeter Brune   PetscFunctionBegin;
11824989b8cSPeter Brune   if (!dmdats->ifunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
11924989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
12024989b8cSPeter Brune 
12124989b8cSPeter Brune   if (dmdats->ijacobianlocal) {
12224989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
12324989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
12424989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
12524989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
12624989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
12724989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr);
12824989b8cSPeter Brune     CHKMEMQ;
12924989b8cSPeter Brune     ierr = (*dmdats->ijacobianlocal)(&info,ptime,x,xdot,shift,*A,*B,mstr,dmdats->ijacobianlocalctx);CHKERRQ(ierr);
13024989b8cSPeter Brune     CHKMEMQ;
13124989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
13224989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr);
13324989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
13424989b8cSPeter Brune   } else {
13524989b8cSPeter Brune     SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()");
13624989b8cSPeter Brune   }
13724989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
13824989b8cSPeter Brune   if (*A != *B) {
13924989b8cSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14024989b8cSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14124989b8cSPeter Brune   }
14224989b8cSPeter Brune   PetscFunctionReturn(0);
14324989b8cSPeter Brune }
14424989b8cSPeter Brune 
14524989b8cSPeter Brune #undef __FUNCT__
14624989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSFunction_DMDA"
14724989b8cSPeter Brune /*
14824989b8cSPeter Brune   This function should eventually replace:
14924989b8cSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
15024989b8cSPeter Brune  */
15124989b8cSPeter Brune static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec F,void *ctx)
15224989b8cSPeter Brune {
15324989b8cSPeter Brune   PetscErrorCode ierr;
15424989b8cSPeter Brune   DM             dm;
15524989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
15624989b8cSPeter Brune   DMDALocalInfo  info;
15724989b8cSPeter Brune   Vec            Xloc;
15824989b8cSPeter Brune   void           *x,*f;
15924989b8cSPeter Brune 
16024989b8cSPeter Brune   PetscFunctionBegin;
16124989b8cSPeter Brune   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
16224989b8cSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
16324989b8cSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
16424989b8cSPeter Brune   if (!dmdats->rhsfunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
16524989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
16624989b8cSPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
16724989b8cSPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
16824989b8cSPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
16924989b8cSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
17024989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
17124989b8cSPeter Brune   switch (dmdats->rhsfunctionlocalimode) {
17224989b8cSPeter Brune   case INSERT_VALUES: {
17324989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
17424989b8cSPeter Brune     CHKMEMQ;
17524989b8cSPeter Brune     ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr);
17624989b8cSPeter Brune     CHKMEMQ;
17724989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
17824989b8cSPeter Brune   } break;
17924989b8cSPeter Brune   case ADD_VALUES: {
18024989b8cSPeter Brune     Vec Floc;
18124989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
18224989b8cSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
18324989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
18424989b8cSPeter Brune     CHKMEMQ;
18524989b8cSPeter Brune     ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr);
18624989b8cSPeter Brune     CHKMEMQ;
18724989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
18824989b8cSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
18924989b8cSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
19024989b8cSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
19124989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
19224989b8cSPeter Brune   } break;
19324989b8cSPeter Brune   default: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->rhsfunctionlocalimode);
19424989b8cSPeter Brune   }
19524989b8cSPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
19624989b8cSPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
19724989b8cSPeter Brune   PetscFunctionReturn(0);
19824989b8cSPeter Brune }
19924989b8cSPeter Brune 
20024989b8cSPeter Brune #undef __FUNCT__
20124989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSJacobian_DMDA"
20224989b8cSPeter Brune /*
20324989b8cSPeter Brune   This function should eventually replace:
20424989b8cSPeter Brune     DMComputeJacobian() and DMDAComputeJacobian1()
20524989b8cSPeter Brune  */
20624989b8cSPeter Brune static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
20724989b8cSPeter Brune {
20824989b8cSPeter Brune   PetscErrorCode ierr;
20924989b8cSPeter Brune   DM             dm;
21024989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
21124989b8cSPeter Brune   DMDALocalInfo  info;
21224989b8cSPeter Brune   Vec            Xloc;
21324989b8cSPeter Brune   void           *x;
21424989b8cSPeter Brune 
21524989b8cSPeter Brune   PetscFunctionBegin;
21624989b8cSPeter Brune   if (!dmdats->rhsfunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
21724989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
21824989b8cSPeter Brune 
21924989b8cSPeter Brune   if (dmdats->rhsjacobianlocal) {
22024989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
22124989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
22224989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
22324989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
22424989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
22524989b8cSPeter Brune     CHKMEMQ;
22624989b8cSPeter Brune     ierr = (*dmdats->rhsjacobianlocal)(&info,ptime,x,*A,*B,mstr,dmdats->rhsjacobianlocalctx);CHKERRQ(ierr);
22724989b8cSPeter Brune     CHKMEMQ;
22824989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
22924989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
23024989b8cSPeter Brune   } else {
23124989b8cSPeter Brune     SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()");
23224989b8cSPeter Brune   }
23324989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
23424989b8cSPeter Brune   if (*A != *B) {
23524989b8cSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
23624989b8cSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
23724989b8cSPeter Brune   }
23824989b8cSPeter Brune   PetscFunctionReturn(0);
23924989b8cSPeter Brune }
24024989b8cSPeter Brune 
24124989b8cSPeter Brune 
24224989b8cSPeter Brune #undef __FUNCT__
24324989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSFunctionLocal"
24424989b8cSPeter Brune /*@C
24524989b8cSPeter Brune    DMDATSSetRHSFunctionLocal - set a local residual evaluation function
24624989b8cSPeter Brune 
24724989b8cSPeter Brune    Logically Collective
24824989b8cSPeter Brune 
24924989b8cSPeter Brune    Input Arguments:
25024989b8cSPeter Brune +  dm - DM to associate callback with
25124989b8cSPeter Brune .  func - local residual evaluation
25224989b8cSPeter Brune -  ctx - optional context for local residual evaluation
25324989b8cSPeter Brune 
25424989b8cSPeter Brune    Calling sequence for func:
25524989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
25624989b8cSPeter Brune .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
25724989b8cSPeter Brune .  x - dimensional pointer to state at which to evaluate residual
25824989b8cSPeter Brune .  f - dimensional pointer to residual, write the residual here
25924989b8cSPeter Brune -  ctx - optional context passed above
26024989b8cSPeter Brune 
26124989b8cSPeter Brune    Level: beginner
26224989b8cSPeter Brune 
26324989b8cSPeter Brune .seealso: DMTSSetFunction(), DMDATSSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
26424989b8cSPeter Brune @*/
265*6c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm,InsertMode imode,DMDATSRHSFunctionLocal func,void *ctx)
26624989b8cSPeter Brune {
26724989b8cSPeter Brune   PetscErrorCode ierr;
26824989b8cSPeter Brune   TSDM         sdm;
26924989b8cSPeter Brune   DM_DA_TS     *dmdats;
27024989b8cSPeter Brune 
27124989b8cSPeter Brune   PetscFunctionBegin;
27224989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
27324989b8cSPeter Brune   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
27424989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
27524989b8cSPeter Brune   dmdats->rhsfunctionlocalimode = imode;
27624989b8cSPeter Brune   dmdats->rhsfunctionlocal = func;
27724989b8cSPeter Brune   dmdats->rhsfunctionlocalctx = ctx;
27824989b8cSPeter Brune   ierr = DMTSSetRHSFunction(dm,TSComputeRHSFunction_DMDA,dmdats);CHKERRQ(ierr);
27924989b8cSPeter Brune   PetscFunctionReturn(0);
28024989b8cSPeter Brune }
28124989b8cSPeter Brune 
28224989b8cSPeter Brune #undef __FUNCT__
28324989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSJacobianLocal"
28424989b8cSPeter Brune /*@C
28524989b8cSPeter Brune    DMDATSSetRHSJacobianLocal - set a local residual evaluation function
28624989b8cSPeter Brune 
28724989b8cSPeter Brune    Logically Collective
28824989b8cSPeter Brune 
28924989b8cSPeter Brune    Input Arguments:
29024989b8cSPeter Brune +  dm - DM to associate callback with
29124989b8cSPeter Brune .  func - local residual evaluation
29224989b8cSPeter Brune -  ctx - optional context for local residual evaluation
29324989b8cSPeter Brune 
29424989b8cSPeter Brune    Calling sequence for func:
29524989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
29624989b8cSPeter Brune .  x - dimensional pointer to state at which to evaluate residual
29724989b8cSPeter Brune .  f - dimensional pointer to residual, write the residual here
29824989b8cSPeter Brune -  ctx - optional context passed above
29924989b8cSPeter Brune 
30024989b8cSPeter Brune    Level: beginner
30124989b8cSPeter Brune 
30224989b8cSPeter Brune .seealso: DMTSSetJacobian(), DMDATSSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
30324989b8cSPeter Brune @*/
304*6c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm,DMDATSRHSJacobianLocal func,void *ctx)
30524989b8cSPeter Brune {
30624989b8cSPeter Brune   PetscErrorCode ierr;
30724989b8cSPeter Brune   TSDM         sdm;
30824989b8cSPeter Brune   DM_DA_TS     *dmdats;
30924989b8cSPeter Brune 
31024989b8cSPeter Brune   PetscFunctionBegin;
31124989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
31224989b8cSPeter Brune   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
31324989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
31424989b8cSPeter Brune   dmdats->rhsjacobianlocal = func;
31524989b8cSPeter Brune   dmdats->rhsjacobianlocalctx = ctx;
31624989b8cSPeter Brune   ierr = DMTSSetRHSJacobian(dm,TSComputeRHSJacobian_DMDA,dmdats);CHKERRQ(ierr);
31724989b8cSPeter Brune   PetscFunctionReturn(0);
31824989b8cSPeter Brune }
31924989b8cSPeter Brune 
32024989b8cSPeter Brune 
32124989b8cSPeter Brune #undef __FUNCT__
32224989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIFunctionLocal"
32324989b8cSPeter Brune /*@C
32424989b8cSPeter Brune    DMDATSSetIFunctionLocal - set a local residual evaluation function
32524989b8cSPeter Brune 
32624989b8cSPeter Brune    Logically Collective
32724989b8cSPeter Brune 
32824989b8cSPeter Brune    Input Arguments:
32924989b8cSPeter Brune +  dm - DM to associate callback with
33024989b8cSPeter Brune .  func - local residual evaluation
33124989b8cSPeter Brune -  ctx - optional context for local residual evaluation
33224989b8cSPeter Brune 
33324989b8cSPeter Brune    Calling sequence for func:
33424989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
33524989b8cSPeter Brune .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
33624989b8cSPeter Brune .  x - dimensional pointer to state at which to evaluate residual
33724989b8cSPeter Brune .  f - dimensional pointer to residual, write the residual here
33824989b8cSPeter Brune -  ctx - optional context passed above
33924989b8cSPeter Brune 
34024989b8cSPeter Brune    Level: beginner
34124989b8cSPeter Brune 
34224989b8cSPeter Brune .seealso: DMTSSetFunction(), DMDATSSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
34324989b8cSPeter Brune @*/
344*6c6b9e74SPeter Brune PetscErrorCode DMDATSSetIFunctionLocal(DM dm,InsertMode imode,DMDATSIFunctionLocal func,void *ctx)
34524989b8cSPeter Brune {
34624989b8cSPeter Brune   PetscErrorCode ierr;
34724989b8cSPeter Brune   TSDM         sdm;
34824989b8cSPeter Brune   DM_DA_TS     *dmdats;
34924989b8cSPeter Brune 
35024989b8cSPeter Brune   PetscFunctionBegin;
35124989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
35224989b8cSPeter Brune   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
35324989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
35424989b8cSPeter Brune   dmdats->ifunctionlocalimode = imode;
35524989b8cSPeter Brune   dmdats->ifunctionlocal = func;
35624989b8cSPeter Brune   dmdats->ifunctionlocalctx = ctx;
35724989b8cSPeter Brune   ierr = DMTSSetIFunction(dm,TSComputeIFunction_DMDA,dmdats);CHKERRQ(ierr);
35824989b8cSPeter Brune   PetscFunctionReturn(0);
35924989b8cSPeter Brune }
36024989b8cSPeter Brune 
36124989b8cSPeter Brune #undef __FUNCT__
36224989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIJacobianLocal"
36324989b8cSPeter Brune /*@C
36424989b8cSPeter Brune    DMDATSSetIJacobianLocal - set a local residual evaluation function
36524989b8cSPeter Brune 
36624989b8cSPeter Brune    Logically Collective
36724989b8cSPeter Brune 
36824989b8cSPeter Brune    Input Arguments:
36924989b8cSPeter Brune +  dm - DM to associate callback with
37024989b8cSPeter Brune .  func - local residual evaluation
37124989b8cSPeter Brune -  ctx - optional context for local residual evaluation
37224989b8cSPeter Brune 
37324989b8cSPeter Brune    Calling sequence for func:
37424989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
37524989b8cSPeter Brune .  x - dimensional pointer to state at which to evaluate residual
37624989b8cSPeter Brune .  f - dimensional pointer to residual, write the residual here
37724989b8cSPeter Brune -  ctx - optional context passed above
37824989b8cSPeter Brune 
37924989b8cSPeter Brune    Level: beginner
38024989b8cSPeter Brune 
38124989b8cSPeter Brune .seealso: DMTSSetJacobian(), DMDATSSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
38224989b8cSPeter Brune @*/
383*6c6b9e74SPeter Brune PetscErrorCode DMDATSSetIJacobianLocal(DM dm,DMDATSIJacobianLocal func,void *ctx)
38424989b8cSPeter Brune {
38524989b8cSPeter Brune   PetscErrorCode ierr;
38624989b8cSPeter Brune   TSDM         sdm;
38724989b8cSPeter Brune   DM_DA_TS     *dmdats;
38824989b8cSPeter Brune 
38924989b8cSPeter Brune   PetscFunctionBegin;
39024989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
39124989b8cSPeter Brune   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
39224989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
39324989b8cSPeter Brune   dmdats->ijacobianlocal = func;
39424989b8cSPeter Brune   dmdats->ijacobianlocalctx = ctx;
39524989b8cSPeter Brune   ierr = DMTSSetIJacobian(dm,TSComputeIJacobian_DMDA,dmdats);CHKERRQ(ierr);
39624989b8cSPeter Brune   PetscFunctionReturn(0);
39724989b8cSPeter Brune }
398