xref: /petsc/src/ts/utils/dmdats.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
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;
17942e3340SBarry Smith } DMTS_DA;
1824989b8cSPeter Brune 
1924989b8cSPeter Brune #undef __FUNCT__
20942e3340SBarry Smith #define __FUNCT__ "DMTSDestroy_DMDA"
21942e3340SBarry Smith static PetscErrorCode DMTSDestroy_DMDA(DMTS 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__
31942e3340SBarry Smith #define __FUNCT__ "DMTSDuplicate_DMDA"
32d74926cbSBarry Smith static PetscErrorCode DMTSDuplicate_DMDA(DMTS oldsdm,DMTS sdm)
330313d0e0SJed Brown {
340313d0e0SJed Brown   PetscErrorCode ierr;
350313d0e0SJed Brown 
360313d0e0SJed Brown   PetscFunctionBegin;
37d74926cbSBarry Smith   ierr = PetscNewLog(sdm,DMTS_DA,&sdm->data);CHKERRQ(ierr);
38942e3340SBarry Smith   if (oldsdm->data) {ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMTS_DA));CHKERRQ(ierr);}
390313d0e0SJed Brown   PetscFunctionReturn(0);
400313d0e0SJed Brown }
410313d0e0SJed Brown 
420313d0e0SJed Brown #undef __FUNCT__
4324989b8cSPeter Brune #define __FUNCT__ "DMDATSGetContext"
44942e3340SBarry Smith static PetscErrorCode DMDATSGetContext(DM dm,DMTS sdm,DMTS_DA **dmdats)
4524989b8cSPeter Brune {
4624989b8cSPeter Brune   PetscErrorCode ierr;
4724989b8cSPeter Brune 
4824989b8cSPeter Brune   PetscFunctionBegin;
490298fd71SBarry Smith   *dmdats = NULL;
5024989b8cSPeter Brune   if (!sdm->data) {
51942e3340SBarry Smith     ierr = PetscNewLog(dm,DMTS_DA,&sdm->data);CHKERRQ(ierr);
52d74926cbSBarry Smith     sdm->ops->destroy   = DMTSDestroy_DMDA;
53d74926cbSBarry Smith     sdm->ops->duplicate = DMTSDuplicate_DMDA;
5424989b8cSPeter Brune   }
55942e3340SBarry Smith   *dmdats = (DMTS_DA*)sdm->data;
5624989b8cSPeter Brune   PetscFunctionReturn(0);
5724989b8cSPeter Brune }
5824989b8cSPeter Brune 
5924989b8cSPeter Brune #undef __FUNCT__
6024989b8cSPeter Brune #define __FUNCT__ "TSComputeIFunction_DMDA"
6124989b8cSPeter Brune /*
6224989b8cSPeter Brune   This function should eventually replace:
6324989b8cSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
6424989b8cSPeter Brune  */
6524989b8cSPeter Brune static PetscErrorCode TSComputeIFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *ctx)
6624989b8cSPeter Brune {
6724989b8cSPeter Brune   PetscErrorCode ierr;
6824989b8cSPeter Brune   DM             dm;
69942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)ctx;
7024989b8cSPeter Brune   DMDALocalInfo  info;
7124989b8cSPeter Brune   Vec            Xloc;
7224989b8cSPeter Brune   void           *x,*f,*xdot;
7324989b8cSPeter Brune 
7424989b8cSPeter Brune   PetscFunctionBegin;
7524989b8cSPeter Brune   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
7624989b8cSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
7724989b8cSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
78*ce94432eSBarry Smith   if (!dmdats->ifunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
7924989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
8024989b8cSPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
8124989b8cSPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
8224989b8cSPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
8324989b8cSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
8424989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
8524989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr);
8624989b8cSPeter Brune   switch (dmdats->ifunctionlocalimode) {
8724989b8cSPeter Brune   case INSERT_VALUES: {
8824989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
8924989b8cSPeter Brune     CHKMEMQ;
9024989b8cSPeter Brune     ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr);
9124989b8cSPeter Brune     CHKMEMQ;
9224989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
9324989b8cSPeter Brune   } break;
9424989b8cSPeter Brune   case ADD_VALUES: {
9524989b8cSPeter Brune     Vec Floc;
9624989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
9724989b8cSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
9824989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
9924989b8cSPeter Brune     CHKMEMQ;
10024989b8cSPeter Brune     ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr);
10124989b8cSPeter Brune     CHKMEMQ;
10224989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
10324989b8cSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
10424989b8cSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
10524989b8cSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
10624989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
10724989b8cSPeter Brune   } break;
108*ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->ifunctionlocalimode);
10924989b8cSPeter Brune   }
11024989b8cSPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
11124989b8cSPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
11211c2e1cfSJed Brown   ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr);
11324989b8cSPeter Brune   PetscFunctionReturn(0);
11424989b8cSPeter Brune }
11524989b8cSPeter Brune 
11624989b8cSPeter Brune #undef __FUNCT__
11724989b8cSPeter Brune #define __FUNCT__ "TSComputeIJacobian_DMDA"
11824989b8cSPeter Brune static PetscErrorCode TSComputeIJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
11924989b8cSPeter Brune {
12024989b8cSPeter Brune   PetscErrorCode ierr;
12124989b8cSPeter Brune   DM             dm;
122942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)ctx;
12324989b8cSPeter Brune   DMDALocalInfo  info;
12424989b8cSPeter Brune   Vec            Xloc;
12524989b8cSPeter Brune   void           *x,*xdot;
12624989b8cSPeter Brune 
12724989b8cSPeter Brune   PetscFunctionBegin;
128*ce94432eSBarry Smith   if (!dmdats->ifunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
12924989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
13024989b8cSPeter Brune 
13124989b8cSPeter Brune   if (dmdats->ijacobianlocal) {
13224989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
13324989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
13424989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
13524989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
13624989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
13724989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr);
13824989b8cSPeter Brune     CHKMEMQ;
13924989b8cSPeter Brune     ierr = (*dmdats->ijacobianlocal)(&info,ptime,x,xdot,shift,*A,*B,mstr,dmdats->ijacobianlocalctx);CHKERRQ(ierr);
14024989b8cSPeter Brune     CHKMEMQ;
14124989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
14224989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr);
14324989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
144*ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()");
14524989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
14624989b8cSPeter Brune   if (*A != *B) {
14724989b8cSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14824989b8cSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14924989b8cSPeter Brune   }
15024989b8cSPeter Brune   PetscFunctionReturn(0);
15124989b8cSPeter Brune }
15224989b8cSPeter Brune 
15324989b8cSPeter Brune #undef __FUNCT__
15424989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSFunction_DMDA"
15524989b8cSPeter Brune /*
15624989b8cSPeter Brune   This function should eventually replace:
15724989b8cSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
15824989b8cSPeter Brune  */
15924989b8cSPeter Brune static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec F,void *ctx)
16024989b8cSPeter Brune {
16124989b8cSPeter Brune   PetscErrorCode ierr;
16224989b8cSPeter Brune   DM             dm;
163942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)ctx;
16424989b8cSPeter Brune   DMDALocalInfo  info;
16524989b8cSPeter Brune   Vec            Xloc;
16624989b8cSPeter Brune   void           *x,*f;
16724989b8cSPeter Brune 
16824989b8cSPeter Brune   PetscFunctionBegin;
16924989b8cSPeter Brune   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
17024989b8cSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
17124989b8cSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
172*ce94432eSBarry Smith   if (!dmdats->rhsfunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
17324989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
17424989b8cSPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
17524989b8cSPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
17624989b8cSPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
17724989b8cSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
17824989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
17924989b8cSPeter Brune   switch (dmdats->rhsfunctionlocalimode) {
18024989b8cSPeter Brune   case INSERT_VALUES: {
18124989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
18224989b8cSPeter Brune     CHKMEMQ;
18324989b8cSPeter Brune     ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr);
18424989b8cSPeter Brune     CHKMEMQ;
18524989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
18624989b8cSPeter Brune   } break;
18724989b8cSPeter Brune   case ADD_VALUES: {
18824989b8cSPeter Brune     Vec Floc;
18924989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
19024989b8cSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
19124989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
19224989b8cSPeter Brune     CHKMEMQ;
19324989b8cSPeter Brune     ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr);
19424989b8cSPeter Brune     CHKMEMQ;
19524989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
19624989b8cSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
19724989b8cSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
19824989b8cSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
19924989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
20024989b8cSPeter Brune   } break;
201*ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->rhsfunctionlocalimode);
20224989b8cSPeter Brune   }
20324989b8cSPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
20424989b8cSPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
20524989b8cSPeter Brune   PetscFunctionReturn(0);
20624989b8cSPeter Brune }
20724989b8cSPeter Brune 
20824989b8cSPeter Brune #undef __FUNCT__
20924989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSJacobian_DMDA"
21024989b8cSPeter Brune static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
21124989b8cSPeter Brune {
21224989b8cSPeter Brune   PetscErrorCode ierr;
21324989b8cSPeter Brune   DM             dm;
214942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)ctx;
21524989b8cSPeter Brune   DMDALocalInfo  info;
21624989b8cSPeter Brune   Vec            Xloc;
21724989b8cSPeter Brune   void           *x;
21824989b8cSPeter Brune 
21924989b8cSPeter Brune   PetscFunctionBegin;
220*ce94432eSBarry Smith   if (!dmdats->rhsfunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
22124989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
22224989b8cSPeter Brune 
22324989b8cSPeter Brune   if (dmdats->rhsjacobianlocal) {
22424989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
22524989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
22624989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
22724989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
22824989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
22924989b8cSPeter Brune     CHKMEMQ;
23024989b8cSPeter Brune     ierr = (*dmdats->rhsjacobianlocal)(&info,ptime,x,*A,*B,mstr,dmdats->rhsjacobianlocalctx);CHKERRQ(ierr);
23124989b8cSPeter Brune     CHKMEMQ;
23224989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
23324989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
234*ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()");
23524989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
23624989b8cSPeter Brune   if (*A != *B) {
23724989b8cSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
23824989b8cSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
23924989b8cSPeter Brune   }
24024989b8cSPeter Brune   PetscFunctionReturn(0);
24124989b8cSPeter Brune }
24224989b8cSPeter Brune 
24324989b8cSPeter Brune 
24424989b8cSPeter Brune #undef __FUNCT__
24524989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSFunctionLocal"
24624989b8cSPeter Brune /*@C
24724989b8cSPeter Brune    DMDATSSetRHSFunctionLocal - set a local residual evaluation function
24824989b8cSPeter Brune 
24924989b8cSPeter Brune    Logically Collective
25024989b8cSPeter Brune 
25124989b8cSPeter Brune    Input Arguments:
25224989b8cSPeter Brune +  dm - DM to associate callback with
253b398ffe5SPeter Brune .  imode - insert mode for the residual
25424989b8cSPeter Brune .  func - local residual evaluation
25524989b8cSPeter Brune -  ctx - optional context for local residual evaluation
25624989b8cSPeter Brune 
25724989b8cSPeter Brune    Calling sequence for func:
258b398ffe5SPeter Brune 
259b398ffe5SPeter Brune $ func(DMDALocalInfo info,PetscReal t,void *x,void *f,void *ctx)
260b398ffe5SPeter Brune 
26124989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
262b398ffe5SPeter Brune .  t - time at which to evaluate residual
263b398ffe5SPeter Brune .  x - array of local state information
264b398ffe5SPeter Brune .  f - output array of local residual information
265b398ffe5SPeter Brune -  ctx - optional user context
26624989b8cSPeter Brune 
26724989b8cSPeter Brune    Level: beginner
26824989b8cSPeter Brune 
2698e625dcfSPeter Brune .seealso: DMTSSetRHSFunction(), DMDATSSetRHSJacobianLocal(), DMDASNESSetFunctionLocal()
27024989b8cSPeter Brune @*/
2716c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm,InsertMode imode,DMDATSRHSFunctionLocal func,void *ctx)
27224989b8cSPeter Brune {
27324989b8cSPeter Brune   PetscErrorCode ierr;
274942e3340SBarry Smith   DMTS           sdm;
275942e3340SBarry Smith   DMTS_DA        *dmdats;
27624989b8cSPeter Brune 
27724989b8cSPeter Brune   PetscFunctionBegin;
27824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
279942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
28024989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
28124989b8cSPeter Brune   dmdats->rhsfunctionlocalimode = imode;
28224989b8cSPeter Brune   dmdats->rhsfunctionlocal      = func;
28324989b8cSPeter Brune   dmdats->rhsfunctionlocalctx   = ctx;
28424989b8cSPeter Brune   ierr = DMTSSetRHSFunction(dm,TSComputeRHSFunction_DMDA,dmdats);CHKERRQ(ierr);
28524989b8cSPeter Brune   PetscFunctionReturn(0);
28624989b8cSPeter Brune }
28724989b8cSPeter Brune 
28824989b8cSPeter Brune #undef __FUNCT__
28924989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSJacobianLocal"
29024989b8cSPeter Brune /*@C
29124989b8cSPeter Brune    DMDATSSetRHSJacobianLocal - set a local residual evaluation function
29224989b8cSPeter Brune 
29324989b8cSPeter Brune    Logically Collective
29424989b8cSPeter Brune 
29524989b8cSPeter Brune    Input Arguments:
29624989b8cSPeter Brune +  dm    - DM to associate callback with
297b398ffe5SPeter Brune .  func  - local RHS Jacobian evaluation routine
298b398ffe5SPeter Brune -  ctx   - optional context for local jacobian evaluation
29924989b8cSPeter Brune 
30024989b8cSPeter Brune    Calling sequence for func:
301b398ffe5SPeter Brune 
302b398ffe5SPeter Brune $ func(DMDALocalInfo* info,PetscReal t,void* x,Mat J,Mat B,MatStructure *flg,void *ctx);
303b398ffe5SPeter Brune 
30424989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
305b398ffe5SPeter Brune .  t    - time at which to evaluate residual
306b398ffe5SPeter Brune .  x    - array of local state information
307b398ffe5SPeter Brune .  J    - Jacobian matrix
308b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
309b398ffe5SPeter Brune .  flg  - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators())
31024989b8cSPeter Brune -  ctx  - optional context passed above
31124989b8cSPeter Brune 
31224989b8cSPeter Brune    Level: beginner
31324989b8cSPeter Brune 
3148e625dcfSPeter Brune .seealso: DMTSSetRHSJacobian(), DMDATSSetRHSFunctionLocal(), DMDASNESSetJacobianLocal()
31524989b8cSPeter Brune @*/
3166c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm,DMDATSRHSJacobianLocal func,void *ctx)
31724989b8cSPeter Brune {
31824989b8cSPeter Brune   PetscErrorCode ierr;
319942e3340SBarry Smith   DMTS           sdm;
320942e3340SBarry Smith   DMTS_DA        *dmdats;
32124989b8cSPeter Brune 
32224989b8cSPeter Brune   PetscFunctionBegin;
32324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
324942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
32524989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
32624989b8cSPeter Brune   dmdats->rhsjacobianlocal    = func;
32724989b8cSPeter Brune   dmdats->rhsjacobianlocalctx = ctx;
32824989b8cSPeter Brune   ierr = DMTSSetRHSJacobian(dm,TSComputeRHSJacobian_DMDA,dmdats);CHKERRQ(ierr);
32924989b8cSPeter Brune   PetscFunctionReturn(0);
33024989b8cSPeter Brune }
33124989b8cSPeter Brune 
33224989b8cSPeter Brune 
33324989b8cSPeter Brune #undef __FUNCT__
33424989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIFunctionLocal"
33524989b8cSPeter Brune /*@C
33624989b8cSPeter Brune    DMDATSSetIFunctionLocal - set a local residual evaluation function
33724989b8cSPeter Brune 
33824989b8cSPeter Brune    Logically Collective
33924989b8cSPeter Brune 
34024989b8cSPeter Brune    Input Arguments:
34124989b8cSPeter Brune +  dm   - DM to associate callback with
34224989b8cSPeter Brune .  func - local residual evaluation
34324989b8cSPeter Brune -  ctx  - optional context for local residual evaluation
34424989b8cSPeter Brune 
34524989b8cSPeter Brune    Calling sequence for func:
34624989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
347b398ffe5SPeter Brune .  t    - time at which to evaluate residual
348b398ffe5SPeter Brune .  x    - array of local state information
349b398ffe5SPeter Brune .  xdot - array of local time derivative information
350b398ffe5SPeter Brune .  f    - output array of local function evaluation information
35124989b8cSPeter Brune -  ctx - optional context passed above
35224989b8cSPeter Brune 
35324989b8cSPeter Brune    Level: beginner
35424989b8cSPeter Brune 
3558e625dcfSPeter Brune .seealso: DMTSSetIFunction(), DMDATSSetIJacobianLocal(), DMDASNESSetFunctionLocal()
35624989b8cSPeter Brune @*/
3576c6b9e74SPeter Brune PetscErrorCode DMDATSSetIFunctionLocal(DM dm,InsertMode imode,DMDATSIFunctionLocal func,void *ctx)
35824989b8cSPeter Brune {
35924989b8cSPeter Brune   PetscErrorCode ierr;
360942e3340SBarry Smith   DMTS           sdm;
361942e3340SBarry Smith   DMTS_DA        *dmdats;
36224989b8cSPeter Brune 
36324989b8cSPeter Brune   PetscFunctionBegin;
36424989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
365942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
36624989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
36724989b8cSPeter Brune   dmdats->ifunctionlocalimode = imode;
36824989b8cSPeter Brune   dmdats->ifunctionlocal      = func;
36924989b8cSPeter Brune   dmdats->ifunctionlocalctx   = ctx;
37024989b8cSPeter Brune   ierr = DMTSSetIFunction(dm,TSComputeIFunction_DMDA,dmdats);CHKERRQ(ierr);
37124989b8cSPeter Brune   PetscFunctionReturn(0);
37224989b8cSPeter Brune }
37324989b8cSPeter Brune 
37424989b8cSPeter Brune #undef __FUNCT__
37524989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIJacobianLocal"
37624989b8cSPeter Brune /*@C
37724989b8cSPeter Brune    DMDATSSetIJacobianLocal - set a local residual evaluation function
37824989b8cSPeter Brune 
37924989b8cSPeter Brune    Logically Collective
38024989b8cSPeter Brune 
38124989b8cSPeter Brune    Input Arguments:
38224989b8cSPeter Brune +  dm   - DM to associate callback with
38324989b8cSPeter Brune .  func - local residual evaluation
38424989b8cSPeter Brune -  ctx   - optional context for local residual evaluation
38524989b8cSPeter Brune 
38624989b8cSPeter Brune    Calling sequence for func:
387b398ffe5SPeter Brune 
388b398ffe5SPeter Brune $ func(DMDALocalInfo* info,PetscReal t,void* x,void *xdot,Mat J,Mat B,MatStructure *flg,void *ctx);
389b398ffe5SPeter Brune 
39024989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
391b398ffe5SPeter Brune .  t    - time at which to evaluate the jacobian
392b398ffe5SPeter Brune .  x    - array of local state information
393b398ffe5SPeter Brune .  xdot - time derivative at this state
394b398ffe5SPeter Brune .  J    - Jacobian matrix
395b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
396b398ffe5SPeter Brune .  flg  - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators())
39724989b8cSPeter Brune -  ctx  - optional context passed above
39824989b8cSPeter Brune 
39924989b8cSPeter Brune    Level: beginner
40024989b8cSPeter Brune 
4018e625dcfSPeter Brune .seealso: DMTSSetJacobian(), DMDATSSetIFunctionLocal(), DMDASNESSetJacobianLocal()
40224989b8cSPeter Brune @*/
4036c6b9e74SPeter Brune PetscErrorCode DMDATSSetIJacobianLocal(DM dm,DMDATSIJacobianLocal func,void *ctx)
40424989b8cSPeter Brune {
40524989b8cSPeter Brune   PetscErrorCode ierr;
406942e3340SBarry Smith   DMTS           sdm;
407942e3340SBarry Smith   DMTS_DA        *dmdats;
40824989b8cSPeter Brune 
40924989b8cSPeter Brune   PetscFunctionBegin;
41024989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
411942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
41224989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
41324989b8cSPeter Brune   dmdats->ijacobianlocal    = func;
41424989b8cSPeter Brune   dmdats->ijacobianlocalctx = ctx;
41524989b8cSPeter Brune   ierr = DMTSSetIJacobian(dm,TSComputeIJacobian_DMDA,dmdats);CHKERRQ(ierr);
41624989b8cSPeter Brune   PetscFunctionReturn(0);
41724989b8cSPeter Brune }
418d1212d36SBarry Smith 
419d1212d36SBarry Smith #undef __FUNCT__
420d1212d36SBarry Smith #define __FUNCT__ "TSMonitorDMDARayDestroy"
421d1212d36SBarry Smith PetscErrorCode TSMonitorDMDARayDestroy(void **mctx)
422d1212d36SBarry Smith {
423d1212d36SBarry Smith   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx*)*mctx;
424d1212d36SBarry Smith   PetscErrorCode      ierr;
425d1212d36SBarry Smith 
426d1212d36SBarry Smith   PetscFunctionBegin;
427d1212d36SBarry Smith   ierr = VecDestroy(&rayctx->ray);CHKERRQ(ierr);
428d1212d36SBarry Smith   ierr = VecScatterDestroy(&rayctx->scatter);CHKERRQ(ierr);
429d1212d36SBarry Smith   ierr = PetscViewerDestroy(&rayctx->viewer);CHKERRQ(ierr);
430d1212d36SBarry Smith   ierr = PetscFree(rayctx);CHKERRQ(ierr);
431d1212d36SBarry Smith   PetscFunctionReturn(0);
432d1212d36SBarry Smith }
433d1212d36SBarry Smith 
434d1212d36SBarry Smith #undef __FUNCT__
435d1212d36SBarry Smith #define __FUNCT__ "TSMonitorDMDARay"
436d1212d36SBarry Smith PetscErrorCode TSMonitorDMDARay(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx)
437d1212d36SBarry Smith {
438d1212d36SBarry Smith   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx*)mctx;
439d1212d36SBarry Smith   Vec                 solution;
440d1212d36SBarry Smith   PetscErrorCode      ierr;
441d1212d36SBarry Smith 
442d1212d36SBarry Smith   PetscFunctionBegin;
443d1212d36SBarry Smith   ierr = TSGetSolution(ts,&solution);CHKERRQ(ierr);
444d1212d36SBarry Smith   ierr = VecScatterBegin(rayctx->scatter,solution,rayctx->ray,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
445d1212d36SBarry Smith   ierr = VecScatterEnd(rayctx->scatter,solution,rayctx->ray,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
446d1212d36SBarry Smith   if (rayctx->viewer) {
447d1212d36SBarry Smith     ierr = VecView(rayctx->ray,rayctx->viewer);CHKERRQ(ierr);
448d1212d36SBarry Smith   }
449d1212d36SBarry Smith   PetscFunctionReturn(0);
450d1212d36SBarry Smith }
451