xref: /petsc/src/ts/utils/dmdats.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
124989b8cSPeter Brune #include <petscdmda.h> /*I "petscdmda.h" I*/
2af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
3af0996ceSBarry Smith #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 
209371c9d4SSatish Balay static PetscErrorCode DMTSDestroy_DMDA(DMTS sdm) {
2124989b8cSPeter Brune   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   PetscCall(PetscFree(sdm->data));
2324989b8cSPeter Brune   PetscFunctionReturn(0);
2424989b8cSPeter Brune }
2524989b8cSPeter Brune 
269371c9d4SSatish Balay static PetscErrorCode DMTSDuplicate_DMDA(DMTS oldsdm, DMTS sdm) {
270313d0e0SJed Brown   PetscFunctionBegin;
28*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew((DMTS_DA **)&sdm->data));
299566063dSJacob Faibussowitsch   if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data, oldsdm->data, sizeof(DMTS_DA)));
300313d0e0SJed Brown   PetscFunctionReturn(0);
310313d0e0SJed Brown }
320313d0e0SJed Brown 
339371c9d4SSatish Balay static PetscErrorCode DMDATSGetContext(DM dm, DMTS sdm, DMTS_DA **dmdats) {
3424989b8cSPeter Brune   PetscFunctionBegin;
350298fd71SBarry Smith   *dmdats = NULL;
3624989b8cSPeter Brune   if (!sdm->data) {
37*4dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMTS_DA **)&sdm->data));
38d74926cbSBarry Smith     sdm->ops->destroy   = DMTSDestroy_DMDA;
39d74926cbSBarry Smith     sdm->ops->duplicate = DMTSDuplicate_DMDA;
4024989b8cSPeter Brune   }
41942e3340SBarry Smith   *dmdats = (DMTS_DA *)sdm->data;
4224989b8cSPeter Brune   PetscFunctionReturn(0);
4324989b8cSPeter Brune }
4424989b8cSPeter Brune 
459371c9d4SSatish Balay static PetscErrorCode TSComputeIFunction_DMDA(TS ts, PetscReal ptime, Vec X, Vec Xdot, Vec F, void *ctx) {
4624989b8cSPeter Brune   DM            dm;
47942e3340SBarry Smith   DMTS_DA      *dmdats = (DMTS_DA *)ctx;
4824989b8cSPeter Brune   DMDALocalInfo info;
49d2fc88d6Ssepideh   Vec           Xloc, Xdotloc;
5024989b8cSPeter Brune   void         *x, *f, *xdot;
5124989b8cSPeter Brune 
5224989b8cSPeter Brune   PetscFunctionBegin;
5324989b8cSPeter Brune   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
54064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(X, VEC_CLASSID, 3);
55064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(F, VEC_CLASSID, 5);
563c633725SBarry Smith   PetscCheck(dmdats->ifunctionlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Corrupt context");
579566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
589566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xdotloc));
599566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, Xdot, INSERT_VALUES, Xdotloc));
609566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, Xdot, INSERT_VALUES, Xdotloc));
619566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
629566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
639566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
649566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
659566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
669566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xdotloc, &xdot));
6724989b8cSPeter Brune   switch (dmdats->ifunctionlocalimode) {
6824989b8cSPeter Brune   case INSERT_VALUES: {
699566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, F, &f));
7024989b8cSPeter Brune     CHKMEMQ;
719566063dSJacob Faibussowitsch     PetscCall((*dmdats->ifunctionlocal)(&info, ptime, x, xdot, f, dmdats->ifunctionlocalctx));
7224989b8cSPeter Brune     CHKMEMQ;
739566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, F, &f));
7424989b8cSPeter Brune   } break;
7524989b8cSPeter Brune   case ADD_VALUES: {
7624989b8cSPeter Brune     Vec Floc;
779566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Floc));
789566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
799566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Floc, &f));
8024989b8cSPeter Brune     CHKMEMQ;
819566063dSJacob Faibussowitsch     PetscCall((*dmdats->ifunctionlocal)(&info, ptime, x, xdot, f, dmdats->ifunctionlocalctx));
8224989b8cSPeter Brune     CHKMEMQ;
839566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
849566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
859566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
869566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
879566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Floc));
8824989b8cSPeter Brune   } break;
8998921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdats->ifunctionlocalimode);
9024989b8cSPeter Brune   }
919566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
929566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
939566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xdotloc, &xdot));
949566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xdotloc));
9524989b8cSPeter Brune   PetscFunctionReturn(0);
9624989b8cSPeter Brune }
9724989b8cSPeter Brune 
989371c9d4SSatish Balay static PetscErrorCode TSComputeIJacobian_DMDA(TS ts, PetscReal ptime, Vec X, Vec Xdot, PetscReal shift, Mat A, Mat B, void *ctx) {
9924989b8cSPeter Brune   DM            dm;
100942e3340SBarry Smith   DMTS_DA      *dmdats = (DMTS_DA *)ctx;
10124989b8cSPeter Brune   DMDALocalInfo info;
10224989b8cSPeter Brune   Vec           Xloc;
10324989b8cSPeter Brune   void         *x, *xdot;
10424989b8cSPeter Brune 
10524989b8cSPeter Brune   PetscFunctionBegin;
1063c633725SBarry Smith   PetscCheck(dmdats->ifunctionlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Corrupt context");
1079566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
10824989b8cSPeter Brune 
10924989b8cSPeter Brune   if (dmdats->ijacobianlocal) {
1109566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Xloc));
1119566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1129566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1139566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(dm, &info));
1149566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Xloc, &x));
1159566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Xdot, &xdot));
11624989b8cSPeter Brune     CHKMEMQ;
1179566063dSJacob Faibussowitsch     PetscCall((*dmdats->ijacobianlocal)(&info, ptime, x, xdot, shift, A, B, dmdats->ijacobianlocalctx));
11824989b8cSPeter Brune     CHKMEMQ;
1199566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1209566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Xdot, &xdot));
1219566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Xloc));
122ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()");
12324989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
12494ab13aaSBarry Smith   if (A != B) {
1259566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1269566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
12724989b8cSPeter Brune   }
12824989b8cSPeter Brune   PetscFunctionReturn(0);
12924989b8cSPeter Brune }
13024989b8cSPeter Brune 
1319371c9d4SSatish Balay static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts, PetscReal ptime, Vec X, Vec F, void *ctx) {
13224989b8cSPeter Brune   DM            dm;
133942e3340SBarry Smith   DMTS_DA      *dmdats = (DMTS_DA *)ctx;
13424989b8cSPeter Brune   DMDALocalInfo info;
13524989b8cSPeter Brune   Vec           Xloc;
13624989b8cSPeter Brune   void         *x, *f;
13724989b8cSPeter Brune 
13824989b8cSPeter Brune   PetscFunctionBegin;
13924989b8cSPeter Brune   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
140064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(X, VEC_CLASSID, 3);
141064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(F, VEC_CLASSID, 4);
1423c633725SBarry Smith   PetscCheck(dmdats->rhsfunctionlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Corrupt context");
1439566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
1449566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
1459566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1469566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1479566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
1489566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
14924989b8cSPeter Brune   switch (dmdats->rhsfunctionlocalimode) {
15024989b8cSPeter Brune   case INSERT_VALUES: {
1519566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, F, &f));
15224989b8cSPeter Brune     CHKMEMQ;
1539566063dSJacob Faibussowitsch     PetscCall((*dmdats->rhsfunctionlocal)(&info, ptime, x, f, dmdats->rhsfunctionlocalctx));
15424989b8cSPeter Brune     CHKMEMQ;
1559566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, F, &f));
15624989b8cSPeter Brune   } break;
15724989b8cSPeter Brune   case ADD_VALUES: {
15824989b8cSPeter Brune     Vec Floc;
1599566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Floc));
1609566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
1619566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Floc, &f));
16224989b8cSPeter Brune     CHKMEMQ;
1639566063dSJacob Faibussowitsch     PetscCall((*dmdats->rhsfunctionlocal)(&info, ptime, x, f, dmdats->rhsfunctionlocalctx));
16424989b8cSPeter Brune     CHKMEMQ;
1659566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
1669566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
1679566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
1689566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
1699566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Floc));
17024989b8cSPeter Brune   } break;
17198921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdats->rhsfunctionlocalimode);
17224989b8cSPeter Brune   }
1739566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1749566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
17524989b8cSPeter Brune   PetscFunctionReturn(0);
17624989b8cSPeter Brune }
17724989b8cSPeter Brune 
1789371c9d4SSatish Balay static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts, PetscReal ptime, Vec X, Mat A, Mat B, void *ctx) {
17924989b8cSPeter Brune   DM            dm;
180942e3340SBarry Smith   DMTS_DA      *dmdats = (DMTS_DA *)ctx;
18124989b8cSPeter Brune   DMDALocalInfo info;
18224989b8cSPeter Brune   Vec           Xloc;
18324989b8cSPeter Brune   void         *x;
18424989b8cSPeter Brune 
18524989b8cSPeter Brune   PetscFunctionBegin;
1863c633725SBarry Smith   PetscCheck(dmdats->rhsfunctionlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Corrupt context");
1879566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
18824989b8cSPeter Brune 
18924989b8cSPeter Brune   if (dmdats->rhsjacobianlocal) {
1909566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Xloc));
1919566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1929566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1939566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(dm, &info));
1949566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Xloc, &x));
19524989b8cSPeter Brune     CHKMEMQ;
1969566063dSJacob Faibussowitsch     PetscCall((*dmdats->rhsjacobianlocal)(&info, ptime, x, A, B, dmdats->rhsjacobianlocalctx));
19724989b8cSPeter Brune     CHKMEMQ;
1989566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1999566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Xloc));
200ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()");
20124989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
20294ab13aaSBarry Smith   if (A != B) {
2039566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2049566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
20524989b8cSPeter Brune   }
20624989b8cSPeter Brune   PetscFunctionReturn(0);
20724989b8cSPeter Brune }
20824989b8cSPeter Brune 
20924989b8cSPeter Brune /*@C
21024989b8cSPeter Brune    DMDATSSetRHSFunctionLocal - set a local residual evaluation function
21124989b8cSPeter Brune 
21224989b8cSPeter Brune    Logically Collective
21324989b8cSPeter Brune 
2144165533cSJose E. Roman    Input Parameters:
21524989b8cSPeter Brune +  dm - DM to associate callback with
216b398ffe5SPeter Brune .  imode - insert mode for the residual
21724989b8cSPeter Brune .  func - local residual evaluation
21824989b8cSPeter Brune -  ctx - optional context for local residual evaluation
21924989b8cSPeter Brune 
22024989b8cSPeter Brune    Calling sequence for func:
221b398ffe5SPeter Brune 
222b398ffe5SPeter Brune $ func(DMDALocalInfo info,PetscReal t,void *x,void *f,void *ctx)
223b398ffe5SPeter Brune 
22424989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
225b398ffe5SPeter Brune .  t - time at which to evaluate residual
226b398ffe5SPeter Brune .  x - array of local state information
227b398ffe5SPeter Brune .  f - output array of local residual information
228b398ffe5SPeter Brune -  ctx - optional user context
22924989b8cSPeter Brune 
23024989b8cSPeter Brune    Level: beginner
23124989b8cSPeter Brune 
232db781477SPatrick Sanan .seealso: `DMTSSetRHSFunction()`, `DMDATSSetRHSJacobianLocal()`, `DMDASNESSetFunctionLocal()`
23324989b8cSPeter Brune @*/
2349371c9d4SSatish Balay PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm, InsertMode imode, DMDATSRHSFunctionLocal func, void *ctx) {
235942e3340SBarry Smith   DMTS     sdm;
236942e3340SBarry Smith   DMTS_DA *dmdats;
23724989b8cSPeter Brune 
23824989b8cSPeter Brune   PetscFunctionBegin;
23924989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2409566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &sdm));
2419566063dSJacob Faibussowitsch   PetscCall(DMDATSGetContext(dm, sdm, &dmdats));
24224989b8cSPeter Brune   dmdats->rhsfunctionlocalimode = imode;
24324989b8cSPeter Brune   dmdats->rhsfunctionlocal      = func;
24424989b8cSPeter Brune   dmdats->rhsfunctionlocalctx   = ctx;
2459566063dSJacob Faibussowitsch   PetscCall(DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMDA, dmdats));
24624989b8cSPeter Brune   PetscFunctionReturn(0);
24724989b8cSPeter Brune }
24824989b8cSPeter Brune 
24924989b8cSPeter Brune /*@C
25024989b8cSPeter Brune    DMDATSSetRHSJacobianLocal - set a local residual evaluation function
25124989b8cSPeter Brune 
25224989b8cSPeter Brune    Logically Collective
25324989b8cSPeter Brune 
2544165533cSJose E. Roman    Input Parameters:
25524989b8cSPeter Brune +  dm    - DM to associate callback with
256b398ffe5SPeter Brune .  func  - local RHS Jacobian evaluation routine
257b398ffe5SPeter Brune -  ctx   - optional context for local jacobian evaluation
25824989b8cSPeter Brune 
25924989b8cSPeter Brune    Calling sequence for func:
260b398ffe5SPeter Brune 
2617f6be190SBarry Smith $ func(DMDALocalInfo* info,PetscReal t,void* x,Mat J,Mat B,void *ctx);
262b398ffe5SPeter Brune 
26324989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
264b398ffe5SPeter Brune .  t    - time at which to evaluate residual
265b398ffe5SPeter Brune .  x    - array of local state information
266b398ffe5SPeter Brune .  J    - Jacobian matrix
267b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
26824989b8cSPeter Brune -  ctx  - optional context passed above
26924989b8cSPeter Brune 
27024989b8cSPeter Brune    Level: beginner
27124989b8cSPeter Brune 
272db781477SPatrick Sanan .seealso: `DMTSSetRHSJacobian()`, `DMDATSSetRHSFunctionLocal()`, `DMDASNESSetJacobianLocal()`
27324989b8cSPeter Brune @*/
2749371c9d4SSatish Balay PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm, DMDATSRHSJacobianLocal func, void *ctx) {
275942e3340SBarry Smith   DMTS     sdm;
276942e3340SBarry Smith   DMTS_DA *dmdats;
27724989b8cSPeter Brune 
27824989b8cSPeter Brune   PetscFunctionBegin;
27924989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2809566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &sdm));
2819566063dSJacob Faibussowitsch   PetscCall(DMDATSGetContext(dm, sdm, &dmdats));
28224989b8cSPeter Brune   dmdats->rhsjacobianlocal    = func;
28324989b8cSPeter Brune   dmdats->rhsjacobianlocalctx = ctx;
2849566063dSJacob Faibussowitsch   PetscCall(DMTSSetRHSJacobian(dm, TSComputeRHSJacobian_DMDA, dmdats));
28524989b8cSPeter Brune   PetscFunctionReturn(0);
28624989b8cSPeter Brune }
28724989b8cSPeter Brune 
28824989b8cSPeter Brune /*@C
28924989b8cSPeter Brune    DMDATSSetIFunctionLocal - set a local residual evaluation function
29024989b8cSPeter Brune 
29124989b8cSPeter Brune    Logically Collective
29224989b8cSPeter Brune 
2934165533cSJose E. Roman    Input Parameters:
29424989b8cSPeter Brune +  dm   - DM to associate callback with
29524989b8cSPeter Brune .  func - local residual evaluation
29624989b8cSPeter Brune -  ctx  - optional context for local residual evaluation
29724989b8cSPeter Brune 
29824989b8cSPeter Brune    Calling sequence for func:
29924989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
300b398ffe5SPeter Brune .  t    - time at which to evaluate residual
301b398ffe5SPeter Brune .  x    - array of local state information
302b398ffe5SPeter Brune .  xdot - array of local time derivative information
303b398ffe5SPeter Brune .  f    - output array of local function evaluation information
30424989b8cSPeter Brune -  ctx - optional context passed above
30524989b8cSPeter Brune 
30624989b8cSPeter Brune    Level: beginner
30724989b8cSPeter Brune 
308db781477SPatrick Sanan .seealso: `DMTSSetIFunction()`, `DMDATSSetIJacobianLocal()`, `DMDASNESSetFunctionLocal()`
30924989b8cSPeter Brune @*/
3109371c9d4SSatish Balay PetscErrorCode DMDATSSetIFunctionLocal(DM dm, InsertMode imode, DMDATSIFunctionLocal func, void *ctx) {
311942e3340SBarry Smith   DMTS     sdm;
312942e3340SBarry Smith   DMTS_DA *dmdats;
31324989b8cSPeter Brune 
31424989b8cSPeter Brune   PetscFunctionBegin;
31524989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3169566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &sdm));
3179566063dSJacob Faibussowitsch   PetscCall(DMDATSGetContext(dm, sdm, &dmdats));
31824989b8cSPeter Brune   dmdats->ifunctionlocalimode = imode;
31924989b8cSPeter Brune   dmdats->ifunctionlocal      = func;
32024989b8cSPeter Brune   dmdats->ifunctionlocalctx   = ctx;
3219566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMDA, dmdats));
32224989b8cSPeter Brune   PetscFunctionReturn(0);
32324989b8cSPeter Brune }
32424989b8cSPeter Brune 
32524989b8cSPeter Brune /*@C
32624989b8cSPeter Brune    DMDATSSetIJacobianLocal - set a local residual evaluation function
32724989b8cSPeter Brune 
32824989b8cSPeter Brune    Logically Collective
32924989b8cSPeter Brune 
3304165533cSJose E. Roman    Input Parameters:
33124989b8cSPeter Brune +  dm   - DM to associate callback with
33224989b8cSPeter Brune .  func - local residual evaluation
33324989b8cSPeter Brune -  ctx   - optional context for local residual evaluation
33424989b8cSPeter Brune 
33524989b8cSPeter Brune    Calling sequence for func:
336b398ffe5SPeter Brune 
337dab4423fSBarry Smith $ func(DMDALocalInfo* info,PetscReal t,void* x,void *xdot,PetscScalar shift,Mat J,Mat B,void *ctx);
338b398ffe5SPeter Brune 
33924989b8cSPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
340b398ffe5SPeter Brune .  t    - time at which to evaluate the jacobian
341b398ffe5SPeter Brune .  x    - array of local state information
342b398ffe5SPeter Brune .  xdot - time derivative at this state
343dab4423fSBarry Smith .  shift - see TSSetIJacobian() for the meaning of this parameter
344b398ffe5SPeter Brune .  J    - Jacobian matrix
345b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
34624989b8cSPeter Brune -  ctx  - optional context passed above
34724989b8cSPeter Brune 
34824989b8cSPeter Brune    Level: beginner
34924989b8cSPeter Brune 
350db781477SPatrick Sanan .seealso: `DMTSSetJacobian()`, `DMDATSSetIFunctionLocal()`, `DMDASNESSetJacobianLocal()`
35124989b8cSPeter Brune @*/
3529371c9d4SSatish Balay PetscErrorCode DMDATSSetIJacobianLocal(DM dm, DMDATSIJacobianLocal func, void *ctx) {
353942e3340SBarry Smith   DMTS     sdm;
354942e3340SBarry Smith   DMTS_DA *dmdats;
35524989b8cSPeter Brune 
35624989b8cSPeter Brune   PetscFunctionBegin;
35724989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3589566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &sdm));
3599566063dSJacob Faibussowitsch   PetscCall(DMDATSGetContext(dm, sdm, &dmdats));
36024989b8cSPeter Brune   dmdats->ijacobianlocal    = func;
36124989b8cSPeter Brune   dmdats->ijacobianlocalctx = ctx;
3629566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMDA, dmdats));
36324989b8cSPeter Brune   PetscFunctionReturn(0);
36424989b8cSPeter Brune }
365d1212d36SBarry Smith 
3669371c9d4SSatish Balay PetscErrorCode TSMonitorDMDARayDestroy(void **mctx) {
367d1212d36SBarry Smith   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *)*mctx;
368d1212d36SBarry Smith 
369d1212d36SBarry Smith   PetscFunctionBegin;
3709566063dSJacob Faibussowitsch   if (rayctx->lgctx) PetscCall(TSMonitorLGCtxDestroy(&rayctx->lgctx));
3719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rayctx->ray));
3729566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&rayctx->scatter));
3739566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&rayctx->viewer));
3749566063dSJacob Faibussowitsch   PetscCall(PetscFree(rayctx));
375d1212d36SBarry Smith   PetscFunctionReturn(0);
376d1212d36SBarry Smith }
377d1212d36SBarry Smith 
3789371c9d4SSatish Balay PetscErrorCode TSMonitorDMDARay(TS ts, PetscInt steps, PetscReal time, Vec u, void *mctx) {
379d1212d36SBarry Smith   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *)mctx;
380d1212d36SBarry Smith   Vec                  solution;
381d1212d36SBarry Smith 
382d1212d36SBarry Smith   PetscFunctionBegin;
3839566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &solution));
3849566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(rayctx->scatter, solution, rayctx->ray, INSERT_VALUES, SCATTER_FORWARD));
3859566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(rayctx->scatter, solution, rayctx->ray, INSERT_VALUES, SCATTER_FORWARD));
3861baa6e33SBarry Smith   if (rayctx->viewer) PetscCall(VecView(rayctx->ray, rayctx->viewer));
387d1212d36SBarry Smith   PetscFunctionReturn(0);
388d1212d36SBarry Smith }
38951b4a12fSMatthew G. Knepley 
3909371c9d4SSatish Balay PetscErrorCode TSMonitorLGDMDARay(TS ts, PetscInt step, PetscReal ptime, Vec u, void *ctx) {
39151b4a12fSMatthew G. Knepley   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *)ctx;
39251b4a12fSMatthew G. Knepley   TSMonitorLGCtx       lgctx  = (TSMonitorLGCtx)rayctx->lgctx;
39351b4a12fSMatthew G. Knepley   Vec                  v      = rayctx->ray;
39451b4a12fSMatthew G. Knepley   const PetscScalar   *a;
39551b4a12fSMatthew G. Knepley   PetscInt             dim;
39651b4a12fSMatthew G. Knepley 
39751b4a12fSMatthew G. Knepley   PetscFunctionBegin;
3989566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD));
3999566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD));
40051b4a12fSMatthew G. Knepley   if (!step) {
40151b4a12fSMatthew G. Knepley     PetscDrawAxis axis;
40251b4a12fSMatthew G. Knepley 
4039566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(lgctx->lg, &axis));
4049566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Solution Ray as function of time", "Time", "Solution"));
4059566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(rayctx->ray, &dim));
4069566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetDimension(lgctx->lg, dim));
4079566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(lgctx->lg));
40851b4a12fSMatthew G. Knepley   }
4099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(v, &a));
41051b4a12fSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX)
41151b4a12fSMatthew G. Knepley   {
41251b4a12fSMatthew G. Knepley     PetscReal *areal;
41351b4a12fSMatthew G. Knepley     PetscInt   i, n;
4149566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v, &n));
4159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &areal));
41651b4a12fSMatthew G. Knepley     for (i = 0; i < n; ++i) areal[i] = PetscRealPart(a[i]);
4179566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(lgctx->lg, ptime, areal));
4189566063dSJacob Faibussowitsch     PetscCall(PetscFree(areal));
41951b4a12fSMatthew G. Knepley   }
42051b4a12fSMatthew G. Knepley #else
4219566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddCommonPoint(lgctx->lg, ptime, a));
42251b4a12fSMatthew G. Knepley #endif
4239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(v, &a));
42451b4a12fSMatthew G. Knepley   if (((lgctx->howoften > 0) && (!(step % lgctx->howoften))) || ((lgctx->howoften == -1) && ts->reason)) {
4259566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(lgctx->lg));
4269566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(lgctx->lg));
42751b4a12fSMatthew G. Knepley   }
42851b4a12fSMatthew G. Knepley   PetscFunctionReturn(0);
42951b4a12fSMatthew G. Knepley }
430