xref: /petsc/src/ts/utils/dmdats.c (revision 8434afd195968570cfdb5bc7b9cfc0a316d974ae)
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));
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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)));
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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;
453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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));
1003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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   }
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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));
1833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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   }
2153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21624989b8cSPeter Brune }
21724989b8cSPeter Brune 
21824989b8cSPeter Brune /*@C
219bcf0153eSBarry 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:
224bcf0153eSBarry Smith + dm    - `DM` to associate callback with
225b398ffe5SPeter Brune . imode - insert mode for the residual
226*8434afd1SBarry Smith . func  - local residual evaluation, see `DMDATSRHSFunctionLocalFn` for the calling sequence
22724989b8cSPeter Brune - ctx   - optional context for local residual evaluation
22824989b8cSPeter Brune 
22924989b8cSPeter Brune   Level: beginner
23024989b8cSPeter Brune 
231*8434afd1SBarry Smith .seealso: [](ch_ts), `DMDA`, `DMDATSRHSFunctionLocalFn`, `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `DMDATSSetRHSJacobianLocal()`, `DMDASNESSetFunctionLocal()`
23224989b8cSPeter Brune @*/
233*8434afd1SBarry Smith PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm, InsertMode imode, DMDATSRHSFunctionLocalFn *func, void *ctx)
234d71ae5a4SJacob Faibussowitsch {
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));
2463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24724989b8cSPeter Brune }
24824989b8cSPeter Brune 
24924989b8cSPeter Brune /*@C
250bcf0153eSBarry Smith   DMDATSSetRHSJacobianLocal - set a local residual evaluation function for use with `DMDA`
25124989b8cSPeter Brune 
25224989b8cSPeter Brune   Logically Collective
25324989b8cSPeter Brune 
2544165533cSJose E. Roman   Input Parameters:
255bcf0153eSBarry Smith + dm   - `DM` to associate callback with
256*8434afd1SBarry Smith . func - local RHS Jacobian evaluation routine, see `DMDATSRHSJacobianLocalFn` for the calling sequence
257b398ffe5SPeter Brune - ctx  - optional context for local jacobian evaluation
25824989b8cSPeter Brune 
25924989b8cSPeter Brune   Level: beginner
26024989b8cSPeter Brune 
261*8434afd1SBarry Smith .seealso: [](ch_ts), `DMDA`, `DMDATSRHSJacobianLocalFn`, `DMTSSetRHSJacobian()`,
26214d0ab18SJacob Faibussowitsch `DMDATSSetRHSFunctionLocal()`, `DMDASNESSetJacobianLocal()`
26324989b8cSPeter Brune @*/
264*8434afd1SBarry Smith PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm, DMDATSRHSJacobianLocalFn *func, void *ctx)
265d71ae5a4SJacob Faibussowitsch {
266942e3340SBarry Smith   DMTS     sdm;
267942e3340SBarry Smith   DMTS_DA *dmdats;
26824989b8cSPeter Brune 
26924989b8cSPeter Brune   PetscFunctionBegin;
27024989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2719566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &sdm));
2729566063dSJacob Faibussowitsch   PetscCall(DMDATSGetContext(dm, sdm, &dmdats));
27324989b8cSPeter Brune   dmdats->rhsjacobianlocal    = func;
27424989b8cSPeter Brune   dmdats->rhsjacobianlocalctx = ctx;
2759566063dSJacob Faibussowitsch   PetscCall(DMTSSetRHSJacobian(dm, TSComputeRHSJacobian_DMDA, dmdats));
2763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27724989b8cSPeter Brune }
27824989b8cSPeter Brune 
27924989b8cSPeter Brune /*@C
280bcf0153eSBarry Smith   DMDATSSetIFunctionLocal - set a local residual evaluation function for use with `DMDA`
28124989b8cSPeter Brune 
28224989b8cSPeter Brune   Logically Collective
28324989b8cSPeter Brune 
2844165533cSJose E. Roman   Input Parameters:
285bcf0153eSBarry Smith + dm    - `DM` to associate callback with
28614d0ab18SJacob Faibussowitsch . imode - the insert mode of the function
287*8434afd1SBarry Smith . func  - local residual evaluation, see `DMDATSIFunctionLocalFn` for the calling sequence
28824989b8cSPeter Brune - ctx   - optional context for local residual evaluation
28924989b8cSPeter Brune 
29024989b8cSPeter Brune   Level: beginner
29124989b8cSPeter Brune 
292*8434afd1SBarry Smith .seealso: [](ch_ts), `DMDA`, `DMDATSIFunctionLocalFn`, `DMTSSetIFunction()`,
29314d0ab18SJacob Faibussowitsch `DMDATSSetIJacobianLocal()`, `DMDASNESSetFunctionLocal()`
29424989b8cSPeter Brune @*/
295*8434afd1SBarry Smith PetscErrorCode DMDATSSetIFunctionLocal(DM dm, InsertMode imode, DMDATSIFunctionLocalFn *func, void *ctx)
296d71ae5a4SJacob Faibussowitsch {
297942e3340SBarry Smith   DMTS     sdm;
298942e3340SBarry Smith   DMTS_DA *dmdats;
29924989b8cSPeter Brune 
30024989b8cSPeter Brune   PetscFunctionBegin;
30124989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3029566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &sdm));
3039566063dSJacob Faibussowitsch   PetscCall(DMDATSGetContext(dm, sdm, &dmdats));
30424989b8cSPeter Brune   dmdats->ifunctionlocalimode = imode;
30524989b8cSPeter Brune   dmdats->ifunctionlocal      = func;
30624989b8cSPeter Brune   dmdats->ifunctionlocalctx   = ctx;
3079566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMDA, dmdats));
3083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30924989b8cSPeter Brune }
31024989b8cSPeter Brune 
31124989b8cSPeter Brune /*@C
312bcf0153eSBarry Smith   DMDATSSetIJacobianLocal - set a local residual evaluation function for use with `DMDA`
31324989b8cSPeter Brune 
31424989b8cSPeter Brune   Logically Collective
31524989b8cSPeter Brune 
3164165533cSJose E. Roman   Input Parameters:
317bcf0153eSBarry Smith + dm   - `DM` to associate callback with
318*8434afd1SBarry Smith . func - local residual evaluation, see `DMDATSIJacobianLocalFn` for the calling sequence
31924989b8cSPeter Brune - ctx  - optional context for local residual evaluation
32024989b8cSPeter Brune 
32124989b8cSPeter Brune   Level: beginner
32224989b8cSPeter Brune 
323*8434afd1SBarry Smith .seealso: [](ch_ts), `DMDA`, `DMDATSIJacobianLocalFn`, `DMTSSetJacobian()`,
32414d0ab18SJacob Faibussowitsch `DMDATSSetIFunctionLocal()`, `DMDASNESSetJacobianLocal()`
32524989b8cSPeter Brune @*/
326*8434afd1SBarry Smith PetscErrorCode DMDATSSetIJacobianLocal(DM dm, DMDATSIJacobianLocalFn *func, void *ctx)
327d71ae5a4SJacob Faibussowitsch {
328942e3340SBarry Smith   DMTS     sdm;
329942e3340SBarry Smith   DMTS_DA *dmdats;
33024989b8cSPeter Brune 
33124989b8cSPeter Brune   PetscFunctionBegin;
33224989b8cSPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3339566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &sdm));
3349566063dSJacob Faibussowitsch   PetscCall(DMDATSGetContext(dm, sdm, &dmdats));
33524989b8cSPeter Brune   dmdats->ijacobianlocal    = func;
33624989b8cSPeter Brune   dmdats->ijacobianlocalctx = ctx;
3379566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMDA, dmdats));
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33924989b8cSPeter Brune }
340d1212d36SBarry Smith 
341d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDMDARayDestroy(void **mctx)
342d71ae5a4SJacob Faibussowitsch {
343d1212d36SBarry Smith   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *)*mctx;
344d1212d36SBarry Smith 
345d1212d36SBarry Smith   PetscFunctionBegin;
3469566063dSJacob Faibussowitsch   if (rayctx->lgctx) PetscCall(TSMonitorLGCtxDestroy(&rayctx->lgctx));
3479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rayctx->ray));
3489566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&rayctx->scatter));
3499566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&rayctx->viewer));
3509566063dSJacob Faibussowitsch   PetscCall(PetscFree(rayctx));
3513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
352d1212d36SBarry Smith }
353d1212d36SBarry Smith 
354d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDMDARay(TS ts, PetscInt steps, PetscReal time, Vec u, void *mctx)
355d71ae5a4SJacob Faibussowitsch {
356d1212d36SBarry Smith   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *)mctx;
357d1212d36SBarry Smith   Vec                  solution;
358d1212d36SBarry Smith 
359d1212d36SBarry Smith   PetscFunctionBegin;
3609566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &solution));
3619566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(rayctx->scatter, solution, rayctx->ray, INSERT_VALUES, SCATTER_FORWARD));
3629566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(rayctx->scatter, solution, rayctx->ray, INSERT_VALUES, SCATTER_FORWARD));
3631baa6e33SBarry Smith   if (rayctx->viewer) PetscCall(VecView(rayctx->ray, rayctx->viewer));
3643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
365d1212d36SBarry Smith }
36651b4a12fSMatthew G. Knepley 
367d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorLGDMDARay(TS ts, PetscInt step, PetscReal ptime, Vec u, void *ctx)
368d71ae5a4SJacob Faibussowitsch {
36951b4a12fSMatthew G. Knepley   TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *)ctx;
37051b4a12fSMatthew G. Knepley   TSMonitorLGCtx       lgctx  = (TSMonitorLGCtx)rayctx->lgctx;
37151b4a12fSMatthew G. Knepley   Vec                  v      = rayctx->ray;
37251b4a12fSMatthew G. Knepley   const PetscScalar   *a;
37351b4a12fSMatthew G. Knepley   PetscInt             dim;
37451b4a12fSMatthew G. Knepley 
37551b4a12fSMatthew G. Knepley   PetscFunctionBegin;
3769566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD));
3779566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD));
37851b4a12fSMatthew G. Knepley   if (!step) {
37951b4a12fSMatthew G. Knepley     PetscDrawAxis axis;
38051b4a12fSMatthew G. Knepley 
3819566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(lgctx->lg, &axis));
3829566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Solution Ray as function of time", "Time", "Solution"));
3839566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(rayctx->ray, &dim));
3849566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetDimension(lgctx->lg, dim));
3859566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGReset(lgctx->lg));
38651b4a12fSMatthew G. Knepley   }
3879566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(v, &a));
38851b4a12fSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX)
38951b4a12fSMatthew G. Knepley   {
39051b4a12fSMatthew G. Knepley     PetscReal *areal;
39151b4a12fSMatthew G. Knepley     PetscInt   i, n;
3929566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v, &n));
3939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &areal));
39451b4a12fSMatthew G. Knepley     for (i = 0; i < n; ++i) areal[i] = PetscRealPart(a[i]);
3959566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddCommonPoint(lgctx->lg, ptime, areal));
3969566063dSJacob Faibussowitsch     PetscCall(PetscFree(areal));
39751b4a12fSMatthew G. Knepley   }
39851b4a12fSMatthew G. Knepley #else
3999566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGAddCommonPoint(lgctx->lg, ptime, a));
40051b4a12fSMatthew G. Knepley #endif
4019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(v, &a));
40251b4a12fSMatthew G. Knepley   if (((lgctx->howoften > 0) && (!(step % lgctx->howoften))) || ((lgctx->howoften == -1) && ts->reason)) {
4039566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(lgctx->lg));
4049566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSave(lgctx->lg));
40551b4a12fSMatthew G. Knepley   }
4063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40751b4a12fSMatthew G. Knepley }
408