xref: /petsc/src/ts/utils/dmdats.c (revision 0313d0e01224e48a64c42069d7ba393b1648c687)
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__
31*0313d0e0SJed Brown #define __FUNCT__ "TSDMDuplicate_DMDA"
32*0313d0e0SJed Brown static PetscErrorCode TSDMDuplicate_DMDA(TSDM oldsdm,DM dm)
33*0313d0e0SJed Brown {
34*0313d0e0SJed Brown   PetscErrorCode ierr;
35*0313d0e0SJed Brown   TSDM           sdm;
36*0313d0e0SJed Brown 
37*0313d0e0SJed Brown   PetscFunctionBegin;
38*0313d0e0SJed Brown   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
39*0313d0e0SJed Brown   ierr = PetscNewLog(dm,DM_DA_TS,&sdm->data);CHKERRQ(ierr);
40*0313d0e0SJed Brown   if (oldsdm->data) {ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DM_DA_TS));CHKERRQ(ierr);}
41*0313d0e0SJed Brown   PetscFunctionReturn(0);
42*0313d0e0SJed Brown }
43*0313d0e0SJed Brown 
44*0313d0e0SJed Brown #undef __FUNCT__
4524989b8cSPeter Brune #define __FUNCT__ "DMDATSGetContext"
4624989b8cSPeter Brune static PetscErrorCode DMDATSGetContext(DM dm,TSDM sdm,DM_DA_TS **dmdats)
4724989b8cSPeter Brune {
4824989b8cSPeter Brune   PetscErrorCode ierr;
4924989b8cSPeter Brune 
5024989b8cSPeter Brune   PetscFunctionBegin;
5124989b8cSPeter Brune   *dmdats = PETSC_NULL;
5224989b8cSPeter Brune   if (!sdm->data) {
5324989b8cSPeter Brune     ierr = PetscNewLog(dm,DM_DA_TS,&sdm->data);CHKERRQ(ierr);
5424989b8cSPeter Brune     sdm->destroy = TSDMDestroy_DMDA;
55*0313d0e0SJed Brown     sdm->duplicate = TSDMDuplicate_DMDA;
5624989b8cSPeter Brune   }
5724989b8cSPeter Brune   *dmdats = (DM_DA_TS*)sdm->data;
5824989b8cSPeter Brune   PetscFunctionReturn(0);
5924989b8cSPeter Brune }
6024989b8cSPeter Brune 
6124989b8cSPeter Brune #undef __FUNCT__
6224989b8cSPeter Brune #define __FUNCT__ "TSComputeIFunction_DMDA"
6324989b8cSPeter Brune /*
6424989b8cSPeter Brune   This function should eventually replace:
6524989b8cSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
6624989b8cSPeter Brune  */
6724989b8cSPeter Brune static PetscErrorCode TSComputeIFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *ctx)
6824989b8cSPeter Brune {
6924989b8cSPeter Brune   PetscErrorCode ierr;
7024989b8cSPeter Brune   DM             dm;
7124989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
7224989b8cSPeter Brune   DMDALocalInfo  info;
7324989b8cSPeter Brune   Vec            Xloc;
7424989b8cSPeter Brune   void           *x,*f,*xdot;
7524989b8cSPeter Brune 
7624989b8cSPeter Brune   PetscFunctionBegin;
7724989b8cSPeter Brune   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
7824989b8cSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
7924989b8cSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
8024989b8cSPeter Brune   if (!dmdats->ifunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
8124989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
8224989b8cSPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
8324989b8cSPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
8424989b8cSPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
8524989b8cSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
8624989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
8724989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr);
8824989b8cSPeter Brune   switch (dmdats->ifunctionlocalimode) {
8924989b8cSPeter Brune   case INSERT_VALUES: {
9024989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
9124989b8cSPeter Brune     CHKMEMQ;
9224989b8cSPeter Brune     ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr);
9324989b8cSPeter Brune     CHKMEMQ;
9424989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
9524989b8cSPeter Brune   } break;
9624989b8cSPeter Brune   case ADD_VALUES: {
9724989b8cSPeter Brune     Vec Floc;
9824989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
9924989b8cSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
10024989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
10124989b8cSPeter Brune     CHKMEMQ;
10224989b8cSPeter Brune     ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr);
10324989b8cSPeter Brune     CHKMEMQ;
10424989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
10524989b8cSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
10624989b8cSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
10724989b8cSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
10824989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
10924989b8cSPeter Brune   } break;
11024989b8cSPeter Brune   default: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->ifunctionlocalimode);
11124989b8cSPeter Brune   }
11224989b8cSPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
11324989b8cSPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
11424989b8cSPeter Brune   PetscFunctionReturn(0);
11524989b8cSPeter Brune }
11624989b8cSPeter Brune 
11724989b8cSPeter Brune #undef __FUNCT__
11824989b8cSPeter Brune #define __FUNCT__ "TSComputeIJacobian_DMDA"
11924989b8cSPeter Brune /*
12024989b8cSPeter Brune   This function should eventually replace:
12124989b8cSPeter Brune     DMComputeJacobian() and DMDAComputeJacobian1()
12224989b8cSPeter Brune  */
12324989b8cSPeter Brune static PetscErrorCode TSComputeIJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
12424989b8cSPeter Brune {
12524989b8cSPeter Brune   PetscErrorCode ierr;
12624989b8cSPeter Brune   DM             dm;
12724989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
12824989b8cSPeter Brune   DMDALocalInfo  info;
12924989b8cSPeter Brune   Vec            Xloc;
13024989b8cSPeter Brune   void           *x,*xdot;
13124989b8cSPeter Brune 
13224989b8cSPeter Brune   PetscFunctionBegin;
13324989b8cSPeter Brune   if (!dmdats->ifunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
13424989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
13524989b8cSPeter Brune 
13624989b8cSPeter Brune   if (dmdats->ijacobianlocal) {
13724989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
13824989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
13924989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
14024989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
14124989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
14224989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr);
14324989b8cSPeter Brune     CHKMEMQ;
14424989b8cSPeter Brune     ierr = (*dmdats->ijacobianlocal)(&info,ptime,x,xdot,shift,*A,*B,mstr,dmdats->ijacobianlocalctx);CHKERRQ(ierr);
14524989b8cSPeter Brune     CHKMEMQ;
14624989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
14724989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr);
14824989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
14924989b8cSPeter Brune   } else {
15024989b8cSPeter Brune     SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()");
15124989b8cSPeter Brune   }
15224989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
15324989b8cSPeter Brune   if (*A != *B) {
15424989b8cSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15524989b8cSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15624989b8cSPeter Brune   }
15724989b8cSPeter Brune   PetscFunctionReturn(0);
15824989b8cSPeter Brune }
15924989b8cSPeter Brune 
16024989b8cSPeter Brune #undef __FUNCT__
16124989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSFunction_DMDA"
16224989b8cSPeter Brune /*
16324989b8cSPeter Brune   This function should eventually replace:
16424989b8cSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
16524989b8cSPeter Brune  */
16624989b8cSPeter Brune static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec F,void *ctx)
16724989b8cSPeter Brune {
16824989b8cSPeter Brune   PetscErrorCode ierr;
16924989b8cSPeter Brune   DM             dm;
17024989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
17124989b8cSPeter Brune   DMDALocalInfo  info;
17224989b8cSPeter Brune   Vec            Xloc;
17324989b8cSPeter Brune   void           *x,*f;
17424989b8cSPeter Brune 
17524989b8cSPeter Brune   PetscFunctionBegin;
17624989b8cSPeter Brune   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
17724989b8cSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
17824989b8cSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
17924989b8cSPeter Brune   if (!dmdats->rhsfunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
18024989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
18124989b8cSPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
18224989b8cSPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
18324989b8cSPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
18424989b8cSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
18524989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
18624989b8cSPeter Brune   switch (dmdats->rhsfunctionlocalimode) {
18724989b8cSPeter Brune   case INSERT_VALUES: {
18824989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
18924989b8cSPeter Brune     CHKMEMQ;
19024989b8cSPeter Brune     ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr);
19124989b8cSPeter Brune     CHKMEMQ;
19224989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
19324989b8cSPeter Brune   } break;
19424989b8cSPeter Brune   case ADD_VALUES: {
19524989b8cSPeter Brune     Vec Floc;
19624989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
19724989b8cSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
19824989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
19924989b8cSPeter Brune     CHKMEMQ;
20024989b8cSPeter Brune     ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr);
20124989b8cSPeter Brune     CHKMEMQ;
20224989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
20324989b8cSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
20424989b8cSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
20524989b8cSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
20624989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
20724989b8cSPeter Brune   } break;
20824989b8cSPeter Brune   default: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->rhsfunctionlocalimode);
20924989b8cSPeter Brune   }
21024989b8cSPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
21124989b8cSPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
21224989b8cSPeter Brune   PetscFunctionReturn(0);
21324989b8cSPeter Brune }
21424989b8cSPeter Brune 
21524989b8cSPeter Brune #undef __FUNCT__
21624989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSJacobian_DMDA"
21724989b8cSPeter Brune /*
21824989b8cSPeter Brune   This function should eventually replace:
21924989b8cSPeter Brune     DMComputeJacobian() and DMDAComputeJacobian1()
22024989b8cSPeter Brune  */
22124989b8cSPeter Brune static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
22224989b8cSPeter Brune {
22324989b8cSPeter Brune   PetscErrorCode ierr;
22424989b8cSPeter Brune   DM             dm;
22524989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
22624989b8cSPeter Brune   DMDALocalInfo  info;
22724989b8cSPeter Brune   Vec            Xloc;
22824989b8cSPeter Brune   void           *x;
22924989b8cSPeter Brune 
23024989b8cSPeter Brune   PetscFunctionBegin;
23124989b8cSPeter Brune   if (!dmdats->rhsfunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
23224989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
23324989b8cSPeter Brune 
23424989b8cSPeter Brune   if (dmdats->rhsjacobianlocal) {
23524989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
23624989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
23724989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
23824989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
23924989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
24024989b8cSPeter Brune     CHKMEMQ;
24124989b8cSPeter Brune     ierr = (*dmdats->rhsjacobianlocal)(&info,ptime,x,*A,*B,mstr,dmdats->rhsjacobianlocalctx);CHKERRQ(ierr);
24224989b8cSPeter Brune     CHKMEMQ;
24324989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
24424989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
24524989b8cSPeter Brune   } else {
24624989b8cSPeter Brune     SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()");
24724989b8cSPeter Brune   }
24824989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
24924989b8cSPeter Brune   if (*A != *B) {
25024989b8cSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
25124989b8cSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
25224989b8cSPeter Brune   }
25324989b8cSPeter Brune   PetscFunctionReturn(0);
25424989b8cSPeter Brune }
25524989b8cSPeter Brune 
25624989b8cSPeter Brune 
25724989b8cSPeter Brune #undef __FUNCT__
25824989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSFunctionLocal"
25924989b8cSPeter Brune /*@C
26024989b8cSPeter Brune    DMDATSSetRHSFunctionLocal - set a local residual evaluation function
26124989b8cSPeter Brune 
26224989b8cSPeter Brune    Logically Collective
26324989b8cSPeter Brune 
26424989b8cSPeter Brune    Input Arguments:
26524989b8cSPeter Brune +  dm - DM to associate callback with
266b398ffe5SPeter Brune .  imode - insert mode for the residual
26724989b8cSPeter Brune .  func - local residual evaluation
26824989b8cSPeter Brune -  ctx - optional context for local residual evaluation
26924989b8cSPeter Brune 
27024989b8cSPeter Brune    Calling sequence for func:
271b398ffe5SPeter Brune 
272b398ffe5SPeter Brune $ func(DMDALocalInfo info,PetscReal t,void *x,void *f,void *ctx)
273b398ffe5SPeter Brune 
27424989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
275b398ffe5SPeter Brune .  t - time at which to evaluate residual
276b398ffe5SPeter Brune .  x - array of local state information
277b398ffe5SPeter Brune .  f - output array of local residual information
278b398ffe5SPeter Brune -  ctx - optional user context
27924989b8cSPeter Brune 
28024989b8cSPeter Brune    Level: beginner
28124989b8cSPeter Brune 
2828e625dcfSPeter Brune .seealso: DMTSSetRHSFunction(), DMDATSSetRHSJacobianLocal(), DMDASNESSetFunctionLocal()
28324989b8cSPeter Brune @*/
2846c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm,InsertMode imode,DMDATSRHSFunctionLocal func,void *ctx)
28524989b8cSPeter Brune {
28624989b8cSPeter Brune   PetscErrorCode ierr;
28724989b8cSPeter Brune   TSDM         sdm;
28824989b8cSPeter Brune   DM_DA_TS     *dmdats;
28924989b8cSPeter Brune 
29024989b8cSPeter Brune   PetscFunctionBegin;
29124989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
292*0313d0e0SJed Brown   ierr = DMTSGetContextWrite(dm,&sdm);CHKERRQ(ierr);
29324989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
29424989b8cSPeter Brune   dmdats->rhsfunctionlocalimode = imode;
29524989b8cSPeter Brune   dmdats->rhsfunctionlocal = func;
29624989b8cSPeter Brune   dmdats->rhsfunctionlocalctx = ctx;
29724989b8cSPeter Brune   ierr = DMTSSetRHSFunction(dm,TSComputeRHSFunction_DMDA,dmdats);CHKERRQ(ierr);
29824989b8cSPeter Brune   PetscFunctionReturn(0);
29924989b8cSPeter Brune }
30024989b8cSPeter Brune 
30124989b8cSPeter Brune #undef __FUNCT__
30224989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSJacobianLocal"
30324989b8cSPeter Brune /*@C
30424989b8cSPeter Brune    DMDATSSetRHSJacobianLocal - set a local residual evaluation function
30524989b8cSPeter Brune 
30624989b8cSPeter Brune    Logically Collective
30724989b8cSPeter Brune 
30824989b8cSPeter Brune    Input Arguments:
30924989b8cSPeter Brune +  dm    - DM to associate callback with
310b398ffe5SPeter Brune .  func  - local RHS Jacobian evaluation routine
311b398ffe5SPeter Brune -  ctx   - optional context for local jacobian evaluation
31224989b8cSPeter Brune 
31324989b8cSPeter Brune    Calling sequence for func:
314b398ffe5SPeter Brune 
315b398ffe5SPeter Brune $ func(DMDALocalInfo* info,PetscReal t,void* x,Mat J,Mat B,MatStructure *flg,void *ctx);
316b398ffe5SPeter Brune 
31724989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
318b398ffe5SPeter Brune .  t    - time at which to evaluate residual
319b398ffe5SPeter Brune .  x    - array of local state information
320b398ffe5SPeter Brune .  J    - Jacobian matrix
321b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
322b398ffe5SPeter Brune .  flg  - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators())
32324989b8cSPeter Brune -  ctx  - optional context passed above
32424989b8cSPeter Brune 
32524989b8cSPeter Brune    Level: beginner
32624989b8cSPeter Brune 
3278e625dcfSPeter Brune .seealso: DMTSSetRHSJacobian(), DMDATSSetRHSFunctionLocal(), DMDASNESSetJacobianLocal()
32824989b8cSPeter Brune @*/
3296c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm,DMDATSRHSJacobianLocal func,void *ctx)
33024989b8cSPeter Brune {
33124989b8cSPeter Brune   PetscErrorCode ierr;
33224989b8cSPeter Brune   TSDM         sdm;
33324989b8cSPeter Brune   DM_DA_TS     *dmdats;
33424989b8cSPeter Brune 
33524989b8cSPeter Brune   PetscFunctionBegin;
33624989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
337*0313d0e0SJed Brown   ierr = DMTSGetContextWrite(dm,&sdm);CHKERRQ(ierr);
33824989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
33924989b8cSPeter Brune   dmdats->rhsjacobianlocal = func;
34024989b8cSPeter Brune   dmdats->rhsjacobianlocalctx = ctx;
34124989b8cSPeter Brune   ierr = DMTSSetRHSJacobian(dm,TSComputeRHSJacobian_DMDA,dmdats);CHKERRQ(ierr);
34224989b8cSPeter Brune   PetscFunctionReturn(0);
34324989b8cSPeter Brune }
34424989b8cSPeter Brune 
34524989b8cSPeter Brune 
34624989b8cSPeter Brune #undef __FUNCT__
34724989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIFunctionLocal"
34824989b8cSPeter Brune /*@C
34924989b8cSPeter Brune    DMDATSSetIFunctionLocal - set a local residual evaluation function
35024989b8cSPeter Brune 
35124989b8cSPeter Brune    Logically Collective
35224989b8cSPeter Brune 
35324989b8cSPeter Brune    Input Arguments:
35424989b8cSPeter Brune +  dm   - DM to associate callback with
35524989b8cSPeter Brune .  func - local residual evaluation
35624989b8cSPeter Brune -  ctx  - optional context for local residual evaluation
35724989b8cSPeter Brune 
35824989b8cSPeter Brune    Calling sequence for func:
35924989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
360b398ffe5SPeter Brune .  t    - time at which to evaluate residual
361b398ffe5SPeter Brune .  x    - array of local state information
362b398ffe5SPeter Brune .  xdot - array of local time derivative information
363b398ffe5SPeter Brune .  f    - output array of local function evaluation information
36424989b8cSPeter Brune -  ctx - optional context passed above
36524989b8cSPeter Brune 
36624989b8cSPeter Brune    Level: beginner
36724989b8cSPeter Brune 
3688e625dcfSPeter Brune .seealso: DMTSSetIFunction(), DMDATSSetIJacobianLocal(), DMDASNESSetFunctionLocal()
36924989b8cSPeter Brune @*/
3706c6b9e74SPeter Brune PetscErrorCode DMDATSSetIFunctionLocal(DM dm,InsertMode imode,DMDATSIFunctionLocal func,void *ctx)
37124989b8cSPeter Brune {
37224989b8cSPeter Brune   PetscErrorCode ierr;
37324989b8cSPeter Brune   TSDM         sdm;
37424989b8cSPeter Brune   DM_DA_TS     *dmdats;
37524989b8cSPeter Brune 
37624989b8cSPeter Brune   PetscFunctionBegin;
37724989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
378*0313d0e0SJed Brown   ierr = DMTSGetContextWrite(dm,&sdm);CHKERRQ(ierr);
37924989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
38024989b8cSPeter Brune   dmdats->ifunctionlocalimode = imode;
38124989b8cSPeter Brune   dmdats->ifunctionlocal = func;
38224989b8cSPeter Brune   dmdats->ifunctionlocalctx = ctx;
38324989b8cSPeter Brune   ierr = DMTSSetIFunction(dm,TSComputeIFunction_DMDA,dmdats);CHKERRQ(ierr);
38424989b8cSPeter Brune   PetscFunctionReturn(0);
38524989b8cSPeter Brune }
38624989b8cSPeter Brune 
38724989b8cSPeter Brune #undef __FUNCT__
38824989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIJacobianLocal"
38924989b8cSPeter Brune /*@C
39024989b8cSPeter Brune    DMDATSSetIJacobianLocal - set a local residual evaluation function
39124989b8cSPeter Brune 
39224989b8cSPeter Brune    Logically Collective
39324989b8cSPeter Brune 
39424989b8cSPeter Brune    Input Arguments:
39524989b8cSPeter Brune +  dm   - DM to associate callback with
39624989b8cSPeter Brune .  func - local residual evaluation
39724989b8cSPeter Brune -  ctx   - optional context for local residual evaluation
39824989b8cSPeter Brune 
39924989b8cSPeter Brune    Calling sequence for func:
400b398ffe5SPeter Brune 
401b398ffe5SPeter Brune $ func(DMDALocalInfo* info,PetscReal t,void* x,void *xdot,Mat J,Mat B,MatStructure *flg,void *ctx);
402b398ffe5SPeter Brune 
40324989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
404b398ffe5SPeter Brune .  t    - time at which to evaluate the jacobian
405b398ffe5SPeter Brune .  x    - array of local state information
406b398ffe5SPeter Brune .  xdot - time derivative at this state
407b398ffe5SPeter Brune .  J    - Jacobian matrix
408b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
409b398ffe5SPeter Brune .  flg  - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators())
41024989b8cSPeter Brune -  ctx  - optional context passed above
41124989b8cSPeter Brune 
41224989b8cSPeter Brune    Level: beginner
41324989b8cSPeter Brune 
4148e625dcfSPeter Brune .seealso: DMTSSetJacobian(), DMDATSSetIFunctionLocal(), DMDASNESSetJacobianLocal()
41524989b8cSPeter Brune @*/
4166c6b9e74SPeter Brune PetscErrorCode DMDATSSetIJacobianLocal(DM dm,DMDATSIJacobianLocal func,void *ctx)
41724989b8cSPeter Brune {
41824989b8cSPeter Brune   PetscErrorCode ierr;
41924989b8cSPeter Brune   TSDM         sdm;
42024989b8cSPeter Brune   DM_DA_TS     *dmdats;
42124989b8cSPeter Brune 
42224989b8cSPeter Brune   PetscFunctionBegin;
42324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
424*0313d0e0SJed Brown   ierr = DMTSGetContextWrite(dm,&sdm);CHKERRQ(ierr);
42524989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
42624989b8cSPeter Brune   dmdats->ijacobianlocal = func;
42724989b8cSPeter Brune   dmdats->ijacobianlocalctx = ctx;
42824989b8cSPeter Brune   ierr = DMTSSetIJacobian(dm,TSComputeIJacobian_DMDA,dmdats);CHKERRQ(ierr);
42924989b8cSPeter Brune   PetscFunctionReturn(0);
43024989b8cSPeter Brune }
431