xref: /petsc/src/ts/utils/dmlocalts.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
3d433e6cbSMatthew G. Knepley 
4d433e6cbSMatthew G. Knepley typedef struct {
5f2b97b97SMatthew G. Knepley   PetscErrorCode (*boundarylocal)(DM, PetscReal, Vec, Vec, void *);
6d433e6cbSMatthew G. Knepley   PetscErrorCode (*ifunctionlocal)(DM, PetscReal, Vec, Vec, Vec, void *);
7d433e6cbSMatthew G. Knepley   PetscErrorCode (*ijacobianlocal)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
8d433e6cbSMatthew G. Knepley   PetscErrorCode (*rhsfunctionlocal)(DM, PetscReal, Vec, Vec, void *);
9f2b97b97SMatthew G. Knepley   void *boundarylocalctx;
10d433e6cbSMatthew G. Knepley   void *ifunctionlocalctx;
11d433e6cbSMatthew G. Knepley   void *ijacobianlocalctx;
12d433e6cbSMatthew G. Knepley   void *rhsfunctionlocalctx;
13b4937a87SMatthew G. Knepley   Vec   lumpedmassinv;
14b4937a87SMatthew G. Knepley   Mat   mass;
15b4937a87SMatthew G. Knepley   KSP   kspmass;
16d433e6cbSMatthew G. Knepley } DMTS_Local;
17d433e6cbSMatthew G. Knepley 
189371c9d4SSatish Balay static PetscErrorCode DMTSDestroy_DMLocal(DMTS tdm) {
19d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(PetscFree(tdm->data));
21d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
22d433e6cbSMatthew G. Knepley }
23d433e6cbSMatthew G. Knepley 
249371c9d4SSatish Balay static PetscErrorCode DMTSDuplicate_DMLocal(DMTS oldtdm, DMTS tdm) {
25d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
26*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew((DMTS_Local **)&tdm->data));
279566063dSJacob Faibussowitsch   if (oldtdm->data) PetscCall(PetscMemcpy(tdm->data, oldtdm->data, sizeof(DMTS_Local)));
28d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
29d433e6cbSMatthew G. Knepley }
30d433e6cbSMatthew G. Knepley 
319371c9d4SSatish Balay static PetscErrorCode DMLocalTSGetContext(DM dm, DMTS tdm, DMTS_Local **dmlocalts) {
32d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
33d433e6cbSMatthew G. Knepley   *dmlocalts = NULL;
34d433e6cbSMatthew G. Knepley   if (!tdm->data) {
35*4dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMTS_Local **)&tdm->data));
36d433e6cbSMatthew G. Knepley 
37d433e6cbSMatthew G. Knepley     tdm->ops->destroy   = DMTSDestroy_DMLocal;
38d433e6cbSMatthew G. Knepley     tdm->ops->duplicate = DMTSDuplicate_DMLocal;
39d433e6cbSMatthew G. Knepley   }
40d433e6cbSMatthew G. Knepley   *dmlocalts = (DMTS_Local *)tdm->data;
41d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
42d433e6cbSMatthew G. Knepley }
43d433e6cbSMatthew G. Knepley 
449371c9d4SSatish Balay static PetscErrorCode TSComputeIFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, Vec F, void *ctx) {
45d433e6cbSMatthew G. Knepley   DM          dm;
46d433e6cbSMatthew G. Knepley   Vec         locX, locX_t, locF;
47d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts = (DMTS_Local *)ctx;
48d433e6cbSMatthew G. Knepley 
49d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
50d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
51d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(X, VEC_CLASSID, 3);
52d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(X_t, VEC_CLASSID, 4);
53d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(F, VEC_CLASSID, 5);
549566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
559566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX));
569566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX_t));
579566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locF));
589566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locX));
599566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locX_t));
609566063dSJacob Faibussowitsch   if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, locX_t, dmlocalts->boundarylocalctx));
619566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
629566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
639566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t));
649566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t));
659566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locF));
66d433e6cbSMatthew G. Knepley   CHKMEMQ;
679566063dSJacob Faibussowitsch   PetscCall((*dmlocalts->ifunctionlocal)(dm, time, locX, locX_t, locF, dmlocalts->ifunctionlocalctx));
68d433e6cbSMatthew G. Knepley   CHKMEMQ;
699566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
709566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F));
719566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F));
729566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX));
739566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX_t));
749566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locF));
75d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
76d433e6cbSMatthew G. Knepley }
77d433e6cbSMatthew G. Knepley 
789371c9d4SSatish Balay static PetscErrorCode TSComputeRHSFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec F, void *ctx) {
79d433e6cbSMatthew G. Knepley   DM          dm;
80b4937a87SMatthew G. Knepley   Vec         locX, locF;
81d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts = (DMTS_Local *)ctx;
82d433e6cbSMatthew G. Knepley 
83d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
84d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
85d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(X, VEC_CLASSID, 3);
86064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(F, VEC_CLASSID, 4);
879566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
889566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX));
899566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locF));
909566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locX));
919566063dSJacob Faibussowitsch   if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, NULL, dmlocalts->boundarylocalctx));
929566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
939566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
949566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locF));
95b4937a87SMatthew G. Knepley   CHKMEMQ;
969566063dSJacob Faibussowitsch   PetscCall((*dmlocalts->rhsfunctionlocal)(dm, time, locX, locF, dmlocalts->rhsfunctionlocalctx));
97b4937a87SMatthew G. Knepley   CHKMEMQ;
989566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
999566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F));
1009566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F));
101b4937a87SMatthew G. Knepley   if (dmlocalts->lumpedmassinv) {
1029566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(F, dmlocalts->lumpedmassinv, F));
103b4937a87SMatthew G. Knepley   } else if (dmlocalts->kspmass) {
104b4937a87SMatthew G. Knepley     Vec tmp;
105b4937a87SMatthew G. Knepley 
1069566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &tmp));
1079566063dSJacob Faibussowitsch     PetscCall(KSPSolve(dmlocalts->kspmass, F, tmp));
1089566063dSJacob Faibussowitsch     PetscCall(VecCopy(tmp, F));
1099566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &tmp));
110b4937a87SMatthew G. Knepley   }
1119566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX));
1129566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locF));
113d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
114d433e6cbSMatthew G. Knepley }
115d433e6cbSMatthew G. Knepley 
1169371c9d4SSatish Balay static PetscErrorCode TSComputeIJacobian_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, PetscReal a, Mat A, Mat B, void *ctx) {
117d433e6cbSMatthew G. Knepley   DM          dm;
118d433e6cbSMatthew G. Knepley   Vec         locX, locX_t;
119d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts = (DMTS_Local *)ctx;
120d433e6cbSMatthew G. Knepley 
121d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
1229566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
123d433e6cbSMatthew G. Knepley   if (dmlocalts->ijacobianlocal) {
1249566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &locX));
1259566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &locX_t));
1269566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(locX));
1279566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(locX_t));
1289566063dSJacob Faibussowitsch     if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, locX_t, dmlocalts->boundarylocalctx));
1299566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
1309566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
1319566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t));
1329566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t));
133d433e6cbSMatthew G. Knepley     CHKMEMQ;
1349566063dSJacob Faibussowitsch     PetscCall((*dmlocalts->ijacobianlocal)(dm, time, locX, locX_t, a, A, B, dmlocalts->ijacobianlocalctx));
135d433e6cbSMatthew G. Knepley     CHKMEMQ;
1369566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &locX));
1379566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &locX_t));
138d433e6cbSMatthew G. Knepley   } else {
139d433e6cbSMatthew G. Knepley     MatFDColoring fdcoloring;
1409566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring));
141d433e6cbSMatthew G. Knepley     if (!fdcoloring) {
142d433e6cbSMatthew G. Knepley       ISColoring coloring;
143d433e6cbSMatthew G. Knepley 
1449566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring));
1459566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring));
1469566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&coloring));
147d433e6cbSMatthew G. Knepley       switch (dm->coloringtype) {
1489371c9d4SSatish Balay       case IS_COLORING_GLOBAL: PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))TSComputeIFunction_DMLocal, dmlocalts)); break;
14998921bdaSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]);
150d433e6cbSMatthew G. Knepley       }
1519566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix));
1529566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1539566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring));
1549566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring));
1559566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)fdcoloring));
156d433e6cbSMatthew G. Knepley 
157d433e6cbSMatthew G. Knepley       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
158d433e6cbSMatthew G. Knepley        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
159d433e6cbSMatthew G. Knepley        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
160d433e6cbSMatthew G. Knepley        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
161d433e6cbSMatthew G. Knepley        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
162d433e6cbSMatthew G. Knepley        */
1639566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dm));
164d433e6cbSMatthew G. Knepley     }
1659566063dSJacob Faibussowitsch     PetscCall(MatFDColoringApply(B, fdcoloring, X, ts));
166d433e6cbSMatthew G. Knepley   }
167d433e6cbSMatthew G. Knepley   /* This will be redundant if the user called both, but it's too common to forget. */
168d433e6cbSMatthew G. Knepley   if (A != B) {
1699566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1709566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
171d433e6cbSMatthew G. Knepley   }
172d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
173d433e6cbSMatthew G. Knepley }
174d433e6cbSMatthew G. Knepley 
1757eb36faeSMatthew G. Knepley /*@C
176f2b97b97SMatthew G. Knepley   DMTSSetBoundaryLocal - set the function for essential boundary data for a local implicit function evaluation.
177a40652d4SToby Isaac     It should set the essential boundary data for the local portion of the solution X, as well its time derivative X_t (if it is not NULL).
178a40652d4SToby Isaac     Vectors are initialized to zero before this function, so it is only needed for non homogeneous data.
179a40652d4SToby Isaac 
180a40652d4SToby Isaac   Note that this function is somewhat optional: boundary data could potentially be inserted by a function passed to
181a40652d4SToby Isaac   DMTSSetIFunctionLocal().  The use case for this function is for discretizations with constraints (see
182a40652d4SToby Isaac   DMGetDefaultConstraints()): this function inserts boundary values before constraint interpolation.
1837eb36faeSMatthew G. Knepley 
1847eb36faeSMatthew G. Knepley   Logically Collective
1857eb36faeSMatthew G. Knepley 
1864165533cSJose E. Roman   Input Parameters:
1877eb36faeSMatthew G. Knepley + dm   - DM to associate callback with
188a40652d4SToby Isaac . func - local function evaluation
189a40652d4SToby Isaac - ctx  - context for function evaluation
1907eb36faeSMatthew G. Knepley 
1917eb36faeSMatthew G. Knepley   Level: intermediate
1927eb36faeSMatthew G. Knepley 
193db781477SPatrick Sanan .seealso: `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()`
1947eb36faeSMatthew G. Knepley @*/
1959371c9d4SSatish Balay PetscErrorCode DMTSSetBoundaryLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx) {
1967eb36faeSMatthew G. Knepley   DMTS        tdm;
1977eb36faeSMatthew G. Knepley   DMTS_Local *dmlocalts;
1987eb36faeSMatthew G. Knepley 
1997eb36faeSMatthew G. Knepley   PetscFunctionBegin;
2007eb36faeSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2019566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
2029566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
2037eb36faeSMatthew G. Knepley 
2047eb36faeSMatthew G. Knepley   dmlocalts->boundarylocal    = func;
2057eb36faeSMatthew G. Knepley   dmlocalts->boundarylocalctx = ctx;
206a40652d4SToby Isaac 
2077eb36faeSMatthew G. Knepley   PetscFunctionReturn(0);
2087eb36faeSMatthew G. Knepley }
2097eb36faeSMatthew G. Knepley 
210d433e6cbSMatthew G. Knepley /*@C
2110c9409e4SMatthew G. Knepley   DMTSGetIFunctionLocal - get the local implicit function evaluation function. This function is called with local vector
2120c9409e4SMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
2130c9409e4SMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
2140c9409e4SMatthew G. Knepley 
2150c9409e4SMatthew G. Knepley   Logically Collective
2160c9409e4SMatthew G. Knepley 
2170c9409e4SMatthew G. Knepley   Input Parameter:
2180c9409e4SMatthew G. Knepley . dm   - DM to associate callback with
2190c9409e4SMatthew G. Knepley 
2200c9409e4SMatthew G. Knepley   Output Parameters:
2210c9409e4SMatthew G. Knepley + func - local function evaluation
2220c9409e4SMatthew G. Knepley - ctx  - context for function evaluation
2230c9409e4SMatthew G. Knepley 
2240c9409e4SMatthew G. Knepley   Level: beginner
2250c9409e4SMatthew G. Knepley 
226db781477SPatrick Sanan .seealso: `DMTSSetIFunctionLocal(()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()`
2270c9409e4SMatthew G. Knepley @*/
2289371c9d4SSatish Balay PetscErrorCode DMTSGetIFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, Vec, void *), void **ctx) {
2290c9409e4SMatthew G. Knepley   DMTS           tdm;
2300c9409e4SMatthew G. Knepley   DMTS_Local    *dmlocalts;
2310c9409e4SMatthew G. Knepley   PetscErrorCode ierr;
2320c9409e4SMatthew G. Knepley 
2330c9409e4SMatthew G. Knepley   PetscFunctionBegin;
2340c9409e4SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2359371c9d4SSatish Balay   ierr = DMGetDMTS(dm, &tdm);
2369371c9d4SSatish Balay   CHKERRQ(ierr);
2379371c9d4SSatish Balay   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);
2389371c9d4SSatish Balay   CHKERRQ(ierr);
2399371c9d4SSatish Balay   if (func) {
2409371c9d4SSatish Balay     PetscValidPointer(func, 2);
2419371c9d4SSatish Balay     *func = dmlocalts->ifunctionlocal;
2429371c9d4SSatish Balay   }
2439371c9d4SSatish Balay   if (ctx) {
2449371c9d4SSatish Balay     PetscValidPointer(ctx, 3);
2459371c9d4SSatish Balay     *ctx = dmlocalts->ifunctionlocalctx;
2469371c9d4SSatish Balay   }
2470c9409e4SMatthew G. Knepley   PetscFunctionReturn(0);
2480c9409e4SMatthew G. Knepley }
2490c9409e4SMatthew G. Knepley 
2500c9409e4SMatthew G. Knepley /*@C
251d433e6cbSMatthew G. Knepley   DMTSSetIFunctionLocal - set a local implicit function evaluation function. This function is called with local vector
252d433e6cbSMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
253d433e6cbSMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
254d433e6cbSMatthew G. Knepley 
255d433e6cbSMatthew G. Knepley   Logically Collective
256d433e6cbSMatthew G. Knepley 
2574165533cSJose E. Roman   Input Parameters:
258d433e6cbSMatthew G. Knepley + dm   - DM to associate callback with
259d433e6cbSMatthew G. Knepley . func - local function evaluation
260d433e6cbSMatthew G. Knepley - ctx  - context for function evaluation
261d433e6cbSMatthew G. Knepley 
262d433e6cbSMatthew G. Knepley   Level: beginner
263d433e6cbSMatthew G. Knepley 
264db781477SPatrick Sanan .seealso: `DMTSGetIFunctionLocal()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()`
265d433e6cbSMatthew G. Knepley @*/
2669371c9d4SSatish Balay PetscErrorCode DMTSSetIFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, Vec, void *), void *ctx) {
267d433e6cbSMatthew G. Knepley   DMTS        tdm;
268d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts;
269d433e6cbSMatthew G. Knepley 
270d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
271d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2729566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
2739566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
274d433e6cbSMatthew G. Knepley 
275d433e6cbSMatthew G. Knepley   dmlocalts->ifunctionlocal    = func;
276d433e6cbSMatthew G. Knepley   dmlocalts->ifunctionlocalctx = ctx;
277d433e6cbSMatthew G. Knepley 
2789566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMLocal, dmlocalts));
279d433e6cbSMatthew G. Knepley   if (!tdm->ops->ijacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
2809566063dSJacob Faibussowitsch     PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts));
281d433e6cbSMatthew G. Knepley   }
282d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
283d433e6cbSMatthew G. Knepley }
284d433e6cbSMatthew G. Knepley 
285d433e6cbSMatthew G. Knepley /*@C
2860c9409e4SMatthew G. Knepley   DMTSGetIJacobianLocal - get a local Jacobian evaluation function
2870c9409e4SMatthew G. Knepley 
2880c9409e4SMatthew G. Knepley   Logically Collective
2890c9409e4SMatthew G. Knepley 
2900c9409e4SMatthew G. Knepley   Input Parameter:
2910c9409e4SMatthew G. Knepley . dm - DM to associate callback with
2920c9409e4SMatthew G. Knepley 
2930c9409e4SMatthew G. Knepley   Output Parameters:
2940c9409e4SMatthew G. Knepley + func - local Jacobian evaluation
2950c9409e4SMatthew G. Knepley - ctx - optional context for local Jacobian evaluation
2960c9409e4SMatthew G. Knepley 
2970c9409e4SMatthew G. Knepley   Level: beginner
2980c9409e4SMatthew G. Knepley 
299db781477SPatrick Sanan .seealso: `DMTSSetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()`
3000c9409e4SMatthew G. Knepley @*/
3019371c9d4SSatish Balay PetscErrorCode DMTSGetIJacobianLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **ctx) {
3020c9409e4SMatthew G. Knepley   DMTS           tdm;
3030c9409e4SMatthew G. Knepley   DMTS_Local    *dmlocalts;
3040c9409e4SMatthew G. Knepley   PetscErrorCode ierr;
3050c9409e4SMatthew G. Knepley 
3060c9409e4SMatthew G. Knepley   PetscFunctionBegin;
3070c9409e4SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3089371c9d4SSatish Balay   ierr = DMGetDMTS(dm, &tdm);
3099371c9d4SSatish Balay   CHKERRQ(ierr);
3109371c9d4SSatish Balay   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);
3119371c9d4SSatish Balay   CHKERRQ(ierr);
3129371c9d4SSatish Balay   if (func) {
3139371c9d4SSatish Balay     PetscValidPointer(func, 2);
3149371c9d4SSatish Balay     *func = dmlocalts->ijacobianlocal;
3159371c9d4SSatish Balay   }
3169371c9d4SSatish Balay   if (ctx) {
3179371c9d4SSatish Balay     PetscValidPointer(ctx, 3);
3189371c9d4SSatish Balay     *ctx = dmlocalts->ijacobianlocalctx;
3199371c9d4SSatish Balay   }
3200c9409e4SMatthew G. Knepley   PetscFunctionReturn(0);
3210c9409e4SMatthew G. Knepley }
3220c9409e4SMatthew G. Knepley 
3230c9409e4SMatthew G. Knepley /*@C
324d433e6cbSMatthew G. Knepley   DMTSSetIJacobianLocal - set a local Jacobian evaluation function
325d433e6cbSMatthew G. Knepley 
326d433e6cbSMatthew G. Knepley   Logically Collective
327d433e6cbSMatthew G. Knepley 
3284165533cSJose E. Roman   Input Parameters:
329d433e6cbSMatthew G. Knepley + dm - DM to associate callback with
330d433e6cbSMatthew G. Knepley . func - local Jacobian evaluation
331d433e6cbSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation
332d433e6cbSMatthew G. Knepley 
333d433e6cbSMatthew G. Knepley   Level: beginner
334d433e6cbSMatthew G. Knepley 
335db781477SPatrick Sanan .seealso: `DMTSGetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()`
336d433e6cbSMatthew G. Knepley @*/
3379371c9d4SSatish Balay PetscErrorCode DMTSSetIJacobianLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *ctx) {
338d433e6cbSMatthew G. Knepley   DMTS        tdm;
339d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts;
340d433e6cbSMatthew G. Knepley 
341d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
342d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3439566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
3449566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
345d433e6cbSMatthew G. Knepley 
346d433e6cbSMatthew G. Knepley   dmlocalts->ijacobianlocal    = func;
347d433e6cbSMatthew G. Knepley   dmlocalts->ijacobianlocalctx = ctx;
348d433e6cbSMatthew G. Knepley 
3499566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts));
350d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
351d433e6cbSMatthew G. Knepley }
352d433e6cbSMatthew G. Knepley 
353d433e6cbSMatthew G. Knepley /*@C
3540c9409e4SMatthew G. Knepley   DMTSGetRHSFunctionLocal - get a local rhs function evaluation function. This function is called with local vector
3550c9409e4SMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
3560c9409e4SMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
3570c9409e4SMatthew G. Knepley 
3580c9409e4SMatthew G. Knepley   Logically Collective
3590c9409e4SMatthew G. Knepley 
3600c9409e4SMatthew G. Knepley   Input Parameter:
3610c9409e4SMatthew G. Knepley . dm   - DM to associate callback with
3620c9409e4SMatthew G. Knepley 
3630c9409e4SMatthew G. Knepley   Output Parameters:
3640c9409e4SMatthew G. Knepley + func - local function evaluation
3650c9409e4SMatthew G. Knepley - ctx  - context for function evaluation
3660c9409e4SMatthew G. Knepley 
3670c9409e4SMatthew G. Knepley   Level: beginner
3680c9409e4SMatthew G. Knepley 
369db781477SPatrick Sanan .seealso: `DMTSSetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()`
3700c9409e4SMatthew G. Knepley @*/
3719371c9d4SSatish Balay PetscErrorCode DMTSGetRHSFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, void *), void **ctx) {
3720c9409e4SMatthew G. Knepley   DMTS           tdm;
3730c9409e4SMatthew G. Knepley   DMTS_Local    *dmlocalts;
3740c9409e4SMatthew G. Knepley   PetscErrorCode ierr;
3750c9409e4SMatthew G. Knepley 
3760c9409e4SMatthew G. Knepley   PetscFunctionBegin;
3770c9409e4SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3789371c9d4SSatish Balay   ierr = DMGetDMTS(dm, &tdm);
3799371c9d4SSatish Balay   CHKERRQ(ierr);
3809371c9d4SSatish Balay   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);
3819371c9d4SSatish Balay   CHKERRQ(ierr);
3829371c9d4SSatish Balay   if (func) {
3839371c9d4SSatish Balay     PetscValidPointer(func, 2);
3849371c9d4SSatish Balay     *func = dmlocalts->rhsfunctionlocal;
3859371c9d4SSatish Balay   }
3869371c9d4SSatish Balay   if (ctx) {
3879371c9d4SSatish Balay     PetscValidPointer(ctx, 3);
3889371c9d4SSatish Balay     *ctx = dmlocalts->rhsfunctionlocalctx;
3899371c9d4SSatish Balay   }
3900c9409e4SMatthew G. Knepley   PetscFunctionReturn(0);
3910c9409e4SMatthew G. Knepley }
3920c9409e4SMatthew G. Knepley 
3930c9409e4SMatthew G. Knepley /*@C
394d433e6cbSMatthew G. Knepley   DMTSSetRHSFunctionLocal - set a local rhs function evaluation function. This function is called with local vector
395d433e6cbSMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
396d433e6cbSMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
397d433e6cbSMatthew G. Knepley 
398d433e6cbSMatthew G. Knepley   Logically Collective
399d433e6cbSMatthew G. Knepley 
4004165533cSJose E. Roman   Input Parameters:
401d433e6cbSMatthew G. Knepley + dm   - DM to associate callback with
402d433e6cbSMatthew G. Knepley . func - local function evaluation
403d433e6cbSMatthew G. Knepley - ctx  - context for function evaluation
404d433e6cbSMatthew G. Knepley 
405d433e6cbSMatthew G. Knepley   Level: beginner
406d433e6cbSMatthew G. Knepley 
407db781477SPatrick Sanan .seealso: `DMTSGetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()`
408d433e6cbSMatthew G. Knepley @*/
4099371c9d4SSatish Balay PetscErrorCode DMTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx) {
410d433e6cbSMatthew G. Knepley   DMTS        tdm;
411d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts;
412d433e6cbSMatthew G. Knepley 
413d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
414d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4159566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
4169566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
417d433e6cbSMatthew G. Knepley 
418d433e6cbSMatthew G. Knepley   dmlocalts->rhsfunctionlocal    = func;
419d433e6cbSMatthew G. Knepley   dmlocalts->rhsfunctionlocalctx = ctx;
420d433e6cbSMatthew G. Knepley 
4219566063dSJacob Faibussowitsch   PetscCall(DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMLocal, dmlocalts));
422d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
423d433e6cbSMatthew G. Knepley }
424b4937a87SMatthew G. Knepley 
425b4937a87SMatthew G. Knepley /*@C
426b4937a87SMatthew G. Knepley   DMTSCreateRHSMassMatrix - This creates the mass matrix associated with the given DM, and a solver to invert it, and stores them in the DMTS context.
427b4937a87SMatthew G. Knepley 
428b4937a87SMatthew G. Knepley   Collective on dm
429b4937a87SMatthew G. Knepley 
430b4937a87SMatthew G. Knepley   Input Parameters:
431b4937a87SMatthew G. Knepley . dm   - DM providing the mass matrix
432b4937a87SMatthew G. Knepley 
433b4937a87SMatthew G. Knepley   Note: The idea here is that an explicit system can be given a mass matrix, based on the DM, which is inverted on the RHS at each step.
434b4937a87SMatthew G. Knepley 
435b4937a87SMatthew G. Knepley   Level: developer
436b4937a87SMatthew G. Knepley 
437db781477SPatrick Sanan .seealso: `DMTSCreateRHSMassMatrixLumped()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS`
438b4937a87SMatthew G. Knepley @*/
4399371c9d4SSatish Balay PetscErrorCode DMTSCreateRHSMassMatrix(DM dm) {
440b4937a87SMatthew G. Knepley   DMTS        tdm;
441b4937a87SMatthew G. Knepley   DMTS_Local *dmlocalts;
442b4937a87SMatthew G. Knepley   const char *prefix;
443b4937a87SMatthew G. Knepley 
444b4937a87SMatthew G. Knepley   PetscFunctionBegin;
445b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4469566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
4479566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
4489566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dm, dm, &dmlocalts->mass));
4499566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm), &dmlocalts->kspmass));
4509566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
4519566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(dmlocalts->kspmass, prefix));
4529566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(dmlocalts->kspmass, "mass_"));
4539566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(dmlocalts->kspmass));
4549566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(dmlocalts->kspmass, dmlocalts->mass, dmlocalts->mass));
455b4937a87SMatthew G. Knepley   PetscFunctionReturn(0);
456b4937a87SMatthew G. Knepley }
457b4937a87SMatthew G. Knepley 
458b4937a87SMatthew G. Knepley /*@C
459b4937a87SMatthew G. Knepley   DMTSCreateRHSMassMatrixLumped - This creates the lumped mass matrix associated with the given DM, and a solver to invert it, and stores them in the DMTS context.
460b4937a87SMatthew G. Knepley 
461b4937a87SMatthew G. Knepley   Collective on dm
462b4937a87SMatthew G. Knepley 
463b4937a87SMatthew G. Knepley   Input Parameters:
464b4937a87SMatthew G. Knepley . dm   - DM providing the mass matrix
465b4937a87SMatthew G. Knepley 
466b4937a87SMatthew G. Knepley   Note: The idea here is that an explicit system can be given a mass matrix, based on the DM, which is inverted on the RHS at each step.
467b4937a87SMatthew G. Knepley   Since the matrix is lumped, inversion is trivial.
468b4937a87SMatthew G. Knepley 
469b4937a87SMatthew G. Knepley   Level: developer
470b4937a87SMatthew G. Knepley 
471db781477SPatrick Sanan .seealso: `DMTSCreateRHSMassMatrix()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS`
472b4937a87SMatthew G. Knepley @*/
4739371c9d4SSatish Balay PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM dm) {
474b4937a87SMatthew G. Knepley   DMTS        tdm;
475b4937a87SMatthew G. Knepley   DMTS_Local *dmlocalts;
476b4937a87SMatthew G. Knepley 
477b4937a87SMatthew G. Knepley   PetscFunctionBegin;
478b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4799566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
4809566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
4819566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrixLumped(dm, &dmlocalts->lumpedmassinv));
4829566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(dmlocalts->lumpedmassinv));
4839566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(dmlocalts->lumpedmassinv, NULL, "-lumped_mass_inv_view"));
484b4937a87SMatthew G. Knepley   PetscFunctionReturn(0);
485b4937a87SMatthew G. Knepley }
486b4937a87SMatthew G. Knepley 
487b4937a87SMatthew G. Knepley /*@C
488b4937a87SMatthew G. Knepley   DMTSDestroyRHSMassMatrix - Destroys the mass matrix and solver stored in the DMTS context, if they exist.
489b4937a87SMatthew G. Knepley 
490b4937a87SMatthew G. Knepley   Logically Collective
491b4937a87SMatthew G. Knepley 
492b4937a87SMatthew G. Knepley   Input Parameters:
493b4937a87SMatthew G. Knepley . dm   - DM providing the mass matrix
494b4937a87SMatthew G. Knepley 
495b4937a87SMatthew G. Knepley   Level: developer
496b4937a87SMatthew G. Knepley 
497db781477SPatrick Sanan .seealso: `DMTSCreateRHSMassMatrixLumped()`, `DMCreateMassMatrix()`, `DMCreateMassMatrix()`, `DMTS`
498b4937a87SMatthew G. Knepley @*/
4999371c9d4SSatish Balay PetscErrorCode DMTSDestroyRHSMassMatrix(DM dm) {
500b4937a87SMatthew G. Knepley   DMTS        tdm;
501b4937a87SMatthew G. Knepley   DMTS_Local *dmlocalts;
502b4937a87SMatthew G. Knepley 
503b4937a87SMatthew G. Knepley   PetscFunctionBegin;
504b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5059566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
5069566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
5079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&dmlocalts->lumpedmassinv));
5089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&dmlocalts->mass));
5099566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&dmlocalts->kspmass));
510b4937a87SMatthew G. Knepley   PetscFunctionReturn(0);
511b4937a87SMatthew G. Knepley }
512