xref: /petsc/src/ts/utils/dmdats.c (revision 24989b8cf3e165ebaa37f4b2d30151131f137c58)
1*24989b8cSPeter Brune #include <petscdmda.h>          /*I "petscdmda.h" I*/
2*24989b8cSPeter Brune #include <petsc-private/dmimpl.h>
3*24989b8cSPeter Brune #include <petsc-private/tsimpl.h>   /*I "petscts.h" I*/
4*24989b8cSPeter Brune 
5*24989b8cSPeter Brune /* This structure holds the user-provided DMDA callbacks */
6*24989b8cSPeter Brune typedef struct {
7*24989b8cSPeter Brune   PetscErrorCode (*ifunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
8*24989b8cSPeter Brune   PetscErrorCode (*rhsfunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
9*24989b8cSPeter Brune   PetscErrorCode (*ijacobianlocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*);
10*24989b8cSPeter Brune   PetscErrorCode (*rhsjacobianlocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*);
11*24989b8cSPeter Brune   void *ifunctionlocalctx;
12*24989b8cSPeter Brune   void *ijacobianlocalctx;
13*24989b8cSPeter Brune   void *rhsfunctionlocalctx;
14*24989b8cSPeter Brune   void *rhsjacobianlocalctx;
15*24989b8cSPeter Brune   InsertMode ifunctionlocalimode;
16*24989b8cSPeter Brune   InsertMode rhsfunctionlocalimode;
17*24989b8cSPeter Brune } DM_DA_TS;
18*24989b8cSPeter Brune 
19*24989b8cSPeter Brune #undef __FUNCT__
20*24989b8cSPeter Brune #define __FUNCT__ "TSDMDestroy_DMDA"
21*24989b8cSPeter Brune static PetscErrorCode TSDMDestroy_DMDA(TSDM sdm)
22*24989b8cSPeter Brune {
23*24989b8cSPeter Brune   PetscErrorCode ierr;
24*24989b8cSPeter Brune 
25*24989b8cSPeter Brune   PetscFunctionBegin;
26*24989b8cSPeter Brune   ierr = PetscFree(sdm->data);CHKERRQ(ierr);
27*24989b8cSPeter Brune   PetscFunctionReturn(0);
28*24989b8cSPeter Brune }
29*24989b8cSPeter Brune 
30*24989b8cSPeter Brune #undef __FUNCT__
31*24989b8cSPeter Brune #define __FUNCT__ "DMDATSGetContext"
32*24989b8cSPeter Brune static PetscErrorCode DMDATSGetContext(DM dm,TSDM sdm,DM_DA_TS **dmdats)
33*24989b8cSPeter Brune {
34*24989b8cSPeter Brune   PetscErrorCode ierr;
35*24989b8cSPeter Brune 
36*24989b8cSPeter Brune   PetscFunctionBegin;
37*24989b8cSPeter Brune   *dmdats = PETSC_NULL;
38*24989b8cSPeter Brune   if (!sdm->data) {
39*24989b8cSPeter Brune     ierr = PetscNewLog(dm,DM_DA_TS,&sdm->data);CHKERRQ(ierr);
40*24989b8cSPeter Brune     sdm->destroy = TSDMDestroy_DMDA;
41*24989b8cSPeter Brune   }
42*24989b8cSPeter Brune   *dmdats = (DM_DA_TS*)sdm->data;
43*24989b8cSPeter Brune   PetscFunctionReturn(0);
44*24989b8cSPeter Brune }
45*24989b8cSPeter Brune 
46*24989b8cSPeter Brune #undef __FUNCT__
47*24989b8cSPeter Brune #define __FUNCT__ "TSComputeIFunction_DMDA"
48*24989b8cSPeter Brune /*
49*24989b8cSPeter Brune   This function should eventually replace:
50*24989b8cSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
51*24989b8cSPeter Brune  */
52*24989b8cSPeter Brune static PetscErrorCode TSComputeIFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *ctx)
53*24989b8cSPeter Brune {
54*24989b8cSPeter Brune   PetscErrorCode ierr;
55*24989b8cSPeter Brune   DM             dm;
56*24989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
57*24989b8cSPeter Brune   DMDALocalInfo  info;
58*24989b8cSPeter Brune   Vec            Xloc;
59*24989b8cSPeter Brune   void           *x,*f,*xdot;
60*24989b8cSPeter Brune 
61*24989b8cSPeter Brune   PetscFunctionBegin;
62*24989b8cSPeter Brune   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
63*24989b8cSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
64*24989b8cSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
65*24989b8cSPeter Brune   if (!dmdats->ifunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
66*24989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
67*24989b8cSPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
68*24989b8cSPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
69*24989b8cSPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
70*24989b8cSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
71*24989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
72*24989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr);
73*24989b8cSPeter Brune   switch (dmdats->ifunctionlocalimode) {
74*24989b8cSPeter Brune   case INSERT_VALUES: {
75*24989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
76*24989b8cSPeter Brune     CHKMEMQ;
77*24989b8cSPeter Brune     ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr);
78*24989b8cSPeter Brune     CHKMEMQ;
79*24989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
80*24989b8cSPeter Brune   } break;
81*24989b8cSPeter Brune   case ADD_VALUES: {
82*24989b8cSPeter Brune     Vec Floc;
83*24989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
84*24989b8cSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
85*24989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
86*24989b8cSPeter Brune     CHKMEMQ;
87*24989b8cSPeter Brune     ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr);
88*24989b8cSPeter Brune     CHKMEMQ;
89*24989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
90*24989b8cSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
91*24989b8cSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
92*24989b8cSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
93*24989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
94*24989b8cSPeter Brune   } break;
95*24989b8cSPeter Brune   default: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->ifunctionlocalimode);
96*24989b8cSPeter Brune   }
97*24989b8cSPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
98*24989b8cSPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
99*24989b8cSPeter Brune   PetscFunctionReturn(0);
100*24989b8cSPeter Brune }
101*24989b8cSPeter Brune 
102*24989b8cSPeter Brune #undef __FUNCT__
103*24989b8cSPeter Brune #define __FUNCT__ "TSComputeIJacobian_DMDA"
104*24989b8cSPeter Brune /*
105*24989b8cSPeter Brune   This function should eventually replace:
106*24989b8cSPeter Brune     DMComputeJacobian() and DMDAComputeJacobian1()
107*24989b8cSPeter Brune  */
108*24989b8cSPeter Brune static PetscErrorCode TSComputeIJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
109*24989b8cSPeter Brune {
110*24989b8cSPeter Brune   PetscErrorCode ierr;
111*24989b8cSPeter Brune   DM             dm;
112*24989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
113*24989b8cSPeter Brune   DMDALocalInfo  info;
114*24989b8cSPeter Brune   Vec            Xloc;
115*24989b8cSPeter Brune   void           *x,*xdot;
116*24989b8cSPeter Brune 
117*24989b8cSPeter Brune   PetscFunctionBegin;
118*24989b8cSPeter Brune   if (!dmdats->ifunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
119*24989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
120*24989b8cSPeter Brune 
121*24989b8cSPeter Brune   if (dmdats->ijacobianlocal) {
122*24989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
123*24989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
124*24989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
125*24989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
126*24989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
127*24989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr);
128*24989b8cSPeter Brune     CHKMEMQ;
129*24989b8cSPeter Brune     ierr = (*dmdats->ijacobianlocal)(&info,ptime,x,xdot,shift,*A,*B,mstr,dmdats->ijacobianlocalctx);CHKERRQ(ierr);
130*24989b8cSPeter Brune     CHKMEMQ;
131*24989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
132*24989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr);
133*24989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
134*24989b8cSPeter Brune   } else {
135*24989b8cSPeter Brune     SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()");
136*24989b8cSPeter Brune   }
137*24989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
138*24989b8cSPeter Brune   if (*A != *B) {
139*24989b8cSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
140*24989b8cSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
141*24989b8cSPeter Brune   }
142*24989b8cSPeter Brune   PetscFunctionReturn(0);
143*24989b8cSPeter Brune }
144*24989b8cSPeter Brune 
145*24989b8cSPeter Brune #undef __FUNCT__
146*24989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSFunction_DMDA"
147*24989b8cSPeter Brune /*
148*24989b8cSPeter Brune   This function should eventually replace:
149*24989b8cSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
150*24989b8cSPeter Brune  */
151*24989b8cSPeter Brune static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec F,void *ctx)
152*24989b8cSPeter Brune {
153*24989b8cSPeter Brune   PetscErrorCode ierr;
154*24989b8cSPeter Brune   DM             dm;
155*24989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
156*24989b8cSPeter Brune   DMDALocalInfo  info;
157*24989b8cSPeter Brune   Vec            Xloc;
158*24989b8cSPeter Brune   void           *x,*f;
159*24989b8cSPeter Brune 
160*24989b8cSPeter Brune   PetscFunctionBegin;
161*24989b8cSPeter Brune   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
162*24989b8cSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
163*24989b8cSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
164*24989b8cSPeter Brune   if (!dmdats->rhsfunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
165*24989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
166*24989b8cSPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
167*24989b8cSPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
168*24989b8cSPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
169*24989b8cSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
170*24989b8cSPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
171*24989b8cSPeter Brune   switch (dmdats->rhsfunctionlocalimode) {
172*24989b8cSPeter Brune   case INSERT_VALUES: {
173*24989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
174*24989b8cSPeter Brune     CHKMEMQ;
175*24989b8cSPeter Brune     ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr);
176*24989b8cSPeter Brune     CHKMEMQ;
177*24989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
178*24989b8cSPeter Brune   } break;
179*24989b8cSPeter Brune   case ADD_VALUES: {
180*24989b8cSPeter Brune     Vec Floc;
181*24989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
182*24989b8cSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
183*24989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
184*24989b8cSPeter Brune     CHKMEMQ;
185*24989b8cSPeter Brune     ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr);
186*24989b8cSPeter Brune     CHKMEMQ;
187*24989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
188*24989b8cSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
189*24989b8cSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
190*24989b8cSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
191*24989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
192*24989b8cSPeter Brune   } break;
193*24989b8cSPeter Brune   default: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->rhsfunctionlocalimode);
194*24989b8cSPeter Brune   }
195*24989b8cSPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
196*24989b8cSPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
197*24989b8cSPeter Brune   PetscFunctionReturn(0);
198*24989b8cSPeter Brune }
199*24989b8cSPeter Brune 
200*24989b8cSPeter Brune #undef __FUNCT__
201*24989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSJacobian_DMDA"
202*24989b8cSPeter Brune /*
203*24989b8cSPeter Brune   This function should eventually replace:
204*24989b8cSPeter Brune     DMComputeJacobian() and DMDAComputeJacobian1()
205*24989b8cSPeter Brune  */
206*24989b8cSPeter Brune static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
207*24989b8cSPeter Brune {
208*24989b8cSPeter Brune   PetscErrorCode ierr;
209*24989b8cSPeter Brune   DM             dm;
210*24989b8cSPeter Brune   DM_DA_TS     *dmdats = (DM_DA_TS*)ctx;
211*24989b8cSPeter Brune   DMDALocalInfo  info;
212*24989b8cSPeter Brune   Vec            Xloc;
213*24989b8cSPeter Brune   void           *x;
214*24989b8cSPeter Brune 
215*24989b8cSPeter Brune   PetscFunctionBegin;
216*24989b8cSPeter Brune   if (!dmdats->rhsfunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context");
217*24989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
218*24989b8cSPeter Brune 
219*24989b8cSPeter Brune   if (dmdats->rhsjacobianlocal) {
220*24989b8cSPeter Brune     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
221*24989b8cSPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
222*24989b8cSPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
223*24989b8cSPeter Brune     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
224*24989b8cSPeter Brune     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
225*24989b8cSPeter Brune     CHKMEMQ;
226*24989b8cSPeter Brune     ierr = (*dmdats->rhsjacobianlocal)(&info,ptime,x,*A,*B,mstr,dmdats->rhsjacobianlocalctx);CHKERRQ(ierr);
227*24989b8cSPeter Brune     CHKMEMQ;
228*24989b8cSPeter Brune     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
229*24989b8cSPeter Brune     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
230*24989b8cSPeter Brune   } else {
231*24989b8cSPeter Brune     SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()");
232*24989b8cSPeter Brune   }
233*24989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
234*24989b8cSPeter Brune   if (*A != *B) {
235*24989b8cSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
236*24989b8cSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
237*24989b8cSPeter Brune   }
238*24989b8cSPeter Brune   PetscFunctionReturn(0);
239*24989b8cSPeter Brune }
240*24989b8cSPeter Brune 
241*24989b8cSPeter Brune 
242*24989b8cSPeter Brune #undef __FUNCT__
243*24989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSFunctionLocal"
244*24989b8cSPeter Brune /*@C
245*24989b8cSPeter Brune    DMDATSSetRHSFunctionLocal - set a local residual evaluation function
246*24989b8cSPeter Brune 
247*24989b8cSPeter Brune    Logically Collective
248*24989b8cSPeter Brune 
249*24989b8cSPeter Brune    Input Arguments:
250*24989b8cSPeter Brune +  dm - DM to associate callback with
251*24989b8cSPeter Brune .  func - local residual evaluation
252*24989b8cSPeter Brune -  ctx - optional context for local residual evaluation
253*24989b8cSPeter Brune 
254*24989b8cSPeter Brune    Calling sequence for func:
255*24989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
256*24989b8cSPeter Brune .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
257*24989b8cSPeter Brune .  x - dimensional pointer to state at which to evaluate residual
258*24989b8cSPeter Brune .  f - dimensional pointer to residual, write the residual here
259*24989b8cSPeter Brune -  ctx - optional context passed above
260*24989b8cSPeter Brune 
261*24989b8cSPeter Brune    Level: beginner
262*24989b8cSPeter Brune 
263*24989b8cSPeter Brune .seealso: DMTSSetFunction(), DMDATSSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
264*24989b8cSPeter Brune @*/
265*24989b8cSPeter Brune PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *ctx)
266*24989b8cSPeter Brune {
267*24989b8cSPeter Brune   PetscErrorCode ierr;
268*24989b8cSPeter Brune   TSDM         sdm;
269*24989b8cSPeter Brune   DM_DA_TS     *dmdats;
270*24989b8cSPeter Brune 
271*24989b8cSPeter Brune   PetscFunctionBegin;
272*24989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
273*24989b8cSPeter Brune   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
274*24989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
275*24989b8cSPeter Brune   dmdats->rhsfunctionlocalimode = imode;
276*24989b8cSPeter Brune   dmdats->rhsfunctionlocal = func;
277*24989b8cSPeter Brune   dmdats->rhsfunctionlocalctx = ctx;
278*24989b8cSPeter Brune   ierr = DMTSSetRHSFunction(dm,TSComputeRHSFunction_DMDA,dmdats);CHKERRQ(ierr);
279*24989b8cSPeter Brune   PetscFunctionReturn(0);
280*24989b8cSPeter Brune }
281*24989b8cSPeter Brune 
282*24989b8cSPeter Brune #undef __FUNCT__
283*24989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSJacobianLocal"
284*24989b8cSPeter Brune /*@C
285*24989b8cSPeter Brune    DMDATSSetRHSJacobianLocal - set a local residual evaluation function
286*24989b8cSPeter Brune 
287*24989b8cSPeter Brune    Logically Collective
288*24989b8cSPeter Brune 
289*24989b8cSPeter Brune    Input Arguments:
290*24989b8cSPeter Brune +  dm - DM to associate callback with
291*24989b8cSPeter Brune .  func - local residual evaluation
292*24989b8cSPeter Brune -  ctx - optional context for local residual evaluation
293*24989b8cSPeter Brune 
294*24989b8cSPeter Brune    Calling sequence for func:
295*24989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
296*24989b8cSPeter Brune .  x - dimensional pointer to state at which to evaluate residual
297*24989b8cSPeter Brune .  f - dimensional pointer to residual, write the residual here
298*24989b8cSPeter Brune -  ctx - optional context passed above
299*24989b8cSPeter Brune 
300*24989b8cSPeter Brune    Level: beginner
301*24989b8cSPeter Brune 
302*24989b8cSPeter Brune .seealso: DMTSSetJacobian(), DMDATSSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
303*24989b8cSPeter Brune @*/
304*24989b8cSPeter Brune PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*),void *ctx)
305*24989b8cSPeter Brune {
306*24989b8cSPeter Brune   PetscErrorCode ierr;
307*24989b8cSPeter Brune   TSDM         sdm;
308*24989b8cSPeter Brune   DM_DA_TS     *dmdats;
309*24989b8cSPeter Brune 
310*24989b8cSPeter Brune   PetscFunctionBegin;
311*24989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
312*24989b8cSPeter Brune   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
313*24989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
314*24989b8cSPeter Brune   dmdats->rhsjacobianlocal = func;
315*24989b8cSPeter Brune   dmdats->rhsjacobianlocalctx = ctx;
316*24989b8cSPeter Brune   ierr = DMTSSetRHSJacobian(dm,TSComputeRHSJacobian_DMDA,dmdats);CHKERRQ(ierr);
317*24989b8cSPeter Brune   PetscFunctionReturn(0);
318*24989b8cSPeter Brune }
319*24989b8cSPeter Brune 
320*24989b8cSPeter Brune 
321*24989b8cSPeter Brune #undef __FUNCT__
322*24989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIFunctionLocal"
323*24989b8cSPeter Brune /*@C
324*24989b8cSPeter Brune    DMDATSSetIFunctionLocal - set a local residual evaluation function
325*24989b8cSPeter Brune 
326*24989b8cSPeter Brune    Logically Collective
327*24989b8cSPeter Brune 
328*24989b8cSPeter Brune    Input Arguments:
329*24989b8cSPeter Brune +  dm - DM to associate callback with
330*24989b8cSPeter Brune .  func - local residual evaluation
331*24989b8cSPeter Brune -  ctx - optional context for local residual evaluation
332*24989b8cSPeter Brune 
333*24989b8cSPeter Brune    Calling sequence for func:
334*24989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
335*24989b8cSPeter Brune .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
336*24989b8cSPeter Brune .  x - dimensional pointer to state at which to evaluate residual
337*24989b8cSPeter Brune .  f - dimensional pointer to residual, write the residual here
338*24989b8cSPeter Brune -  ctx - optional context passed above
339*24989b8cSPeter Brune 
340*24989b8cSPeter Brune    Level: beginner
341*24989b8cSPeter Brune 
342*24989b8cSPeter Brune .seealso: DMTSSetFunction(), DMDATSSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
343*24989b8cSPeter Brune @*/
344*24989b8cSPeter Brune PetscErrorCode DMDATSSetIFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *ctx)
345*24989b8cSPeter Brune {
346*24989b8cSPeter Brune   PetscErrorCode ierr;
347*24989b8cSPeter Brune   TSDM         sdm;
348*24989b8cSPeter Brune   DM_DA_TS     *dmdats;
349*24989b8cSPeter Brune 
350*24989b8cSPeter Brune   PetscFunctionBegin;
351*24989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
352*24989b8cSPeter Brune   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
353*24989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
354*24989b8cSPeter Brune   dmdats->ifunctionlocalimode = imode;
355*24989b8cSPeter Brune   dmdats->ifunctionlocal = func;
356*24989b8cSPeter Brune   dmdats->ifunctionlocalctx = ctx;
357*24989b8cSPeter Brune   ierr = DMTSSetIFunction(dm,TSComputeIFunction_DMDA,dmdats);CHKERRQ(ierr);
358*24989b8cSPeter Brune   PetscFunctionReturn(0);
359*24989b8cSPeter Brune }
360*24989b8cSPeter Brune 
361*24989b8cSPeter Brune #undef __FUNCT__
362*24989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIJacobianLocal"
363*24989b8cSPeter Brune /*@C
364*24989b8cSPeter Brune    DMDATSSetIJacobianLocal - set a local residual evaluation function
365*24989b8cSPeter Brune 
366*24989b8cSPeter Brune    Logically Collective
367*24989b8cSPeter Brune 
368*24989b8cSPeter Brune    Input Arguments:
369*24989b8cSPeter Brune +  dm - DM to associate callback with
370*24989b8cSPeter Brune .  func - local residual evaluation
371*24989b8cSPeter Brune -  ctx - optional context for local residual evaluation
372*24989b8cSPeter Brune 
373*24989b8cSPeter Brune    Calling sequence for func:
374*24989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
375*24989b8cSPeter Brune .  x - dimensional pointer to state at which to evaluate residual
376*24989b8cSPeter Brune .  f - dimensional pointer to residual, write the residual here
377*24989b8cSPeter Brune -  ctx - optional context passed above
378*24989b8cSPeter Brune 
379*24989b8cSPeter Brune    Level: beginner
380*24989b8cSPeter Brune 
381*24989b8cSPeter Brune .seealso: DMTSSetJacobian(), DMDATSSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
382*24989b8cSPeter Brune @*/
383*24989b8cSPeter Brune PetscErrorCode DMDATSSetIJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*),void *ctx)
384*24989b8cSPeter Brune {
385*24989b8cSPeter Brune   PetscErrorCode ierr;
386*24989b8cSPeter Brune   TSDM         sdm;
387*24989b8cSPeter Brune   DM_DA_TS     *dmdats;
388*24989b8cSPeter Brune 
389*24989b8cSPeter Brune   PetscFunctionBegin;
390*24989b8cSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
391*24989b8cSPeter Brune   ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr);
392*24989b8cSPeter Brune   ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr);
393*24989b8cSPeter Brune   dmdats->ijacobianlocal = func;
394*24989b8cSPeter Brune   dmdats->ijacobianlocalctx = ctx;
395*24989b8cSPeter Brune   ierr = DMTSSetIJacobian(dm,TSComputeIJacobian_DMDA,dmdats);CHKERRQ(ierr);
396*24989b8cSPeter Brune   PetscFunctionReturn(0);
397*24989b8cSPeter Brune }
398