xref: /petsc/src/ts/utils/dmlocalts.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
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 
18d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSDestroy_DMLocal(DMTS tdm)
19d71ae5a4SJacob Faibussowitsch {
20d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(PetscFree(tdm->data));
223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23d433e6cbSMatthew G. Knepley }
24d433e6cbSMatthew G. Knepley 
25d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSDuplicate_DMLocal(DMTS oldtdm, DMTS tdm)
26d71ae5a4SJacob Faibussowitsch {
27d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
284dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew((DMTS_Local **)&tdm->data));
299566063dSJacob Faibussowitsch   if (oldtdm->data) PetscCall(PetscMemcpy(tdm->data, oldtdm->data, sizeof(DMTS_Local)));
303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31d433e6cbSMatthew G. Knepley }
32d433e6cbSMatthew G. Knepley 
33d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMLocalTSGetContext(DM dm, DMTS tdm, DMTS_Local **dmlocalts)
34d71ae5a4SJacob Faibussowitsch {
35d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
36d433e6cbSMatthew G. Knepley   *dmlocalts = NULL;
37d433e6cbSMatthew G. Knepley   if (!tdm->data) {
384dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMTS_Local **)&tdm->data));
39d433e6cbSMatthew G. Knepley 
40d433e6cbSMatthew G. Knepley     tdm->ops->destroy   = DMTSDestroy_DMLocal;
41d433e6cbSMatthew G. Knepley     tdm->ops->duplicate = DMTSDuplicate_DMLocal;
42d433e6cbSMatthew G. Knepley   }
43d433e6cbSMatthew G. Knepley   *dmlocalts = (DMTS_Local *)tdm->data;
443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45d433e6cbSMatthew G. Knepley }
46d433e6cbSMatthew G. Knepley 
47d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeIFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, Vec F, void *ctx)
48d71ae5a4SJacob Faibussowitsch {
49d433e6cbSMatthew G. Knepley   DM          dm;
50d433e6cbSMatthew G. Knepley   Vec         locX, locX_t, locF;
51d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts = (DMTS_Local *)ctx;
52d433e6cbSMatthew G. Knepley 
53d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
54d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
55d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(X, VEC_CLASSID, 3);
56d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(X_t, VEC_CLASSID, 4);
57d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(F, VEC_CLASSID, 5);
589566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
599566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX));
609566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX_t));
619566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locF));
629566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locX));
639566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locX_t));
649566063dSJacob Faibussowitsch   if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, locX_t, dmlocalts->boundarylocalctx));
659566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
669566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
679566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t));
689566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t));
699566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locF));
70d433e6cbSMatthew G. Knepley   CHKMEMQ;
719566063dSJacob Faibussowitsch   PetscCall((*dmlocalts->ifunctionlocal)(dm, time, locX, locX_t, locF, dmlocalts->ifunctionlocalctx));
72d433e6cbSMatthew G. Knepley   CHKMEMQ;
739566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
749566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F));
759566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F));
769566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX));
779566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX_t));
789566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locF));
793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80d433e6cbSMatthew G. Knepley }
81d433e6cbSMatthew G. Knepley 
82d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeRHSFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec F, void *ctx)
83d71ae5a4SJacob Faibussowitsch {
84d433e6cbSMatthew G. Knepley   DM          dm;
85b4937a87SMatthew G. Knepley   Vec         locX, locF;
86d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts = (DMTS_Local *)ctx;
87d433e6cbSMatthew G. Knepley 
88d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
89d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
90d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(X, VEC_CLASSID, 3);
91064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(F, VEC_CLASSID, 4);
929566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
939566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX));
949566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locF));
959566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locX));
969566063dSJacob Faibussowitsch   if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, NULL, dmlocalts->boundarylocalctx));
979566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
989566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
999566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locF));
100b4937a87SMatthew G. Knepley   CHKMEMQ;
1019566063dSJacob Faibussowitsch   PetscCall((*dmlocalts->rhsfunctionlocal)(dm, time, locX, locF, dmlocalts->rhsfunctionlocalctx));
102b4937a87SMatthew G. Knepley   CHKMEMQ;
1039566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
1049566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F));
1059566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F));
106b4937a87SMatthew G. Knepley   if (dmlocalts->lumpedmassinv) {
1079566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(F, dmlocalts->lumpedmassinv, F));
108b4937a87SMatthew G. Knepley   } else if (dmlocalts->kspmass) {
109b4937a87SMatthew G. Knepley     Vec tmp;
110b4937a87SMatthew G. Knepley 
1119566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &tmp));
1129566063dSJacob Faibussowitsch     PetscCall(KSPSolve(dmlocalts->kspmass, F, tmp));
1139566063dSJacob Faibussowitsch     PetscCall(VecCopy(tmp, F));
1149566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &tmp));
115b4937a87SMatthew G. Knepley   }
1169566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX));
1179566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locF));
1183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119d433e6cbSMatthew G. Knepley }
120d433e6cbSMatthew G. Knepley 
121d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeIJacobian_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, PetscReal a, Mat A, Mat B, void *ctx)
122d71ae5a4SJacob Faibussowitsch {
123d433e6cbSMatthew G. Knepley   DM          dm;
124d433e6cbSMatthew G. Knepley   Vec         locX, locX_t;
125d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts = (DMTS_Local *)ctx;
126d433e6cbSMatthew G. Knepley 
127d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
1289566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
129d433e6cbSMatthew G. Knepley   if (dmlocalts->ijacobianlocal) {
1309566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &locX));
1319566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &locX_t));
1329566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(locX));
1339566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(locX_t));
1349566063dSJacob Faibussowitsch     if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, locX_t, dmlocalts->boundarylocalctx));
1359566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
1369566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
1379566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t));
1389566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t));
139d433e6cbSMatthew G. Knepley     CHKMEMQ;
1409566063dSJacob Faibussowitsch     PetscCall((*dmlocalts->ijacobianlocal)(dm, time, locX, locX_t, a, A, B, dmlocalts->ijacobianlocalctx));
141d433e6cbSMatthew G. Knepley     CHKMEMQ;
1429566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &locX));
1439566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &locX_t));
144d433e6cbSMatthew G. Knepley   } else {
145d433e6cbSMatthew G. Knepley     MatFDColoring fdcoloring;
1469566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring));
147d433e6cbSMatthew G. Knepley     if (!fdcoloring) {
148d433e6cbSMatthew G. Knepley       ISColoring coloring;
149d433e6cbSMatthew G. Knepley 
1509566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring));
1519566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring));
1529566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&coloring));
153d433e6cbSMatthew G. Knepley       switch (dm->coloringtype) {
154d71ae5a4SJacob Faibussowitsch       case IS_COLORING_GLOBAL:
155d71ae5a4SJacob Faibussowitsch         PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))TSComputeIFunction_DMLocal, dmlocalts));
156d71ae5a4SJacob Faibussowitsch         break;
157d71ae5a4SJacob Faibussowitsch       default:
158d71ae5a4SJacob Faibussowitsch         SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]);
159d433e6cbSMatthew G. Knepley       }
1609566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix));
1619566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1629566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring));
1639566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring));
1649566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)fdcoloring));
165d433e6cbSMatthew G. Knepley 
166d433e6cbSMatthew G. Knepley       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
167d433e6cbSMatthew G. Knepley        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
168d433e6cbSMatthew G. Knepley        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
169d433e6cbSMatthew G. Knepley        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
170d433e6cbSMatthew G. Knepley        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
171d433e6cbSMatthew G. Knepley        */
1729566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dm));
173d433e6cbSMatthew G. Knepley     }
1749566063dSJacob Faibussowitsch     PetscCall(MatFDColoringApply(B, fdcoloring, X, ts));
175d433e6cbSMatthew G. Knepley   }
176d433e6cbSMatthew G. Knepley   /* This will be redundant if the user called both, but it's too common to forget. */
177d433e6cbSMatthew G. Knepley   if (A != B) {
1789566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1799566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
180d433e6cbSMatthew G. Knepley   }
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
182d433e6cbSMatthew G. Knepley }
183d433e6cbSMatthew G. Knepley 
1847eb36faeSMatthew G. Knepley /*@C
185f2b97b97SMatthew G. Knepley   DMTSSetBoundaryLocal - set the function for essential boundary data for a local implicit function evaluation.
186a40652d4SToby 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).
187a40652d4SToby Isaac     Vectors are initialized to zero before this function, so it is only needed for non homogeneous data.
188a40652d4SToby Isaac 
1897eb36faeSMatthew G. Knepley   Logically Collective
1907eb36faeSMatthew G. Knepley 
1914165533cSJose E. Roman   Input Parameters:
192bcf0153eSBarry Smith + dm   - `DM` to associate callback with
193a40652d4SToby Isaac . func - local function evaluation
194a40652d4SToby Isaac - ctx  - context for function evaluation
1957eb36faeSMatthew G. Knepley 
1967eb36faeSMatthew G. Knepley   Level: intermediate
1977eb36faeSMatthew G. Knepley 
1982fe279fdSBarry Smith   Note:
1992fe279fdSBarry Smith   This function is somewhat optional: boundary data could potentially be inserted by a function passed to
2002fe279fdSBarry Smith   `DMTSSetIFunctionLocal()`.  The use case for this function is for discretizations with constraints (see
2012fe279fdSBarry Smith   `DMGetDefaultConstraints()`): this function inserts boundary values before constraint interpolation.
2022fe279fdSBarry Smith 
203*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()`
2047eb36faeSMatthew G. Knepley @*/
205d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetBoundaryLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx)
206d71ae5a4SJacob Faibussowitsch {
2077eb36faeSMatthew G. Knepley   DMTS        tdm;
2087eb36faeSMatthew G. Knepley   DMTS_Local *dmlocalts;
2097eb36faeSMatthew G. Knepley 
2107eb36faeSMatthew G. Knepley   PetscFunctionBegin;
2117eb36faeSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2129566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
2139566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
2147eb36faeSMatthew G. Knepley 
2157eb36faeSMatthew G. Knepley   dmlocalts->boundarylocal    = func;
2167eb36faeSMatthew G. Knepley   dmlocalts->boundarylocalctx = ctx;
217a40652d4SToby Isaac 
2183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2197eb36faeSMatthew G. Knepley }
2207eb36faeSMatthew G. Knepley 
221d433e6cbSMatthew G. Knepley /*@C
2220c9409e4SMatthew G. Knepley   DMTSGetIFunctionLocal - get the local implicit function evaluation function. This function is called with local vector
2230c9409e4SMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
224bcf0153eSBarry Smith       elements and `DM` will automatically accumulate the overlapping values.
2250c9409e4SMatthew G. Knepley 
2260c9409e4SMatthew G. Knepley   Logically Collective
2270c9409e4SMatthew G. Knepley 
2280c9409e4SMatthew G. Knepley   Input Parameter:
229bcf0153eSBarry Smith . dm   - `DM` to associate callback with
2300c9409e4SMatthew G. Knepley 
2310c9409e4SMatthew G. Knepley   Output Parameters:
2320c9409e4SMatthew G. Knepley + func - local function evaluation
2330c9409e4SMatthew G. Knepley - ctx  - context for function evaluation
2340c9409e4SMatthew G. Knepley 
2350c9409e4SMatthew G. Knepley   Level: beginner
2360c9409e4SMatthew G. Knepley 
237*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetIFunctionLocal()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()`
2380c9409e4SMatthew G. Knepley @*/
239d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetIFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, Vec, void *), void **ctx)
240d71ae5a4SJacob Faibussowitsch {
2410c9409e4SMatthew G. Knepley   DMTS           tdm;
2420c9409e4SMatthew G. Knepley   DMTS_Local    *dmlocalts;
2430c9409e4SMatthew G. Knepley   PetscErrorCode ierr;
2440c9409e4SMatthew G. Knepley 
2450c9409e4SMatthew G. Knepley   PetscFunctionBegin;
2460c9409e4SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2479371c9d4SSatish Balay   ierr = DMGetDMTS(dm, &tdm);
2489371c9d4SSatish Balay   CHKERRQ(ierr);
2499371c9d4SSatish Balay   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);
2509371c9d4SSatish Balay   CHKERRQ(ierr);
2519371c9d4SSatish Balay   if (func) {
2529371c9d4SSatish Balay     PetscValidPointer(func, 2);
2539371c9d4SSatish Balay     *func = dmlocalts->ifunctionlocal;
2549371c9d4SSatish Balay   }
2559371c9d4SSatish Balay   if (ctx) {
2569371c9d4SSatish Balay     PetscValidPointer(ctx, 3);
2579371c9d4SSatish Balay     *ctx = dmlocalts->ifunctionlocalctx;
2589371c9d4SSatish Balay   }
2593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2600c9409e4SMatthew G. Knepley }
2610c9409e4SMatthew G. Knepley 
2620c9409e4SMatthew G. Knepley /*@C
263d433e6cbSMatthew G. Knepley   DMTSSetIFunctionLocal - set a local implicit function evaluation function. This function is called with local vector
264d433e6cbSMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
265bcf0153eSBarry Smith       elements and `DM` will automatically accumulate the overlapping values.
266d433e6cbSMatthew G. Knepley 
267d433e6cbSMatthew G. Knepley   Logically Collective
268d433e6cbSMatthew G. Knepley 
2694165533cSJose E. Roman   Input Parameters:
270bcf0153eSBarry Smith + dm   - `DM` to associate callback with
271d433e6cbSMatthew G. Knepley . func - local function evaluation
272d433e6cbSMatthew G. Knepley - ctx  - context for function evaluation
273d433e6cbSMatthew G. Knepley 
274d433e6cbSMatthew G. Knepley   Level: beginner
275d433e6cbSMatthew G. Knepley 
276*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetIFunctionLocal()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()`
277d433e6cbSMatthew G. Knepley @*/
278d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, Vec, void *), void *ctx)
279d71ae5a4SJacob Faibussowitsch {
280d433e6cbSMatthew G. Knepley   DMTS        tdm;
281d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts;
282d433e6cbSMatthew G. Knepley 
283d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
284d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2859566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
2869566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
287d433e6cbSMatthew G. Knepley 
288d433e6cbSMatthew G. Knepley   dmlocalts->ifunctionlocal    = func;
289d433e6cbSMatthew G. Knepley   dmlocalts->ifunctionlocalctx = ctx;
290d433e6cbSMatthew G. Knepley 
2919566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMLocal, dmlocalts));
292d433e6cbSMatthew G. Knepley   if (!tdm->ops->ijacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
2939566063dSJacob Faibussowitsch     PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts));
294d433e6cbSMatthew G. Knepley   }
2953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
296d433e6cbSMatthew G. Knepley }
297d433e6cbSMatthew G. Knepley 
298d433e6cbSMatthew G. Knepley /*@C
2990c9409e4SMatthew G. Knepley   DMTSGetIJacobianLocal - get a local Jacobian evaluation function
3000c9409e4SMatthew G. Knepley 
3010c9409e4SMatthew G. Knepley   Logically Collective
3020c9409e4SMatthew G. Knepley 
3030c9409e4SMatthew G. Knepley   Input Parameter:
304bcf0153eSBarry Smith . dm - `DM` to associate callback with
3050c9409e4SMatthew G. Knepley 
3060c9409e4SMatthew G. Knepley   Output Parameters:
3070c9409e4SMatthew G. Knepley + func - local Jacobian evaluation
3080c9409e4SMatthew G. Knepley - ctx - optional context for local Jacobian evaluation
3090c9409e4SMatthew G. Knepley 
3100c9409e4SMatthew G. Knepley   Level: beginner
3110c9409e4SMatthew G. Knepley 
312*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()`
3130c9409e4SMatthew G. Knepley @*/
314d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetIJacobianLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **ctx)
315d71ae5a4SJacob Faibussowitsch {
3160c9409e4SMatthew G. Knepley   DMTS           tdm;
3170c9409e4SMatthew G. Knepley   DMTS_Local    *dmlocalts;
3180c9409e4SMatthew G. Knepley   PetscErrorCode ierr;
3190c9409e4SMatthew G. Knepley 
3200c9409e4SMatthew G. Knepley   PetscFunctionBegin;
3210c9409e4SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3229371c9d4SSatish Balay   ierr = DMGetDMTS(dm, &tdm);
3239371c9d4SSatish Balay   CHKERRQ(ierr);
3249371c9d4SSatish Balay   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);
3259371c9d4SSatish Balay   CHKERRQ(ierr);
3269371c9d4SSatish Balay   if (func) {
3279371c9d4SSatish Balay     PetscValidPointer(func, 2);
3289371c9d4SSatish Balay     *func = dmlocalts->ijacobianlocal;
3299371c9d4SSatish Balay   }
3309371c9d4SSatish Balay   if (ctx) {
3319371c9d4SSatish Balay     PetscValidPointer(ctx, 3);
3329371c9d4SSatish Balay     *ctx = dmlocalts->ijacobianlocalctx;
3339371c9d4SSatish Balay   }
3343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3350c9409e4SMatthew G. Knepley }
3360c9409e4SMatthew G. Knepley 
3370c9409e4SMatthew G. Knepley /*@C
338d433e6cbSMatthew G. Knepley   DMTSSetIJacobianLocal - set a local Jacobian evaluation function
339d433e6cbSMatthew G. Knepley 
340d433e6cbSMatthew G. Knepley   Logically Collective
341d433e6cbSMatthew G. Knepley 
3424165533cSJose E. Roman   Input Parameters:
343bcf0153eSBarry Smith + dm - `DM` to associate callback with
344d433e6cbSMatthew G. Knepley . func - local Jacobian evaluation
345d433e6cbSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation
346d433e6cbSMatthew G. Knepley 
347d433e6cbSMatthew G. Knepley   Level: beginner
348d433e6cbSMatthew G. Knepley 
349*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()`
350d433e6cbSMatthew G. Knepley @*/
351d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIJacobianLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *ctx)
352d71ae5a4SJacob Faibussowitsch {
353d433e6cbSMatthew G. Knepley   DMTS        tdm;
354d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts;
355d433e6cbSMatthew G. Knepley 
356d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
357d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3589566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
3599566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
360d433e6cbSMatthew G. Knepley 
361d433e6cbSMatthew G. Knepley   dmlocalts->ijacobianlocal    = func;
362d433e6cbSMatthew G. Knepley   dmlocalts->ijacobianlocalctx = ctx;
363d433e6cbSMatthew G. Knepley 
3649566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts));
3653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
366d433e6cbSMatthew G. Knepley }
367d433e6cbSMatthew G. Knepley 
368d433e6cbSMatthew G. Knepley /*@C
3690c9409e4SMatthew G. Knepley   DMTSGetRHSFunctionLocal - get a local rhs function evaluation function. This function is called with local vector
3700c9409e4SMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
371bcf0153eSBarry Smith       elements and `DM` will automatically accumulate the overlapping values.
3720c9409e4SMatthew G. Knepley 
3730c9409e4SMatthew G. Knepley   Logically Collective
3740c9409e4SMatthew G. Knepley 
3750c9409e4SMatthew G. Knepley   Input Parameter:
376bcf0153eSBarry Smith . dm   - `DM` to associate callback with
3770c9409e4SMatthew G. Knepley 
3780c9409e4SMatthew G. Knepley   Output Parameters:
3790c9409e4SMatthew G. Knepley + func - local function evaluation
3800c9409e4SMatthew G. Knepley - ctx  - context for function evaluation
3810c9409e4SMatthew G. Knepley 
3820c9409e4SMatthew G. Knepley   Level: beginner
3830c9409e4SMatthew G. Knepley 
384*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()`
3850c9409e4SMatthew G. Knepley @*/
386d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetRHSFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, void *), void **ctx)
387d71ae5a4SJacob Faibussowitsch {
3880c9409e4SMatthew G. Knepley   DMTS           tdm;
3890c9409e4SMatthew G. Knepley   DMTS_Local    *dmlocalts;
3900c9409e4SMatthew G. Knepley   PetscErrorCode ierr;
3910c9409e4SMatthew G. Knepley 
3920c9409e4SMatthew G. Knepley   PetscFunctionBegin;
3930c9409e4SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3949371c9d4SSatish Balay   ierr = DMGetDMTS(dm, &tdm);
3959371c9d4SSatish Balay   CHKERRQ(ierr);
3969371c9d4SSatish Balay   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);
3979371c9d4SSatish Balay   CHKERRQ(ierr);
3989371c9d4SSatish Balay   if (func) {
3999371c9d4SSatish Balay     PetscValidPointer(func, 2);
4009371c9d4SSatish Balay     *func = dmlocalts->rhsfunctionlocal;
4019371c9d4SSatish Balay   }
4029371c9d4SSatish Balay   if (ctx) {
4039371c9d4SSatish Balay     PetscValidPointer(ctx, 3);
4049371c9d4SSatish Balay     *ctx = dmlocalts->rhsfunctionlocalctx;
4059371c9d4SSatish Balay   }
4063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4070c9409e4SMatthew G. Knepley }
4080c9409e4SMatthew G. Knepley 
4090c9409e4SMatthew G. Knepley /*@C
410d433e6cbSMatthew G. Knepley   DMTSSetRHSFunctionLocal - set a local rhs function evaluation function. This function is called with local vector
411d433e6cbSMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
412bcf0153eSBarry Smith       elements and `DM` will automatically accumulate the overlapping values.
413d433e6cbSMatthew G. Knepley 
414d433e6cbSMatthew G. Knepley   Logically Collective
415d433e6cbSMatthew G. Knepley 
4164165533cSJose E. Roman   Input Parameters:
417bcf0153eSBarry Smith + dm   - `DM` to associate callback with
418d433e6cbSMatthew G. Knepley . func - local function evaluation
419d433e6cbSMatthew G. Knepley - ctx  - context for function evaluation
420d433e6cbSMatthew G. Knepley 
421d433e6cbSMatthew G. Knepley   Level: beginner
422d433e6cbSMatthew G. Knepley 
423*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()`
424d433e6cbSMatthew G. Knepley @*/
425d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx)
426d71ae5a4SJacob Faibussowitsch {
427d433e6cbSMatthew G. Knepley   DMTS        tdm;
428d433e6cbSMatthew G. Knepley   DMTS_Local *dmlocalts;
429d433e6cbSMatthew G. Knepley 
430d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
431d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4329566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
4339566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
434d433e6cbSMatthew G. Knepley 
435d433e6cbSMatthew G. Knepley   dmlocalts->rhsfunctionlocal    = func;
436d433e6cbSMatthew G. Knepley   dmlocalts->rhsfunctionlocalctx = ctx;
437d433e6cbSMatthew G. Knepley 
4389566063dSJacob Faibussowitsch   PetscCall(DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMLocal, dmlocalts));
4393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
440d433e6cbSMatthew G. Knepley }
441b4937a87SMatthew G. Knepley 
442b4937a87SMatthew G. Knepley /*@C
443bcf0153eSBarry Smith   DMTSCreateRHSMassMatrix - This creates the mass matrix associated with the given `DM`, and a solver to invert it, and stores them in the `DM` context.
444b4937a87SMatthew G. Knepley 
445c3339decSBarry Smith   Collective
446b4937a87SMatthew G. Knepley 
4472fe279fdSBarry Smith   Input Parameter:
448bcf0153eSBarry Smith . dm   - `DM` providing the mass matrix
449b4937a87SMatthew G. Knepley 
450b4937a87SMatthew G. Knepley   Level: developer
451b4937a87SMatthew G. Knepley 
452bcf0153eSBarry Smith   Note:
453bcf0153eSBarry Smith   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.
454bcf0153eSBarry Smith 
455*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrixLumped()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS`
456b4937a87SMatthew G. Knepley @*/
457d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCreateRHSMassMatrix(DM dm)
458d71ae5a4SJacob Faibussowitsch {
459b4937a87SMatthew G. Knepley   DMTS        tdm;
460b4937a87SMatthew G. Knepley   DMTS_Local *dmlocalts;
461b4937a87SMatthew G. Knepley   const char *prefix;
462b4937a87SMatthew G. Knepley 
463b4937a87SMatthew G. Knepley   PetscFunctionBegin;
464b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4659566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
4669566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
4679566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dm, dm, &dmlocalts->mass));
4689566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm), &dmlocalts->kspmass));
4699566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
4709566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(dmlocalts->kspmass, prefix));
4719566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(dmlocalts->kspmass, "mass_"));
4729566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(dmlocalts->kspmass));
4739566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(dmlocalts->kspmass, dmlocalts->mass, dmlocalts->mass));
4743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
475b4937a87SMatthew G. Knepley }
476b4937a87SMatthew G. Knepley 
477b4937a87SMatthew G. Knepley /*@C
478bcf0153eSBarry Smith   DMTSCreateRHSMassMatrixLumped - This creates the lumped mass matrix associated with the given `DM`, and a solver to invert it, and stores them in the `DM` context.
479b4937a87SMatthew G. Knepley 
480c3339decSBarry Smith   Collective
481b4937a87SMatthew G. Knepley 
4822fe279fdSBarry Smith   Input Parameter:
483bcf0153eSBarry Smith . dm   - `DM` providing the mass matrix
484b4937a87SMatthew G. Knepley 
485b4937a87SMatthew G. Knepley   Level: developer
486b4937a87SMatthew G. Knepley 
487bcf0153eSBarry Smith   Note:
488bcf0153eSBarry Smith   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.
489bcf0153eSBarry Smith   Since the matrix is lumped, inversion is trivial.
490bcf0153eSBarry Smith 
491*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrix()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS`
492b4937a87SMatthew G. Knepley @*/
493d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM dm)
494d71ae5a4SJacob Faibussowitsch {
495b4937a87SMatthew G. Knepley   DMTS        tdm;
496b4937a87SMatthew G. Knepley   DMTS_Local *dmlocalts;
497b4937a87SMatthew G. Knepley 
498b4937a87SMatthew G. Knepley   PetscFunctionBegin;
499b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5009566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
5019566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
5029566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrixLumped(dm, &dmlocalts->lumpedmassinv));
5039566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(dmlocalts->lumpedmassinv));
5049566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(dmlocalts->lumpedmassinv, NULL, "-lumped_mass_inv_view"));
5053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
506b4937a87SMatthew G. Knepley }
507b4937a87SMatthew G. Knepley 
508b4937a87SMatthew G. Knepley /*@C
509bcf0153eSBarry Smith   DMTSDestroyRHSMassMatrix - Destroys the mass matrix and solver stored in the `DM` context, if they exist.
510b4937a87SMatthew G. Knepley 
511b4937a87SMatthew G. Knepley   Logically Collective
512b4937a87SMatthew G. Knepley 
5132fe279fdSBarry Smith   Input Parameter:
514bcf0153eSBarry Smith . dm   - `DM` providing the mass matrix
515b4937a87SMatthew G. Knepley 
516b4937a87SMatthew G. Knepley   Level: developer
517b4937a87SMatthew G. Knepley 
518*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrixLumped()`, `DMCreateMassMatrix()`, `DMCreateMassMatrix()`, `DMTS`
519b4937a87SMatthew G. Knepley @*/
520d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSDestroyRHSMassMatrix(DM dm)
521d71ae5a4SJacob Faibussowitsch {
522b4937a87SMatthew G. Knepley   DMTS        tdm;
523b4937a87SMatthew G. Knepley   DMTS_Local *dmlocalts;
524b4937a87SMatthew G. Knepley 
525b4937a87SMatthew G. Knepley   PetscFunctionBegin;
526b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5279566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
5289566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
5299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&dmlocalts->lumpedmassinv));
5309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&dmlocalts->mass));
5319566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&dmlocalts->kspmass));
5323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
533b4937a87SMatthew G. Knepley }
534