xref: /petsc/src/ts/utils/dmdats.c (revision dab4423f8acb25e3442fb80c59fbbbd498575e6e)
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*/
451b4a12fSMatthew G. Knepley #include <petscdraw.h>
524989b8cSPeter Brune 
624989b8cSPeter Brune /* This structure holds the user-provided DMDA callbacks */
724989b8cSPeter Brune typedef struct {
824989b8cSPeter Brune   PetscErrorCode (*ifunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
924989b8cSPeter Brune   PetscErrorCode (*rhsfunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
10d1e9a80fSBarry Smith   PetscErrorCode (*ijacobianlocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*);
11d1e9a80fSBarry Smith   PetscErrorCode (*rhsjacobianlocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*);
1224989b8cSPeter Brune   void       *ifunctionlocalctx;
1324989b8cSPeter Brune   void       *ijacobianlocalctx;
1424989b8cSPeter Brune   void       *rhsfunctionlocalctx;
1524989b8cSPeter Brune   void       *rhsjacobianlocalctx;
1624989b8cSPeter Brune   InsertMode ifunctionlocalimode;
1724989b8cSPeter Brune   InsertMode rhsfunctionlocalimode;
18942e3340SBarry Smith } DMTS_DA;
1924989b8cSPeter Brune 
2024989b8cSPeter Brune #undef __FUNCT__
21942e3340SBarry Smith #define __FUNCT__ "DMTSDestroy_DMDA"
22942e3340SBarry Smith static PetscErrorCode DMTSDestroy_DMDA(DMTS sdm)
2324989b8cSPeter Brune {
2424989b8cSPeter Brune   PetscErrorCode ierr;
2524989b8cSPeter Brune 
2624989b8cSPeter Brune   PetscFunctionBegin;
2724989b8cSPeter Brune   ierr = PetscFree(sdm->data);CHKERRQ(ierr);
2824989b8cSPeter Brune   PetscFunctionReturn(0);
2924989b8cSPeter Brune }
3024989b8cSPeter Brune 
3124989b8cSPeter Brune #undef __FUNCT__
32942e3340SBarry Smith #define __FUNCT__ "DMTSDuplicate_DMDA"
33d74926cbSBarry Smith static PetscErrorCode DMTSDuplicate_DMDA(DMTS oldsdm,DMTS sdm)
340313d0e0SJed Brown {
350313d0e0SJed Brown   PetscErrorCode ierr;
360313d0e0SJed Brown 
370313d0e0SJed Brown   PetscFunctionBegin;
38b00a9115SJed Brown   ierr = PetscNewLog(sdm,(DMTS_DA**)&sdm->data);CHKERRQ(ierr);
39942e3340SBarry Smith   if (oldsdm->data) {ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMTS_DA));CHKERRQ(ierr);}
400313d0e0SJed Brown   PetscFunctionReturn(0);
410313d0e0SJed Brown }
420313d0e0SJed Brown 
430313d0e0SJed Brown #undef __FUNCT__
4424989b8cSPeter Brune #define __FUNCT__ "DMDATSGetContext"
45942e3340SBarry Smith static PetscErrorCode DMDATSGetContext(DM dm,DMTS sdm,DMTS_DA **dmdats)
4624989b8cSPeter Brune {
4724989b8cSPeter Brune   PetscErrorCode ierr;
4824989b8cSPeter Brune 
4924989b8cSPeter Brune   PetscFunctionBegin;
500298fd71SBarry Smith   *dmdats = NULL;
5124989b8cSPeter Brune   if (!sdm->data) {
52b00a9115SJed Brown     ierr = PetscNewLog(dm,(DMTS_DA**)&sdm->data);CHKERRQ(ierr);
53d74926cbSBarry Smith     sdm->ops->destroy   = DMTSDestroy_DMDA;
54d74926cbSBarry Smith     sdm->ops->duplicate = DMTSDuplicate_DMDA;
5524989b8cSPeter Brune   }
56942e3340SBarry Smith   *dmdats = (DMTS_DA*)sdm->data;
5724989b8cSPeter Brune   PetscFunctionReturn(0);
5824989b8cSPeter Brune }
5924989b8cSPeter Brune 
6024989b8cSPeter Brune #undef __FUNCT__
6124989b8cSPeter Brune #define __FUNCT__ "TSComputeIFunction_DMDA"
6224989b8cSPeter Brune static PetscErrorCode TSComputeIFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *ctx)
6324989b8cSPeter Brune {
6424989b8cSPeter Brune   PetscErrorCode ierr;
6524989b8cSPeter Brune   DM             dm;
66942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)ctx;
6724989b8cSPeter Brune   DMDALocalInfo  info;
6824989b8cSPeter Brune   Vec            Xloc;
6924989b8cSPeter Brune   void           *x,*f,*xdot;
7024989b8cSPeter Brune 
7124989b8cSPeter Brune   PetscFunctionBegin;
7224989b8cSPeter Brune   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
7324989b8cSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
7424989b8cSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
75ce94432eSBarry Smith   if (!dmdats->ifunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
7624989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
7724989b8cSPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
7824989b8cSPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
7924989b8cSPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
8024989b8cSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
8124989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
8224989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr);
8324989b8cSPeter Brune   switch (dmdats->ifunctionlocalimode) {
8424989b8cSPeter Brune   case INSERT_VALUES: {
8524989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,F,&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,F,&f);CHKERRQ(ierr);
9024989b8cSPeter Brune   } break;
9124989b8cSPeter Brune   case ADD_VALUES: {
9224989b8cSPeter Brune     Vec Floc;
9324989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
9424989b8cSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
9524989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
9624989b8cSPeter Brune     CHKMEMQ;
9724989b8cSPeter Brune     ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr);
9824989b8cSPeter Brune     CHKMEMQ;
9924989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
10024989b8cSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
10124989b8cSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
10224989b8cSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
10324989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
10424989b8cSPeter Brune   } break;
105ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->ifunctionlocalimode);
10624989b8cSPeter Brune   }
10724989b8cSPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
10824989b8cSPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
10911c2e1cfSJed Brown   ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr);
11024989b8cSPeter Brune   PetscFunctionReturn(0);
11124989b8cSPeter Brune }
11224989b8cSPeter Brune 
11324989b8cSPeter Brune #undef __FUNCT__
11424989b8cSPeter Brune #define __FUNCT__ "TSComputeIJacobian_DMDA"
115d1e9a80fSBarry Smith static PetscErrorCode TSComputeIJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,PetscReal shift,Mat A,Mat B,void *ctx)
11624989b8cSPeter Brune {
11724989b8cSPeter Brune   PetscErrorCode ierr;
11824989b8cSPeter Brune   DM             dm;
119942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)ctx;
12024989b8cSPeter Brune   DMDALocalInfo  info;
12124989b8cSPeter Brune   Vec            Xloc;
12224989b8cSPeter Brune   void           *x,*xdot;
12324989b8cSPeter Brune 
12424989b8cSPeter Brune   PetscFunctionBegin;
125ce94432eSBarry Smith   if (!dmdats->ifunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
12624989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
12724989b8cSPeter Brune 
12824989b8cSPeter Brune   if (dmdats->ijacobianlocal) {
12924989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
13024989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
13124989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
13224989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
13324989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
13424989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr);
13524989b8cSPeter Brune     CHKMEMQ;
136d1e9a80fSBarry Smith     ierr = (*dmdats->ijacobianlocal)(&info,ptime,x,xdot,shift,A,B,dmdats->ijacobianlocalctx);CHKERRQ(ierr);
13724989b8cSPeter Brune     CHKMEMQ;
13824989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
13924989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr);
14024989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
141ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()");
14224989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
14394ab13aaSBarry Smith   if (A != B) {
14494ab13aaSBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14594ab13aaSBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14624989b8cSPeter Brune   }
14724989b8cSPeter Brune   PetscFunctionReturn(0);
14824989b8cSPeter Brune }
14924989b8cSPeter Brune 
15024989b8cSPeter Brune #undef __FUNCT__
15124989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSFunction_DMDA"
15224989b8cSPeter Brune static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec F,void *ctx)
15324989b8cSPeter Brune {
15424989b8cSPeter Brune   PetscErrorCode ierr;
15524989b8cSPeter Brune   DM             dm;
156942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)ctx;
15724989b8cSPeter Brune   DMDALocalInfo  info;
15824989b8cSPeter Brune   Vec            Xloc;
15924989b8cSPeter Brune   void           *x,*f;
16024989b8cSPeter Brune 
16124989b8cSPeter Brune   PetscFunctionBegin;
16224989b8cSPeter Brune   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
16324989b8cSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
16424989b8cSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
165ce94432eSBarry Smith   if (!dmdats->rhsfunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
16624989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
16724989b8cSPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
16824989b8cSPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
16924989b8cSPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
17024989b8cSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
17124989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
17224989b8cSPeter Brune   switch (dmdats->rhsfunctionlocalimode) {
17324989b8cSPeter Brune   case INSERT_VALUES: {
17424989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
17524989b8cSPeter Brune     CHKMEMQ;
17624989b8cSPeter Brune     ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr);
17724989b8cSPeter Brune     CHKMEMQ;
17824989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
17924989b8cSPeter Brune   } break;
18024989b8cSPeter Brune   case ADD_VALUES: {
18124989b8cSPeter Brune     Vec Floc;
18224989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
18324989b8cSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
18424989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&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,Floc,&f);CHKERRQ(ierr);
18924989b8cSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
19024989b8cSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
19124989b8cSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
19224989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
19324989b8cSPeter Brune   } break;
194ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->rhsfunctionlocalimode);
19524989b8cSPeter Brune   }
19624989b8cSPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
19724989b8cSPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
19824989b8cSPeter Brune   PetscFunctionReturn(0);
19924989b8cSPeter Brune }
20024989b8cSPeter Brune 
20124989b8cSPeter Brune #undef __FUNCT__
20224989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSJacobian_DMDA"
203d1e9a80fSBarry Smith static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Mat A,Mat B,void *ctx)
20424989b8cSPeter Brune {
20524989b8cSPeter Brune   PetscErrorCode ierr;
20624989b8cSPeter Brune   DM             dm;
207942e3340SBarry Smith   DMTS_DA        *dmdats = (DMTS_DA*)ctx;
20824989b8cSPeter Brune   DMDALocalInfo  info;
20924989b8cSPeter Brune   Vec            Xloc;
21024989b8cSPeter Brune   void           *x;
21124989b8cSPeter Brune 
21224989b8cSPeter Brune   PetscFunctionBegin;
213ce94432eSBarry Smith   if (!dmdats->rhsfunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
21424989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
21524989b8cSPeter Brune 
21624989b8cSPeter Brune   if (dmdats->rhsjacobianlocal) {
21724989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
21824989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
21924989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
22024989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
22124989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
22224989b8cSPeter Brune     CHKMEMQ;
223d1e9a80fSBarry Smith     ierr = (*dmdats->rhsjacobianlocal)(&info,ptime,x,A,B,dmdats->rhsjacobianlocalctx);CHKERRQ(ierr);
22424989b8cSPeter Brune     CHKMEMQ;
22524989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
22624989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
227ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()");
22824989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
22994ab13aaSBarry Smith   if (A != B) {
23094ab13aaSBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
23194ab13aaSBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
23224989b8cSPeter Brune   }
23324989b8cSPeter Brune   PetscFunctionReturn(0);
23424989b8cSPeter Brune }
23524989b8cSPeter Brune 
23624989b8cSPeter Brune 
23724989b8cSPeter Brune #undef __FUNCT__
23824989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSFunctionLocal"
23924989b8cSPeter Brune /*@C
24024989b8cSPeter Brune    DMDATSSetRHSFunctionLocal - set a local residual evaluation function
24124989b8cSPeter Brune 
24224989b8cSPeter Brune    Logically Collective
24324989b8cSPeter Brune 
24424989b8cSPeter Brune    Input Arguments:
24524989b8cSPeter Brune +  dm - DM to associate callback with
246b398ffe5SPeter Brune .  imode - insert mode for the residual
24724989b8cSPeter Brune .  func - local residual evaluation
24824989b8cSPeter Brune -  ctx - optional context for local residual evaluation
24924989b8cSPeter Brune 
25024989b8cSPeter Brune    Calling sequence for func:
251b398ffe5SPeter Brune 
252b398ffe5SPeter Brune $ func(DMDALocalInfo info,PetscReal t,void *x,void *f,void *ctx)
253b398ffe5SPeter Brune 
25424989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
255b398ffe5SPeter Brune .  t - time at which to evaluate residual
256b398ffe5SPeter Brune .  x - array of local state information
257b398ffe5SPeter Brune .  f - output array of local residual information
258b398ffe5SPeter Brune -  ctx - optional user context
25924989b8cSPeter Brune 
26024989b8cSPeter Brune    Level: beginner
26124989b8cSPeter Brune 
2628e625dcfSPeter Brune .seealso: DMTSSetRHSFunction(), DMDATSSetRHSJacobianLocal(), DMDASNESSetFunctionLocal()
26324989b8cSPeter Brune @*/
2646c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm,InsertMode imode,DMDATSRHSFunctionLocal func,void *ctx)
26524989b8cSPeter Brune {
26624989b8cSPeter Brune   PetscErrorCode ierr;
267942e3340SBarry Smith   DMTS           sdm;
268942e3340SBarry Smith   DMTS_DA        *dmdats;
26924989b8cSPeter Brune 
27024989b8cSPeter Brune   PetscFunctionBegin;
27124989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
272942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
27324989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
27424989b8cSPeter Brune   dmdats->rhsfunctionlocalimode = imode;
27524989b8cSPeter Brune   dmdats->rhsfunctionlocal      = func;
27624989b8cSPeter Brune   dmdats->rhsfunctionlocalctx   = ctx;
27724989b8cSPeter Brune   ierr = DMTSSetRHSFunction(dm,TSComputeRHSFunction_DMDA,dmdats);CHKERRQ(ierr);
27824989b8cSPeter Brune   PetscFunctionReturn(0);
27924989b8cSPeter Brune }
28024989b8cSPeter Brune 
28124989b8cSPeter Brune #undef __FUNCT__
28224989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSJacobianLocal"
28324989b8cSPeter Brune /*@C
28424989b8cSPeter Brune    DMDATSSetRHSJacobianLocal - set a local residual evaluation function
28524989b8cSPeter Brune 
28624989b8cSPeter Brune    Logically Collective
28724989b8cSPeter Brune 
28824989b8cSPeter Brune    Input Arguments:
28924989b8cSPeter Brune +  dm    - DM to associate callback with
290b398ffe5SPeter Brune .  func  - local RHS Jacobian evaluation routine
291b398ffe5SPeter Brune -  ctx   - optional context for local jacobian evaluation
29224989b8cSPeter Brune 
29324989b8cSPeter Brune    Calling sequence for func:
294b398ffe5SPeter Brune 
2957f6be190SBarry Smith $ func(DMDALocalInfo* info,PetscReal t,void* x,Mat J,Mat B,void *ctx);
296b398ffe5SPeter Brune 
29724989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
298b398ffe5SPeter Brune .  t    - time at which to evaluate residual
299b398ffe5SPeter Brune .  x    - array of local state information
300b398ffe5SPeter Brune .  J    - Jacobian matrix
301b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
30224989b8cSPeter Brune -  ctx  - optional context passed above
30324989b8cSPeter Brune 
30424989b8cSPeter Brune    Level: beginner
30524989b8cSPeter Brune 
3068e625dcfSPeter Brune .seealso: DMTSSetRHSJacobian(), DMDATSSetRHSFunctionLocal(), DMDASNESSetJacobianLocal()
30724989b8cSPeter Brune @*/
3086c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm,DMDATSRHSJacobianLocal func,void *ctx)
30924989b8cSPeter Brune {
31024989b8cSPeter Brune   PetscErrorCode ierr;
311942e3340SBarry Smith   DMTS           sdm;
312942e3340SBarry Smith   DMTS_DA        *dmdats;
31324989b8cSPeter Brune 
31424989b8cSPeter Brune   PetscFunctionBegin;
31524989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
316942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
31724989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
31824989b8cSPeter Brune   dmdats->rhsjacobianlocal    = func;
31924989b8cSPeter Brune   dmdats->rhsjacobianlocalctx = ctx;
32024989b8cSPeter Brune   ierr = DMTSSetRHSJacobian(dm,TSComputeRHSJacobian_DMDA,dmdats);CHKERRQ(ierr);
32124989b8cSPeter Brune   PetscFunctionReturn(0);
32224989b8cSPeter Brune }
32324989b8cSPeter Brune 
32424989b8cSPeter Brune 
32524989b8cSPeter Brune #undef __FUNCT__
32624989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIFunctionLocal"
32724989b8cSPeter Brune /*@C
32824989b8cSPeter Brune    DMDATSSetIFunctionLocal - set a local residual evaluation function
32924989b8cSPeter Brune 
33024989b8cSPeter Brune    Logically Collective
33124989b8cSPeter Brune 
33224989b8cSPeter Brune    Input Arguments:
33324989b8cSPeter Brune +  dm   - DM to associate callback with
33424989b8cSPeter Brune .  func - local residual evaluation
33524989b8cSPeter Brune -  ctx  - optional context for local residual evaluation
33624989b8cSPeter Brune 
33724989b8cSPeter Brune    Calling sequence for func:
33824989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
339b398ffe5SPeter Brune .  t    - time at which to evaluate residual
340b398ffe5SPeter Brune .  x    - array of local state information
341b398ffe5SPeter Brune .  xdot - array of local time derivative information
342b398ffe5SPeter Brune .  f    - output array of local function evaluation information
34324989b8cSPeter Brune -  ctx - optional context passed above
34424989b8cSPeter Brune 
34524989b8cSPeter Brune    Level: beginner
34624989b8cSPeter Brune 
3478e625dcfSPeter Brune .seealso: DMTSSetIFunction(), DMDATSSetIJacobianLocal(), DMDASNESSetFunctionLocal()
34824989b8cSPeter Brune @*/
3496c6b9e74SPeter Brune PetscErrorCode DMDATSSetIFunctionLocal(DM dm,InsertMode imode,DMDATSIFunctionLocal func,void *ctx)
35024989b8cSPeter Brune {
35124989b8cSPeter Brune   PetscErrorCode ierr;
352942e3340SBarry Smith   DMTS           sdm;
353942e3340SBarry Smith   DMTS_DA        *dmdats;
35424989b8cSPeter Brune 
35524989b8cSPeter Brune   PetscFunctionBegin;
35624989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
357942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
35824989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
35924989b8cSPeter Brune   dmdats->ifunctionlocalimode = imode;
36024989b8cSPeter Brune   dmdats->ifunctionlocal      = func;
36124989b8cSPeter Brune   dmdats->ifunctionlocalctx   = ctx;
36224989b8cSPeter Brune   ierr = DMTSSetIFunction(dm,TSComputeIFunction_DMDA,dmdats);CHKERRQ(ierr);
36324989b8cSPeter Brune   PetscFunctionReturn(0);
36424989b8cSPeter Brune }
36524989b8cSPeter Brune 
36624989b8cSPeter Brune #undef __FUNCT__
36724989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIJacobianLocal"
36824989b8cSPeter Brune /*@C
36924989b8cSPeter Brune    DMDATSSetIJacobianLocal - set a local residual evaluation function
37024989b8cSPeter Brune 
37124989b8cSPeter Brune    Logically Collective
37224989b8cSPeter Brune 
37324989b8cSPeter Brune    Input Arguments:
37424989b8cSPeter Brune +  dm   - DM to associate callback with
37524989b8cSPeter Brune .  func - local residual evaluation
37624989b8cSPeter Brune -  ctx   - optional context for local residual evaluation
37724989b8cSPeter Brune 
37824989b8cSPeter Brune    Calling sequence for func:
379b398ffe5SPeter Brune 
380*dab4423fSBarry Smith $ func(DMDALocalInfo* info,PetscReal t,void* x,void *xdot,PetscScalar shift,Mat J,Mat B,void *ctx);
381b398ffe5SPeter Brune 
38224989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
383b398ffe5SPeter Brune .  t    - time at which to evaluate the jacobian
384b398ffe5SPeter Brune .  x    - array of local state information
385b398ffe5SPeter Brune .  xdot - time derivative at this state
386*dab4423fSBarry Smith .  shift - see TSSetIJacobian() for the meaning of this parameter
387b398ffe5SPeter Brune .  J    - Jacobian matrix
388b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
38924989b8cSPeter Brune -  ctx  - optional context passed above
39024989b8cSPeter Brune 
39124989b8cSPeter Brune    Level: beginner
39224989b8cSPeter Brune 
3938e625dcfSPeter Brune .seealso: DMTSSetJacobian(), DMDATSSetIFunctionLocal(), DMDASNESSetJacobianLocal()
39424989b8cSPeter Brune @*/
3956c6b9e74SPeter Brune PetscErrorCode DMDATSSetIJacobianLocal(DM dm,DMDATSIJacobianLocal func,void *ctx)
39624989b8cSPeter Brune {
39724989b8cSPeter Brune   PetscErrorCode ierr;
398942e3340SBarry Smith   DMTS           sdm;
399942e3340SBarry Smith   DMTS_DA        *dmdats;
40024989b8cSPeter Brune 
40124989b8cSPeter Brune   PetscFunctionBegin;
40224989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
403942e3340SBarry Smith   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
40424989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
40524989b8cSPeter Brune   dmdats->ijacobianlocal    = func;
40624989b8cSPeter Brune   dmdats->ijacobianlocalctx = ctx;
40724989b8cSPeter Brune   ierr = DMTSSetIJacobian(dm,TSComputeIJacobian_DMDA,dmdats);CHKERRQ(ierr);
40824989b8cSPeter Brune   PetscFunctionReturn(0);
40924989b8cSPeter Brune }
410d1212d36SBarry Smith 
411d1212d36SBarry Smith #undef __FUNCT__
412d1212d36SBarry Smith #define __FUNCT__ "TSMonitorDMDARayDestroy"
413d1212d36SBarry Smith PetscErrorCode TSMonitorDMDARayDestroy(void **mctx)
414d1212d36SBarry Smith {
415d1212d36SBarry Smith   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *) *mctx;
416d1212d36SBarry Smith   PetscErrorCode       ierr;
417d1212d36SBarry Smith 
418d1212d36SBarry Smith   PetscFunctionBegin;
41951b4a12fSMatthew G. Knepley   if (rayctx->lgctx) {ierr = TSMonitorLGCtxDestroy(&rayctx->lgctx);CHKERRQ(ierr);}
420d1212d36SBarry Smith   ierr = VecDestroy(&rayctx->ray);CHKERRQ(ierr);
421d1212d36SBarry Smith   ierr = VecScatterDestroy(&rayctx->scatter);CHKERRQ(ierr);
422d1212d36SBarry Smith   ierr = PetscViewerDestroy(&rayctx->viewer);CHKERRQ(ierr);
423d1212d36SBarry Smith   ierr = PetscFree(rayctx);CHKERRQ(ierr);
424d1212d36SBarry Smith   PetscFunctionReturn(0);
425d1212d36SBarry Smith }
426d1212d36SBarry Smith 
427d1212d36SBarry Smith #undef __FUNCT__
428d1212d36SBarry Smith #define __FUNCT__ "TSMonitorDMDARay"
429d1212d36SBarry Smith PetscErrorCode TSMonitorDMDARay(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx)
430d1212d36SBarry Smith {
431d1212d36SBarry Smith   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx*)mctx;
432d1212d36SBarry Smith   Vec                 solution;
433d1212d36SBarry Smith   PetscErrorCode      ierr;
434d1212d36SBarry Smith 
435d1212d36SBarry Smith   PetscFunctionBegin;
436d1212d36SBarry Smith   ierr = TSGetSolution(ts,&solution);CHKERRQ(ierr);
437d1212d36SBarry Smith   ierr = VecScatterBegin(rayctx->scatter,solution,rayctx->ray,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
438d1212d36SBarry Smith   ierr = VecScatterEnd(rayctx->scatter,solution,rayctx->ray,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
439d1212d36SBarry Smith   if (rayctx->viewer) {
440d1212d36SBarry Smith     ierr = VecView(rayctx->ray,rayctx->viewer);CHKERRQ(ierr);
441d1212d36SBarry Smith   }
442d1212d36SBarry Smith   PetscFunctionReturn(0);
443d1212d36SBarry Smith }
44451b4a12fSMatthew G. Knepley 
44551b4a12fSMatthew G. Knepley #undef __FUNCT__
44651b4a12fSMatthew G. Knepley #define __FUNCT__ "TSMonitorLGDMDARay"
44751b4a12fSMatthew G. Knepley PetscErrorCode  TSMonitorLGDMDARay(TS ts, PetscInt step, PetscReal ptime, Vec u, void *ctx)
44851b4a12fSMatthew G. Knepley {
44951b4a12fSMatthew G. Knepley   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *) ctx;
45051b4a12fSMatthew G. Knepley   TSMonitorLGCtx       lgctx  = (TSMonitorLGCtx) rayctx->lgctx;
45151b4a12fSMatthew G. Knepley   Vec                  v      = rayctx->ray;
45251b4a12fSMatthew G. Knepley   const PetscScalar   *a;
45351b4a12fSMatthew G. Knepley   PetscInt             dim;
45451b4a12fSMatthew G. Knepley   PetscErrorCode       ierr;
45551b4a12fSMatthew G. Knepley 
45651b4a12fSMatthew G. Knepley   PetscFunctionBegin;
45751b4a12fSMatthew G. Knepley   ierr = VecScatterBegin(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
45851b4a12fSMatthew G. Knepley   ierr = VecScatterEnd(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
45951b4a12fSMatthew G. Knepley   if (!step) {
46051b4a12fSMatthew G. Knepley     PetscDrawAxis axis;
46151b4a12fSMatthew G. Knepley 
46251b4a12fSMatthew G. Knepley     ierr = PetscDrawLGGetAxis(lgctx->lg, &axis);CHKERRQ(ierr);
46351b4a12fSMatthew G. Knepley     ierr = PetscDrawAxisSetLabels(axis, "Solution Ray as function of time", "Time", "Solution");CHKERRQ(ierr);
46451b4a12fSMatthew G. Knepley     ierr = VecGetLocalSize(rayctx->ray, &dim);CHKERRQ(ierr);
46551b4a12fSMatthew G. Knepley     ierr = PetscDrawLGSetDimension(lgctx->lg, dim);CHKERRQ(ierr);
46651b4a12fSMatthew G. Knepley     ierr = PetscDrawLGReset(lgctx->lg);CHKERRQ(ierr);
46751b4a12fSMatthew G. Knepley   }
46851b4a12fSMatthew G. Knepley   ierr = VecGetArrayRead(v, &a);CHKERRQ(ierr);
46951b4a12fSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX)
47051b4a12fSMatthew G. Knepley   {
47151b4a12fSMatthew G. Knepley     PetscReal *areal;
47251b4a12fSMatthew G. Knepley     PetscInt   i,n;
47351b4a12fSMatthew G. Knepley     ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
474785e854fSJed Brown     ierr = PetscMalloc1(n, &areal);CHKERRQ(ierr);
47551b4a12fSMatthew G. Knepley     for (i = 0; i < n; ++i) areal[i] = PetscRealPart(a[i]);
47651b4a12fSMatthew G. Knepley     ierr = PetscDrawLGAddCommonPoint(lgctx->lg, ptime, areal);CHKERRQ(ierr);
47751b4a12fSMatthew G. Knepley     ierr = PetscFree(areal);CHKERRQ(ierr);
47851b4a12fSMatthew G. Knepley   }
47951b4a12fSMatthew G. Knepley #else
48051b4a12fSMatthew G. Knepley   ierr = PetscDrawLGAddCommonPoint(lgctx->lg, ptime, a);CHKERRQ(ierr);
48151b4a12fSMatthew G. Knepley #endif
48251b4a12fSMatthew G. Knepley   ierr = VecRestoreArrayRead(v, &a);CHKERRQ(ierr);
48351b4a12fSMatthew G. Knepley   if (((lgctx->howoften > 0) && (!(step % lgctx->howoften))) || ((lgctx->howoften == -1) && ts->reason)) {
48451b4a12fSMatthew G. Knepley     ierr = PetscDrawLGDraw(lgctx->lg);CHKERRQ(ierr);
48551b4a12fSMatthew G. Knepley   }
48651b4a12fSMatthew G. Knepley   PetscFunctionReturn(0);
48751b4a12fSMatthew G. Knepley }
488