xref: /petsc/src/ts/utils/dmdats.c (revision 942e334039a1a47ff909c032afdf0581bd880e9f)
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;
17*942e3340SBarry Smith } DMTS_DA;
1824989b8cSPeter Brune 
1924989b8cSPeter Brune #undef __FUNCT__
20*942e3340SBarry Smith #define __FUNCT__ "DMTSDestroy_DMDA"
21*942e3340SBarry 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__
31*942e3340SBarry Smith #define __FUNCT__ "DMTSDuplicate_DMDA"
32*942e3340SBarry Smith static PetscErrorCode DMTSDuplicate_DMDA(DMTS oldsdm,DM dm)
330313d0e0SJed Brown {
340313d0e0SJed Brown   PetscErrorCode ierr;
35*942e3340SBarry Smith   DMTS           sdm;
360313d0e0SJed Brown 
370313d0e0SJed Brown   PetscFunctionBegin;
38*942e3340SBarry Smith   ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr);
39*942e3340SBarry Smith   ierr = PetscNewLog(dm,DMTS_DA,&sdm->data);CHKERRQ(ierr);
40*942e3340SBarry Smith   if (oldsdm->data) {ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMTS_DA));CHKERRQ(ierr);}
410313d0e0SJed Brown   PetscFunctionReturn(0);
420313d0e0SJed Brown }
430313d0e0SJed Brown 
440313d0e0SJed Brown #undef __FUNCT__
4524989b8cSPeter Brune #define __FUNCT__ "DMDATSGetContext"
46*942e3340SBarry Smith static PetscErrorCode DMDATSGetContext(DM dm,DMTS sdm,DMTS_DA **dmdats)
4724989b8cSPeter Brune {
4824989b8cSPeter Brune   PetscErrorCode ierr;
4924989b8cSPeter Brune 
5024989b8cSPeter Brune   PetscFunctionBegin;
5124989b8cSPeter Brune   *dmdats = PETSC_NULL;
5224989b8cSPeter Brune   if (!sdm->data) {
53*942e3340SBarry Smith     ierr = PetscNewLog(dm,DMTS_DA,&sdm->data);CHKERRQ(ierr);
54*942e3340SBarry Smith     sdm->destroy = DMTSDestroy_DMDA;
55*942e3340SBarry Smith     sdm->duplicate = DMTSDuplicate_DMDA;
5624989b8cSPeter Brune   }
57*942e3340SBarry Smith   *dmdats = (DMTS_DA*)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;
71*942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)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 static PetscErrorCode TSComputeIJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
12024989b8cSPeter Brune {
12124989b8cSPeter Brune   PetscErrorCode ierr;
12224989b8cSPeter Brune   DM             dm;
123*942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)ctx;
12424989b8cSPeter Brune   DMDALocalInfo  info;
12524989b8cSPeter Brune   Vec            Xloc;
12624989b8cSPeter Brune   void           *x,*xdot;
12724989b8cSPeter Brune 
12824989b8cSPeter Brune   PetscFunctionBegin;
12924989b8cSPeter Brune   if (!dmdats->ifunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
13024989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
13124989b8cSPeter Brune 
13224989b8cSPeter Brune   if (dmdats->ijacobianlocal) {
13324989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
13424989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
13524989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
13624989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
13724989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
13824989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr);
13924989b8cSPeter Brune     CHKMEMQ;
14024989b8cSPeter Brune     ierr = (*dmdats->ijacobianlocal)(&info,ptime,x,xdot,shift,*A,*B,mstr,dmdats->ijacobianlocalctx);CHKERRQ(ierr);
14124989b8cSPeter Brune     CHKMEMQ;
14224989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
14324989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr);
14424989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
14524989b8cSPeter Brune   } else {
14624989b8cSPeter Brune     SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()");
14724989b8cSPeter Brune   }
14824989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
14924989b8cSPeter Brune   if (*A != *B) {
15024989b8cSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15124989b8cSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15224989b8cSPeter Brune   }
15324989b8cSPeter Brune   PetscFunctionReturn(0);
15424989b8cSPeter Brune }
15524989b8cSPeter Brune 
15624989b8cSPeter Brune #undef __FUNCT__
15724989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSFunction_DMDA"
15824989b8cSPeter Brune /*
15924989b8cSPeter Brune   This function should eventually replace:
16024989b8cSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
16124989b8cSPeter Brune  */
16224989b8cSPeter Brune static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec F,void *ctx)
16324989b8cSPeter Brune {
16424989b8cSPeter Brune   PetscErrorCode ierr;
16524989b8cSPeter Brune   DM             dm;
166*942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)ctx;
16724989b8cSPeter Brune   DMDALocalInfo  info;
16824989b8cSPeter Brune   Vec            Xloc;
16924989b8cSPeter Brune   void           *x,*f;
17024989b8cSPeter Brune 
17124989b8cSPeter Brune   PetscFunctionBegin;
17224989b8cSPeter Brune   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
17324989b8cSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
17424989b8cSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
17524989b8cSPeter Brune   if (!dmdats->rhsfunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
17624989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
17724989b8cSPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
17824989b8cSPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
17924989b8cSPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
18024989b8cSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
18124989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
18224989b8cSPeter Brune   switch (dmdats->rhsfunctionlocalimode) {
18324989b8cSPeter Brune   case INSERT_VALUES: {
18424989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
18524989b8cSPeter Brune     CHKMEMQ;
18624989b8cSPeter Brune     ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr);
18724989b8cSPeter Brune     CHKMEMQ;
18824989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
18924989b8cSPeter Brune   } break;
19024989b8cSPeter Brune   case ADD_VALUES: {
19124989b8cSPeter Brune     Vec Floc;
19224989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
19324989b8cSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
19424989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
19524989b8cSPeter Brune     CHKMEMQ;
19624989b8cSPeter Brune     ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr);
19724989b8cSPeter Brune     CHKMEMQ;
19824989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
19924989b8cSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
20024989b8cSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
20124989b8cSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
20224989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
20324989b8cSPeter Brune   } break;
20424989b8cSPeter Brune   default: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->rhsfunctionlocalimode);
20524989b8cSPeter Brune   }
20624989b8cSPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
20724989b8cSPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
20824989b8cSPeter Brune   PetscFunctionReturn(0);
20924989b8cSPeter Brune }
21024989b8cSPeter Brune 
21124989b8cSPeter Brune #undef __FUNCT__
21224989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSJacobian_DMDA"
21324989b8cSPeter Brune static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
21424989b8cSPeter Brune {
21524989b8cSPeter Brune   PetscErrorCode ierr;
21624989b8cSPeter Brune   DM             dm;
217*942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)ctx;
21824989b8cSPeter Brune   DMDALocalInfo  info;
21924989b8cSPeter Brune   Vec            Xloc;
22024989b8cSPeter Brune   void           *x;
22124989b8cSPeter Brune 
22224989b8cSPeter Brune   PetscFunctionBegin;
22324989b8cSPeter Brune   if (!dmdats->rhsfunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
22424989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
22524989b8cSPeter Brune 
22624989b8cSPeter Brune   if (dmdats->rhsjacobianlocal) {
22724989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
22824989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
22924989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
23024989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
23124989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
23224989b8cSPeter Brune     CHKMEMQ;
23324989b8cSPeter Brune     ierr = (*dmdats->rhsjacobianlocal)(&info,ptime,x,*A,*B,mstr,dmdats->rhsjacobianlocalctx);CHKERRQ(ierr);
23424989b8cSPeter Brune     CHKMEMQ;
23524989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
23624989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
23724989b8cSPeter Brune   } else {
23824989b8cSPeter Brune     SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()");
23924989b8cSPeter Brune   }
24024989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
24124989b8cSPeter Brune   if (*A != *B) {
24224989b8cSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
24324989b8cSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
24424989b8cSPeter Brune   }
24524989b8cSPeter Brune   PetscFunctionReturn(0);
24624989b8cSPeter Brune }
24724989b8cSPeter Brune 
24824989b8cSPeter Brune 
24924989b8cSPeter Brune #undef __FUNCT__
25024989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSFunctionLocal"
25124989b8cSPeter Brune /*@C
25224989b8cSPeter Brune    DMDATSSetRHSFunctionLocal - set a local residual evaluation function
25324989b8cSPeter Brune 
25424989b8cSPeter Brune    Logically Collective
25524989b8cSPeter Brune 
25624989b8cSPeter Brune    Input Arguments:
25724989b8cSPeter Brune +  dm - DM to associate callback with
258b398ffe5SPeter Brune .  imode - insert mode for the residual
25924989b8cSPeter Brune .  func - local residual evaluation
26024989b8cSPeter Brune -  ctx - optional context for local residual evaluation
26124989b8cSPeter Brune 
26224989b8cSPeter Brune    Calling sequence for func:
263b398ffe5SPeter Brune 
264b398ffe5SPeter Brune $ func(DMDALocalInfo info,PetscReal t,void *x,void *f,void *ctx)
265b398ffe5SPeter Brune 
26624989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
267b398ffe5SPeter Brune .  t - time at which to evaluate residual
268b398ffe5SPeter Brune .  x - array of local state information
269b398ffe5SPeter Brune .  f - output array of local residual information
270b398ffe5SPeter Brune -  ctx - optional user context
27124989b8cSPeter Brune 
27224989b8cSPeter Brune    Level: beginner
27324989b8cSPeter Brune 
2748e625dcfSPeter Brune .seealso: DMTSSetRHSFunction(), DMDATSSetRHSJacobianLocal(), DMDASNESSetFunctionLocal()
27524989b8cSPeter Brune @*/
2766c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm,InsertMode imode,DMDATSRHSFunctionLocal func,void *ctx)
27724989b8cSPeter Brune {
27824989b8cSPeter Brune   PetscErrorCode ierr;
279*942e3340SBarry Smith   DMTS           sdm;
280*942e3340SBarry Smith   DMTS_DA        *dmdats;
28124989b8cSPeter Brune 
28224989b8cSPeter Brune   PetscFunctionBegin;
28324989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
284*942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
28524989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
28624989b8cSPeter Brune   dmdats->rhsfunctionlocalimode = imode;
28724989b8cSPeter Brune   dmdats->rhsfunctionlocal = func;
28824989b8cSPeter Brune   dmdats->rhsfunctionlocalctx = ctx;
28924989b8cSPeter Brune   ierr = DMTSSetRHSFunction(dm,TSComputeRHSFunction_DMDA,dmdats);CHKERRQ(ierr);
29024989b8cSPeter Brune   PetscFunctionReturn(0);
29124989b8cSPeter Brune }
29224989b8cSPeter Brune 
29324989b8cSPeter Brune #undef __FUNCT__
29424989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSJacobianLocal"
29524989b8cSPeter Brune /*@C
29624989b8cSPeter Brune    DMDATSSetRHSJacobianLocal - set a local residual evaluation function
29724989b8cSPeter Brune 
29824989b8cSPeter Brune    Logically Collective
29924989b8cSPeter Brune 
30024989b8cSPeter Brune    Input Arguments:
30124989b8cSPeter Brune +  dm    - DM to associate callback with
302b398ffe5SPeter Brune .  func  - local RHS Jacobian evaluation routine
303b398ffe5SPeter Brune -  ctx   - optional context for local jacobian evaluation
30424989b8cSPeter Brune 
30524989b8cSPeter Brune    Calling sequence for func:
306b398ffe5SPeter Brune 
307b398ffe5SPeter Brune $ func(DMDALocalInfo* info,PetscReal t,void* x,Mat J,Mat B,MatStructure *flg,void *ctx);
308b398ffe5SPeter Brune 
30924989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
310b398ffe5SPeter Brune .  t    - time at which to evaluate residual
311b398ffe5SPeter Brune .  x    - array of local state information
312b398ffe5SPeter Brune .  J    - Jacobian matrix
313b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
314b398ffe5SPeter Brune .  flg  - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators())
31524989b8cSPeter Brune -  ctx  - optional context passed above
31624989b8cSPeter Brune 
31724989b8cSPeter Brune    Level: beginner
31824989b8cSPeter Brune 
3198e625dcfSPeter Brune .seealso: DMTSSetRHSJacobian(), DMDATSSetRHSFunctionLocal(), DMDASNESSetJacobianLocal()
32024989b8cSPeter Brune @*/
3216c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm,DMDATSRHSJacobianLocal func,void *ctx)
32224989b8cSPeter Brune {
32324989b8cSPeter Brune   PetscErrorCode ierr;
324*942e3340SBarry Smith   DMTS           sdm;
325*942e3340SBarry Smith   DMTS_DA        *dmdats;
32624989b8cSPeter Brune 
32724989b8cSPeter Brune   PetscFunctionBegin;
32824989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
329*942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
33024989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
33124989b8cSPeter Brune   dmdats->rhsjacobianlocal = func;
33224989b8cSPeter Brune   dmdats->rhsjacobianlocalctx = ctx;
33324989b8cSPeter Brune   ierr = DMTSSetRHSJacobian(dm,TSComputeRHSJacobian_DMDA,dmdats);CHKERRQ(ierr);
33424989b8cSPeter Brune   PetscFunctionReturn(0);
33524989b8cSPeter Brune }
33624989b8cSPeter Brune 
33724989b8cSPeter Brune 
33824989b8cSPeter Brune #undef __FUNCT__
33924989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIFunctionLocal"
34024989b8cSPeter Brune /*@C
34124989b8cSPeter Brune    DMDATSSetIFunctionLocal - set a local residual evaluation function
34224989b8cSPeter Brune 
34324989b8cSPeter Brune    Logically Collective
34424989b8cSPeter Brune 
34524989b8cSPeter Brune    Input Arguments:
34624989b8cSPeter Brune +  dm   - DM to associate callback with
34724989b8cSPeter Brune .  func - local residual evaluation
34824989b8cSPeter Brune -  ctx  - optional context for local residual evaluation
34924989b8cSPeter Brune 
35024989b8cSPeter Brune    Calling sequence for func:
35124989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
352b398ffe5SPeter Brune .  t    - time at which to evaluate residual
353b398ffe5SPeter Brune .  x    - array of local state information
354b398ffe5SPeter Brune .  xdot - array of local time derivative information
355b398ffe5SPeter Brune .  f    - output array of local function evaluation information
35624989b8cSPeter Brune -  ctx - optional context passed above
35724989b8cSPeter Brune 
35824989b8cSPeter Brune    Level: beginner
35924989b8cSPeter Brune 
3608e625dcfSPeter Brune .seealso: DMTSSetIFunction(), DMDATSSetIJacobianLocal(), DMDASNESSetFunctionLocal()
36124989b8cSPeter Brune @*/
3626c6b9e74SPeter Brune PetscErrorCode DMDATSSetIFunctionLocal(DM dm,InsertMode imode,DMDATSIFunctionLocal func,void *ctx)
36324989b8cSPeter Brune {
36424989b8cSPeter Brune   PetscErrorCode ierr;
365*942e3340SBarry Smith   DMTS           sdm;
366*942e3340SBarry Smith   DMTS_DA        *dmdats;
36724989b8cSPeter Brune 
36824989b8cSPeter Brune   PetscFunctionBegin;
36924989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
370*942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
37124989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
37224989b8cSPeter Brune   dmdats->ifunctionlocalimode = imode;
37324989b8cSPeter Brune   dmdats->ifunctionlocal = func;
37424989b8cSPeter Brune   dmdats->ifunctionlocalctx = ctx;
37524989b8cSPeter Brune   ierr = DMTSSetIFunction(dm,TSComputeIFunction_DMDA,dmdats);CHKERRQ(ierr);
37624989b8cSPeter Brune   PetscFunctionReturn(0);
37724989b8cSPeter Brune }
37824989b8cSPeter Brune 
37924989b8cSPeter Brune #undef __FUNCT__
38024989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIJacobianLocal"
38124989b8cSPeter Brune /*@C
38224989b8cSPeter Brune    DMDATSSetIJacobianLocal - set a local residual evaluation function
38324989b8cSPeter Brune 
38424989b8cSPeter Brune    Logically Collective
38524989b8cSPeter Brune 
38624989b8cSPeter Brune    Input Arguments:
38724989b8cSPeter Brune +  dm   - DM to associate callback with
38824989b8cSPeter Brune .  func - local residual evaluation
38924989b8cSPeter Brune -  ctx   - optional context for local residual evaluation
39024989b8cSPeter Brune 
39124989b8cSPeter Brune    Calling sequence for func:
392b398ffe5SPeter Brune 
393b398ffe5SPeter Brune $ func(DMDALocalInfo* info,PetscReal t,void* x,void *xdot,Mat J,Mat B,MatStructure *flg,void *ctx);
394b398ffe5SPeter Brune 
39524989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
396b398ffe5SPeter Brune .  t    - time at which to evaluate the jacobian
397b398ffe5SPeter Brune .  x    - array of local state information
398b398ffe5SPeter Brune .  xdot - time derivative at this state
399b398ffe5SPeter Brune .  J    - Jacobian matrix
400b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
401b398ffe5SPeter Brune .  flg  - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators())
40224989b8cSPeter Brune -  ctx  - optional context passed above
40324989b8cSPeter Brune 
40424989b8cSPeter Brune    Level: beginner
40524989b8cSPeter Brune 
4068e625dcfSPeter Brune .seealso: DMTSSetJacobian(), DMDATSSetIFunctionLocal(), DMDASNESSetJacobianLocal()
40724989b8cSPeter Brune @*/
4086c6b9e74SPeter Brune PetscErrorCode DMDATSSetIJacobianLocal(DM dm,DMDATSIJacobianLocal func,void *ctx)
40924989b8cSPeter Brune {
41024989b8cSPeter Brune   PetscErrorCode ierr;
411*942e3340SBarry Smith   DMTS           sdm;
412*942e3340SBarry Smith   DMTS_DA        *dmdats;
41324989b8cSPeter Brune 
41424989b8cSPeter Brune   PetscFunctionBegin;
41524989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
416*942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
41724989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
41824989b8cSPeter Brune   dmdats->ijacobianlocal = func;
41924989b8cSPeter Brune   dmdats->ijacobianlocalctx = ctx;
42024989b8cSPeter Brune   ierr = DMTSSetIJacobian(dm,TSComputeIJacobian_DMDA,dmdats);CHKERRQ(ierr);
42124989b8cSPeter Brune   PetscFunctionReturn(0);
42224989b8cSPeter Brune }
423