xref: /petsc/src/ts/utils/dmdats.c (revision bcf0153e883cfed9568ef4557dcc209048fb58f7)
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 
20d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSDestroy_DMDA(DMTS sdm)
21d71ae5a4SJacob Faibussowitsch {
2224989b8cSPeter Brune   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(PetscFree(sdm->data));
2424989b8cSPeter Brune   PetscFunctionReturn(0);
2524989b8cSPeter Brune }
2624989b8cSPeter Brune 
27d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSDuplicate_DMDA(DMTS oldsdm, DMTS sdm)
28d71ae5a4SJacob Faibussowitsch {
290313d0e0SJed Brown   PetscFunctionBegin;
304dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew((DMTS_DA **)&sdm->data));
319566063dSJacob Faibussowitsch   if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data, oldsdm->data, sizeof(DMTS_DA)));
320313d0e0SJed Brown   PetscFunctionReturn(0);
330313d0e0SJed Brown }
340313d0e0SJed Brown 
35d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDATSGetContext(DM dm, DMTS sdm, DMTS_DA **dmdats)
36d71ae5a4SJacob Faibussowitsch {
3724989b8cSPeter Brune   PetscFunctionBegin;
380298fd71SBarry Smith   *dmdats = NULL;
3924989b8cSPeter Brune   if (!sdm->data) {
404dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMTS_DA **)&sdm->data));
41d74926cbSBarry Smith     sdm->ops->destroy   = DMTSDestroy_DMDA;
42d74926cbSBarry Smith     sdm->ops->duplicate = DMTSDuplicate_DMDA;
4324989b8cSPeter Brune   }
44942e3340SBarry Smith   *dmdats = (DMTS_DA *)sdm->data;
4524989b8cSPeter Brune   PetscFunctionReturn(0);
4624989b8cSPeter Brune }
4724989b8cSPeter Brune 
48d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeIFunction_DMDA(TS ts, PetscReal ptime, Vec X, Vec Xdot, Vec F, void *ctx)
49d71ae5a4SJacob Faibussowitsch {
5024989b8cSPeter Brune   DM            dm;
51942e3340SBarry Smith   DMTS_DA      *dmdats = (DMTS_DA *)ctx;
5224989b8cSPeter Brune   DMDALocalInfo info;
53d2fc88d6Ssepideh   Vec           Xloc, Xdotloc;
5424989b8cSPeter Brune   void         *x, *f, *xdot;
5524989b8cSPeter Brune 
5624989b8cSPeter Brune   PetscFunctionBegin;
5724989b8cSPeter Brune   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
58064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(X, VEC_CLASSID, 3);
59064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(F, VEC_CLASSID, 5);
603c633725SBarry Smith   PetscCheck(dmdats->ifunctionlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Corrupt context");
619566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
629566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xdotloc));
639566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, Xdot, INSERT_VALUES, Xdotloc));
649566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, Xdot, INSERT_VALUES, Xdotloc));
659566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
669566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
679566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
689566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
699566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
709566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xdotloc, &xdot));
7124989b8cSPeter Brune   switch (dmdats->ifunctionlocalimode) {
7224989b8cSPeter Brune   case INSERT_VALUES: {
739566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, F, &f));
7424989b8cSPeter Brune     CHKMEMQ;
759566063dSJacob Faibussowitsch     PetscCall((*dmdats->ifunctionlocal)(&info, ptime, x, xdot, f, dmdats->ifunctionlocalctx));
7624989b8cSPeter Brune     CHKMEMQ;
779566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, F, &f));
7824989b8cSPeter Brune   } break;
7924989b8cSPeter Brune   case ADD_VALUES: {
8024989b8cSPeter Brune     Vec Floc;
819566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Floc));
829566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
839566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Floc, &f));
8424989b8cSPeter Brune     CHKMEMQ;
859566063dSJacob Faibussowitsch     PetscCall((*dmdats->ifunctionlocal)(&info, ptime, x, xdot, f, dmdats->ifunctionlocalctx));
8624989b8cSPeter Brune     CHKMEMQ;
879566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
889566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
899566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
909566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
919566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Floc));
9224989b8cSPeter Brune   } break;
93d71ae5a4SJacob Faibussowitsch   default:
94d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdats->ifunctionlocalimode);
9524989b8cSPeter Brune   }
969566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
979566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
989566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xdotloc, &xdot));
999566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xdotloc));
10024989b8cSPeter Brune   PetscFunctionReturn(0);
10124989b8cSPeter Brune }
10224989b8cSPeter Brune 
103d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeIJacobian_DMDA(TS ts, PetscReal ptime, Vec X, Vec Xdot, PetscReal shift, Mat A, Mat B, void *ctx)
104d71ae5a4SJacob Faibussowitsch {
10524989b8cSPeter Brune   DM            dm;
106942e3340SBarry Smith   DMTS_DA      *dmdats = (DMTS_DA *)ctx;
10724989b8cSPeter Brune   DMDALocalInfo info;
10824989b8cSPeter Brune   Vec           Xloc;
10924989b8cSPeter Brune   void         *x, *xdot;
11024989b8cSPeter Brune 
11124989b8cSPeter Brune   PetscFunctionBegin;
1123c633725SBarry Smith   PetscCheck(dmdats->ifunctionlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Corrupt context");
1139566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
11424989b8cSPeter Brune 
11524989b8cSPeter Brune   if (dmdats->ijacobianlocal) {
1169566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Xloc));
1179566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1189566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1199566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(dm, &info));
1209566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Xloc, &x));
1219566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Xdot, &xdot));
12224989b8cSPeter Brune     CHKMEMQ;
1239566063dSJacob Faibussowitsch     PetscCall((*dmdats->ijacobianlocal)(&info, ptime, x, xdot, shift, A, B, dmdats->ijacobianlocalctx));
12424989b8cSPeter Brune     CHKMEMQ;
1259566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1269566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Xdot, &xdot));
1279566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Xloc));
128ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()");
12924989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
13094ab13aaSBarry Smith   if (A != B) {
1319566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1329566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
13324989b8cSPeter Brune   }
13424989b8cSPeter Brune   PetscFunctionReturn(0);
13524989b8cSPeter Brune }
13624989b8cSPeter Brune 
137d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts, PetscReal ptime, Vec X, Vec F, void *ctx)
138d71ae5a4SJacob Faibussowitsch {
13924989b8cSPeter Brune   DM            dm;
140942e3340SBarry Smith   DMTS_DA      *dmdats = (DMTS_DA *)ctx;
14124989b8cSPeter Brune   DMDALocalInfo info;
14224989b8cSPeter Brune   Vec           Xloc;
14324989b8cSPeter Brune   void         *x, *f;
14424989b8cSPeter Brune 
14524989b8cSPeter Brune   PetscFunctionBegin;
14624989b8cSPeter Brune   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
147064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(X, VEC_CLASSID, 3);
148064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(F, VEC_CLASSID, 4);
1493c633725SBarry Smith   PetscCheck(dmdats->rhsfunctionlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Corrupt context");
1509566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
1519566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
1529566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1539566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1549566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
1559566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
15624989b8cSPeter Brune   switch (dmdats->rhsfunctionlocalimode) {
15724989b8cSPeter Brune   case INSERT_VALUES: {
1589566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, F, &f));
15924989b8cSPeter Brune     CHKMEMQ;
1609566063dSJacob Faibussowitsch     PetscCall((*dmdats->rhsfunctionlocal)(&info, ptime, x, f, dmdats->rhsfunctionlocalctx));
16124989b8cSPeter Brune     CHKMEMQ;
1629566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, F, &f));
16324989b8cSPeter Brune   } break;
16424989b8cSPeter Brune   case ADD_VALUES: {
16524989b8cSPeter Brune     Vec Floc;
1669566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Floc));
1679566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
1689566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Floc, &f));
16924989b8cSPeter Brune     CHKMEMQ;
1709566063dSJacob Faibussowitsch     PetscCall((*dmdats->rhsfunctionlocal)(&info, ptime, x, f, dmdats->rhsfunctionlocalctx));
17124989b8cSPeter Brune     CHKMEMQ;
1729566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
1739566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
1749566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
1759566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
1769566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Floc));
17724989b8cSPeter Brune   } break;
178d71ae5a4SJacob Faibussowitsch   default:
179d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdats->rhsfunctionlocalimode);
18024989b8cSPeter Brune   }
1819566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1829566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
18324989b8cSPeter Brune   PetscFunctionReturn(0);
18424989b8cSPeter Brune }
18524989b8cSPeter Brune 
186d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts, PetscReal ptime, Vec X, Mat A, Mat B, void *ctx)
187d71ae5a4SJacob Faibussowitsch {
18824989b8cSPeter Brune   DM            dm;
189942e3340SBarry Smith   DMTS_DA      *dmdats = (DMTS_DA *)ctx;
19024989b8cSPeter Brune   DMDALocalInfo info;
19124989b8cSPeter Brune   Vec           Xloc;
19224989b8cSPeter Brune   void         *x;
19324989b8cSPeter Brune 
19424989b8cSPeter Brune   PetscFunctionBegin;
1953c633725SBarry Smith   PetscCheck(dmdats->rhsfunctionlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Corrupt context");
1969566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
19724989b8cSPeter Brune 
19824989b8cSPeter Brune   if (dmdats->rhsjacobianlocal) {
1999566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Xloc));
2009566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
2019566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
2029566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(dm, &info));
2039566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Xloc, &x));
20424989b8cSPeter Brune     CHKMEMQ;
2059566063dSJacob Faibussowitsch     PetscCall((*dmdats->rhsjacobianlocal)(&info, ptime, x, A, B, dmdats->rhsjacobianlocalctx));
20624989b8cSPeter Brune     CHKMEMQ;
2079566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
2089566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Xloc));
209ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()");
21024989b8cSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
21194ab13aaSBarry Smith   if (A != B) {
2129566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2139566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
21424989b8cSPeter Brune   }
21524989b8cSPeter Brune   PetscFunctionReturn(0);
21624989b8cSPeter Brune }
21724989b8cSPeter Brune 
21824989b8cSPeter Brune /*@C
219*bcf0153eSBarry Smith    DMDATSSetRHSFunctionLocal - set a local residual evaluation function for use with `DMDA`
22024989b8cSPeter Brune 
22124989b8cSPeter Brune    Logically Collective
22224989b8cSPeter Brune 
2234165533cSJose E. Roman    Input Parameters:
224*bcf0153eSBarry Smith +  dm - `DM` to associate callback with
225b398ffe5SPeter Brune .  imode - insert mode for the residual
22624989b8cSPeter Brune .  func - local residual evaluation
22724989b8cSPeter Brune -  ctx - optional context for local residual evaluation
22824989b8cSPeter Brune 
22924989b8cSPeter Brune    Calling sequence for func:
230b398ffe5SPeter Brune 
231b398ffe5SPeter Brune $ func(DMDALocalInfo info,PetscReal t,void *x,void *f,void *ctx)
232b398ffe5SPeter Brune 
233*bcf0153eSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
234b398ffe5SPeter Brune .  t - time at which to evaluate residual
235b398ffe5SPeter Brune .  x - array of local state information
236b398ffe5SPeter Brune .  f - output array of local residual information
237b398ffe5SPeter Brune -  ctx - optional user context
23824989b8cSPeter Brune 
23924989b8cSPeter Brune    Level: beginner
24024989b8cSPeter Brune 
241*bcf0153eSBarry Smith .seealso: [](chapter_ts), `DMDA`, `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `DMDATSSetRHSJacobianLocal()`, `DMDASNESSetFunctionLocal()`
24224989b8cSPeter Brune @*/
243d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm, InsertMode imode, DMDATSRHSFunctionLocal func, void *ctx)
244d71ae5a4SJacob Faibussowitsch {
245942e3340SBarry Smith   DMTS     sdm;
246942e3340SBarry Smith   DMTS_DA *dmdats;
24724989b8cSPeter Brune 
24824989b8cSPeter Brune   PetscFunctionBegin;
24924989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2509566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &sdm));
2519566063dSJacob Faibussowitsch   PetscCall(DMDATSGetContext(dm, sdm, &dmdats));
25224989b8cSPeter Brune   dmdats->rhsfunctionlocalimode = imode;
25324989b8cSPeter Brune   dmdats->rhsfunctionlocal      = func;
25424989b8cSPeter Brune   dmdats->rhsfunctionlocalctx   = ctx;
2559566063dSJacob Faibussowitsch   PetscCall(DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMDA, dmdats));
25624989b8cSPeter Brune   PetscFunctionReturn(0);
25724989b8cSPeter Brune }
25824989b8cSPeter Brune 
25924989b8cSPeter Brune /*@C
260*bcf0153eSBarry Smith    DMDATSSetRHSJacobianLocal - set a local residual evaluation function for use with `DMDA`
26124989b8cSPeter Brune 
26224989b8cSPeter Brune    Logically Collective
26324989b8cSPeter Brune 
2644165533cSJose E. Roman    Input Parameters:
265*bcf0153eSBarry Smith +  dm    - `DM` to associate callback with
266b398ffe5SPeter Brune .  func  - local RHS Jacobian evaluation routine
267b398ffe5SPeter Brune -  ctx   - optional context for local jacobian evaluation
26824989b8cSPeter Brune 
26924989b8cSPeter Brune    Calling sequence for func:
270b398ffe5SPeter Brune 
2717f6be190SBarry Smith $ func(DMDALocalInfo* info,PetscReal t,void* x,Mat J,Mat B,void *ctx);
272b398ffe5SPeter Brune 
273*bcf0153eSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
274b398ffe5SPeter Brune .  t    - time at which to evaluate residual
275b398ffe5SPeter Brune .  x    - array of local state information
276b398ffe5SPeter Brune .  J    - Jacobian matrix
277b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
27824989b8cSPeter Brune -  ctx  - optional context passed above
27924989b8cSPeter Brune 
28024989b8cSPeter Brune    Level: beginner
28124989b8cSPeter Brune 
282*bcf0153eSBarry Smith .seealso: [](chapter_ts), `DMDA`, `DMTSSetRHSJacobian()`, `DMDATSSetRHSFunctionLocal()`, `DMDASNESSetJacobianLocal()`
28324989b8cSPeter Brune @*/
284d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm, DMDATSRHSJacobianLocal func, void *ctx)
285d71ae5a4SJacob Faibussowitsch {
286942e3340SBarry Smith   DMTS     sdm;
287942e3340SBarry Smith   DMTS_DA *dmdats;
28824989b8cSPeter Brune 
28924989b8cSPeter Brune   PetscFunctionBegin;
29024989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2919566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &sdm));
2929566063dSJacob Faibussowitsch   PetscCall(DMDATSGetContext(dm, sdm, &dmdats));
29324989b8cSPeter Brune   dmdats->rhsjacobianlocal    = func;
29424989b8cSPeter Brune   dmdats->rhsjacobianlocalctx = ctx;
2959566063dSJacob Faibussowitsch   PetscCall(DMTSSetRHSJacobian(dm, TSComputeRHSJacobian_DMDA, dmdats));
29624989b8cSPeter Brune   PetscFunctionReturn(0);
29724989b8cSPeter Brune }
29824989b8cSPeter Brune 
29924989b8cSPeter Brune /*@C
300*bcf0153eSBarry Smith    DMDATSSetIFunctionLocal - set a local residual evaluation function for use with `DMDA`
30124989b8cSPeter Brune 
30224989b8cSPeter Brune    Logically Collective
30324989b8cSPeter Brune 
3044165533cSJose E. Roman    Input Parameters:
305*bcf0153eSBarry Smith +  dm   - `DM` to associate callback with
30624989b8cSPeter Brune .  func - local residual evaluation
30724989b8cSPeter Brune -  ctx  - optional context for local residual evaluation
30824989b8cSPeter Brune 
30924989b8cSPeter Brune    Calling sequence for func:
310*bcf0153eSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
311b398ffe5SPeter Brune .  t    - time at which to evaluate residual
312b398ffe5SPeter Brune .  x    - array of local state information
313b398ffe5SPeter Brune .  xdot - array of local time derivative information
314b398ffe5SPeter Brune .  f    - output array of local function evaluation information
31524989b8cSPeter Brune -  ctx - optional context passed above
31624989b8cSPeter Brune 
31724989b8cSPeter Brune    Level: beginner
31824989b8cSPeter Brune 
319*bcf0153eSBarry Smith .seealso: [](chapter_ts), `DMDA`, `DMTSSetIFunction()`, `DMDATSSetIJacobianLocal()`, `DMDASNESSetFunctionLocal()`
32024989b8cSPeter Brune @*/
321d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDATSSetIFunctionLocal(DM dm, InsertMode imode, DMDATSIFunctionLocal func, void *ctx)
322d71ae5a4SJacob Faibussowitsch {
323942e3340SBarry Smith   DMTS     sdm;
324942e3340SBarry Smith   DMTS_DA *dmdats;
32524989b8cSPeter Brune 
32624989b8cSPeter Brune   PetscFunctionBegin;
32724989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3289566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &sdm));
3299566063dSJacob Faibussowitsch   PetscCall(DMDATSGetContext(dm, sdm, &dmdats));
33024989b8cSPeter Brune   dmdats->ifunctionlocalimode = imode;
33124989b8cSPeter Brune   dmdats->ifunctionlocal      = func;
33224989b8cSPeter Brune   dmdats->ifunctionlocalctx   = ctx;
3339566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMDA, dmdats));
33424989b8cSPeter Brune   PetscFunctionReturn(0);
33524989b8cSPeter Brune }
33624989b8cSPeter Brune 
33724989b8cSPeter Brune /*@C
338*bcf0153eSBarry Smith    DMDATSSetIJacobianLocal - set a local residual evaluation function for use with `DMDA`
33924989b8cSPeter Brune 
34024989b8cSPeter Brune    Logically Collective
34124989b8cSPeter Brune 
3424165533cSJose E. Roman    Input Parameters:
343*bcf0153eSBarry Smith +  dm   - `DM` to associate callback with
34424989b8cSPeter Brune .  func - local residual evaluation
34524989b8cSPeter Brune -  ctx   - optional context for local residual evaluation
34624989b8cSPeter Brune 
34724989b8cSPeter Brune    Calling sequence for func:
348b398ffe5SPeter Brune 
349dab4423fSBarry Smith $ func(DMDALocalInfo* info,PetscReal t,void* x,void *xdot,PetscScalar shift,Mat J,Mat B,void *ctx);
350b398ffe5SPeter Brune 
351*bcf0153eSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
352b398ffe5SPeter Brune .  t    - time at which to evaluate the jacobian
353b398ffe5SPeter Brune .  x    - array of local state information
354b398ffe5SPeter Brune .  xdot - time derivative at this state
355dab4423fSBarry Smith .  shift - see TSSetIJacobian() for the meaning of this parameter
356b398ffe5SPeter Brune .  J    - Jacobian matrix
357b398ffe5SPeter Brune .  B    - preconditioner matrix; often same as J
35824989b8cSPeter Brune -  ctx  - optional context passed above
35924989b8cSPeter Brune 
36024989b8cSPeter Brune    Level: beginner
36124989b8cSPeter Brune 
362*bcf0153eSBarry Smith .seealso: [](chapter_ts), `DMDA`, `DMTSSetJacobian()`, `DMDATSSetIFunctionLocal()`, `DMDASNESSetJacobianLocal()`
36324989b8cSPeter Brune @*/
364d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDATSSetIJacobianLocal(DM dm, DMDATSIJacobianLocal func, void *ctx)
365d71ae5a4SJacob Faibussowitsch {
366942e3340SBarry Smith   DMTS     sdm;
367942e3340SBarry Smith   DMTS_DA *dmdats;
36824989b8cSPeter Brune 
36924989b8cSPeter Brune   PetscFunctionBegin;
37024989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3719566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &sdm));
3729566063dSJacob Faibussowitsch   PetscCall(DMDATSGetContext(dm, sdm, &dmdats));
37324989b8cSPeter Brune   dmdats->ijacobianlocal    = func;
37424989b8cSPeter Brune   dmdats->ijacobianlocalctx = ctx;
3759566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMDA, dmdats));
37624989b8cSPeter Brune   PetscFunctionReturn(0);
37724989b8cSPeter Brune }
378d1212d36SBarry Smith 
379d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDMDARayDestroy(void **mctx)
380d71ae5a4SJacob Faibussowitsch {
381d1212d36SBarry Smith   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *)*mctx;
382d1212d36SBarry Smith 
383d1212d36SBarry Smith   PetscFunctionBegin;
3849566063dSJacob Faibussowitsch   if (rayctx->lgctx) PetscCall(TSMonitorLGCtxDestroy(&rayctx->lgctx));
3859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rayctx->ray));
3869566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&rayctx->scatter));
3879566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&rayctx->viewer));
3889566063dSJacob Faibussowitsch   PetscCall(PetscFree(rayctx));
389d1212d36SBarry Smith   PetscFunctionReturn(0);
390d1212d36SBarry Smith }
391d1212d36SBarry Smith 
392d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDMDARay(TS ts, PetscInt steps, PetscReal time, Vec u, void *mctx)
393d71ae5a4SJacob Faibussowitsch {
394d1212d36SBarry Smith   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *)mctx;
395d1212d36SBarry Smith   Vec                  solution;
396d1212d36SBarry Smith 
397d1212d36SBarry Smith   PetscFunctionBegin;
3989566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &solution));
3999566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(rayctx->scatter, solution, rayctx->ray, INSERT_VALUES, SCATTER_FORWARD));
4009566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(rayctx->scatter, solution, rayctx->ray, INSERT_VALUES, SCATTER_FORWARD));
4011baa6e33SBarry Smith   if (rayctx->viewer) PetscCall(VecView(rayctx->ray, rayctx->viewer));
402d1212d36SBarry Smith   PetscFunctionReturn(0);
403d1212d36SBarry Smith }
40451b4a12fSMatthew G. Knepley 
405d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGDMDARay(TS ts, PetscInt step, PetscReal ptime, Vec u, void *ctx)
406d71ae5a4SJacob Faibussowitsch {
40751b4a12fSMatthew G. Knepley   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *)ctx;
40851b4a12fSMatthew G. Knepley   TSMonitorLGCtx       lgctx  = (TSMonitorLGCtx)rayctx->lgctx;
40951b4a12fSMatthew G. Knepley   Vec                  v      = rayctx->ray;
41051b4a12fSMatthew G. Knepley   const PetscScalar   *a;
41151b4a12fSMatthew G. Knepley   PetscInt             dim;
41251b4a12fSMatthew G. Knepley 
41351b4a12fSMatthew G. Knepley   PetscFunctionBegin;
4149566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD));
4159566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD));
41651b4a12fSMatthew G. Knepley   if (!step) {
41751b4a12fSMatthew G. Knepley     PetscDrawAxis axis;
41851b4a12fSMatthew G. Knepley 
4199566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(lgctx->lg, &axis));
4209566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Solution Ray as function of time", "Time", "Solution"));
4219566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(rayctx->ray, &dim));
4229566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetDimension(lgctx->lg, dim));
4239566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(lgctx->lg));
42451b4a12fSMatthew G. Knepley   }
4259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(v, &a));
42651b4a12fSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX)
42751b4a12fSMatthew G. Knepley   {
42851b4a12fSMatthew G. Knepley     PetscReal *areal;
42951b4a12fSMatthew G. Knepley     PetscInt   i, n;
4309566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v, &n));
4319566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &areal));
43251b4a12fSMatthew G. Knepley     for (i = 0; i < n; ++i) areal[i] = PetscRealPart(a[i]);
4339566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(lgctx->lg, ptime, areal));
4349566063dSJacob Faibussowitsch     PetscCall(PetscFree(areal));
43551b4a12fSMatthew G. Knepley   }
43651b4a12fSMatthew G. Knepley #else
4379566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddCommonPoint(lgctx->lg, ptime, a));
43851b4a12fSMatthew G. Knepley #endif
4399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(v, &a));
44051b4a12fSMatthew G. Knepley   if (((lgctx->howoften > 0) && (!(step % lgctx->howoften))) || ((lgctx->howoften == -1) && ts->reason)) {
4419566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(lgctx->lg));
4429566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(lgctx->lg));
44351b4a12fSMatthew G. Knepley   }
44451b4a12fSMatthew G. Knepley   PetscFunctionReturn(0);
44551b4a12fSMatthew G. Knepley }
446