xref: /petsc/src/ts/utils/dmlocalts.c (revision 9566063d113dddea24716c546802770db7481bc0)
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 
18d433e6cbSMatthew G. Knepley static PetscErrorCode DMTSDestroy_DMLocal(DMTS tdm)
19d433e6cbSMatthew G. Knepley {
20d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
21*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(tdm->data));
22d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
23d433e6cbSMatthew G. Knepley }
24d433e6cbSMatthew G. Knepley 
25d433e6cbSMatthew G. Knepley static PetscErrorCode DMTSDuplicate_DMLocal(DMTS oldtdm, DMTS tdm)
26d433e6cbSMatthew G. Knepley {
27d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
28*9566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tdm, (DMTS_Local **) &tdm->data));
29*9566063dSJacob Faibussowitsch   if (oldtdm->data) PetscCall(PetscMemcpy(tdm->data, oldtdm->data, sizeof(DMTS_Local)));
30d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
31d433e6cbSMatthew G. Knepley }
32d433e6cbSMatthew G. Knepley 
33d433e6cbSMatthew G. Knepley static PetscErrorCode DMLocalTSGetContext(DM dm, DMTS tdm, DMTS_Local **dmlocalts)
34d433e6cbSMatthew G. Knepley {
35d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
36d433e6cbSMatthew G. Knepley   *dmlocalts = NULL;
37d433e6cbSMatthew G. Knepley   if (!tdm->data) {
38*9566063dSJacob Faibussowitsch     PetscCall(PetscNewLog(dm, (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;
44d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
45d433e6cbSMatthew G. Knepley }
46d433e6cbSMatthew G. Knepley 
47d433e6cbSMatthew G. Knepley static PetscErrorCode TSComputeIFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, Vec F, void *ctx)
48d433e6cbSMatthew G. Knepley {
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);
58*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
59*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX));
60*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX_t));
61*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locF));
62*9566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locX));
63*9566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locX_t));
64*9566063dSJacob Faibussowitsch   if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, locX_t,dmlocalts->boundarylocalctx));
65*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
66*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
67*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t));
68*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t));
69*9566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locF));
70d433e6cbSMatthew G. Knepley   CHKMEMQ;
71*9566063dSJacob Faibussowitsch   PetscCall((*dmlocalts->ifunctionlocal)(dm, time, locX, locX_t, locF, dmlocalts->ifunctionlocalctx));
72d433e6cbSMatthew G. Knepley   CHKMEMQ;
73*9566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
74*9566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F));
75*9566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F));
76*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX));
77*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX_t));
78*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locF));
79d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
80d433e6cbSMatthew G. Knepley }
81d433e6cbSMatthew G. Knepley 
82d433e6cbSMatthew G. Knepley static PetscErrorCode TSComputeRHSFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec F, void *ctx)
83d433e6cbSMatthew G. Knepley {
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);
92*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
93*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locX));
94*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &locF));
95*9566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locX));
96*9566063dSJacob Faibussowitsch   if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm,time,locX,NULL,dmlocalts->boundarylocalctx));
97*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
98*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
99*9566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locF));
100b4937a87SMatthew G. Knepley   CHKMEMQ;
101*9566063dSJacob Faibussowitsch   PetscCall((*dmlocalts->rhsfunctionlocal)(dm, time, locX, locF, dmlocalts->rhsfunctionlocalctx));
102b4937a87SMatthew G. Knepley   CHKMEMQ;
103*9566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
104*9566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F));
105*9566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F));
106b4937a87SMatthew G. Knepley   if (dmlocalts->lumpedmassinv) {
107*9566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(F, dmlocalts->lumpedmassinv, F));
108b4937a87SMatthew G. Knepley   } else if (dmlocalts->kspmass) {
109b4937a87SMatthew G. Knepley     Vec tmp;
110b4937a87SMatthew G. Knepley 
111*9566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &tmp));
112*9566063dSJacob Faibussowitsch     PetscCall(KSPSolve(dmlocalts->kspmass, F, tmp));
113*9566063dSJacob Faibussowitsch     PetscCall(VecCopy(tmp, F));
114*9566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &tmp));
115b4937a87SMatthew G. Knepley   }
116*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locX));
117*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &locF));
118d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
119d433e6cbSMatthew G. Knepley }
120d433e6cbSMatthew G. Knepley 
121d433e6cbSMatthew G. Knepley static PetscErrorCode TSComputeIJacobian_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, PetscReal a, Mat A, Mat B, void *ctx)
122d433e6cbSMatthew G. Knepley {
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;
128*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
129d433e6cbSMatthew G. Knepley   if (dmlocalts->ijacobianlocal) {
130*9566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &locX));
131*9566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &locX_t));
132*9566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(locX));
133*9566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(locX_t));
134*9566063dSJacob Faibussowitsch     if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm,time,locX,locX_t,dmlocalts->boundarylocalctx));
135*9566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
136*9566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
137*9566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t));
138*9566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t));
139d433e6cbSMatthew G. Knepley     CHKMEMQ;
140*9566063dSJacob Faibussowitsch     PetscCall((*dmlocalts->ijacobianlocal)(dm, time, locX, locX_t, a, A, B, dmlocalts->ijacobianlocalctx));
141d433e6cbSMatthew G. Knepley     CHKMEMQ;
142*9566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &locX));
143*9566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &locX_t));
144d433e6cbSMatthew G. Knepley   } else {
145d433e6cbSMatthew G. Knepley     MatFDColoring fdcoloring;
146*9566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject) dm, "DMDASNES_FDCOLORING", (PetscObject *) &fdcoloring));
147d433e6cbSMatthew G. Knepley     if (!fdcoloring) {
148d433e6cbSMatthew G. Knepley       ISColoring coloring;
149d433e6cbSMatthew G. Knepley 
150*9566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring));
151*9566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring));
152*9566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&coloring));
153d433e6cbSMatthew G. Knepley       switch (dm->coloringtype) {
154d433e6cbSMatthew G. Knepley       case IS_COLORING_GLOBAL:
155*9566063dSJacob Faibussowitsch         PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode (*)(void)) TSComputeIFunction_DMLocal, dmlocalts));
156d433e6cbSMatthew G. Knepley         break;
15798921bdaSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject) ts), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]);
158d433e6cbSMatthew G. Knepley       }
159*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject) fdcoloring, ((PetscObject) dm)->prefix));
160*9566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
161*9566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring));
162*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject) dm, "DMDASNES_FDCOLORING", (PetscObject) fdcoloring));
163*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject) fdcoloring));
164d433e6cbSMatthew G. Knepley 
165d433e6cbSMatthew G. Knepley       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
166d433e6cbSMatthew G. Knepley        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
167d433e6cbSMatthew G. Knepley        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
168d433e6cbSMatthew G. Knepley        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
169d433e6cbSMatthew G. Knepley        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
170d433e6cbSMatthew G. Knepley        */
171*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject) dm));
172d433e6cbSMatthew G. Knepley     }
173*9566063dSJacob Faibussowitsch     PetscCall(MatFDColoringApply(B, fdcoloring, X, ts));
174d433e6cbSMatthew G. Knepley   }
175d433e6cbSMatthew G. Knepley   /* This will be redundant if the user called both, but it's too common to forget. */
176d433e6cbSMatthew G. Knepley   if (A != B) {
177*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
178*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
179d433e6cbSMatthew G. Knepley   }
180d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
181d433e6cbSMatthew G. Knepley }
182d433e6cbSMatthew G. Knepley 
1837eb36faeSMatthew G. Knepley /*@C
184f2b97b97SMatthew G. Knepley   DMTSSetBoundaryLocal - set the function for essential boundary data for a local implicit function evaluation.
185a40652d4SToby 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).
186a40652d4SToby Isaac     Vectors are initialized to zero before this function, so it is only needed for non homogeneous data.
187a40652d4SToby Isaac 
188a40652d4SToby Isaac   Note that this function is somewhat optional: boundary data could potentially be inserted by a function passed to
189a40652d4SToby Isaac   DMTSSetIFunctionLocal().  The use case for this function is for discretizations with constraints (see
190a40652d4SToby Isaac   DMGetDefaultConstraints()): this function inserts boundary values before constraint interpolation.
1917eb36faeSMatthew G. Knepley 
1927eb36faeSMatthew G. Knepley   Logically Collective
1937eb36faeSMatthew G. Knepley 
1944165533cSJose E. Roman   Input Parameters:
1957eb36faeSMatthew G. Knepley + dm   - DM to associate callback with
196a40652d4SToby Isaac . func - local function evaluation
197a40652d4SToby Isaac - ctx  - context for function evaluation
1987eb36faeSMatthew G. Knepley 
1997eb36faeSMatthew G. Knepley   Level: intermediate
2007eb36faeSMatthew G. Knepley 
201a40652d4SToby Isaac .seealso: DMTSSetIFunction(), DMTSSetIJacobianLocal()
2027eb36faeSMatthew G. Knepley @*/
203f2b97b97SMatthew G. Knepley PetscErrorCode DMTSSetBoundaryLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx)
2047eb36faeSMatthew G. Knepley {
2057eb36faeSMatthew G. Knepley   DMTS           tdm;
2067eb36faeSMatthew G. Knepley   DMTS_Local    *dmlocalts;
2077eb36faeSMatthew G. Knepley 
2087eb36faeSMatthew G. Knepley   PetscFunctionBegin;
2097eb36faeSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
210*9566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
211*9566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
2127eb36faeSMatthew G. Knepley 
2137eb36faeSMatthew G. Knepley   dmlocalts->boundarylocal    = func;
2147eb36faeSMatthew G. Knepley   dmlocalts->boundarylocalctx = ctx;
215a40652d4SToby Isaac 
2167eb36faeSMatthew G. Knepley   PetscFunctionReturn(0);
2177eb36faeSMatthew G. Knepley }
2187eb36faeSMatthew G. Knepley 
219d433e6cbSMatthew G. Knepley /*@C
220d433e6cbSMatthew G. Knepley   DMTSSetIFunctionLocal - set a local implicit function evaluation function. This function is called with local vector
221d433e6cbSMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
222d433e6cbSMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
223d433e6cbSMatthew G. Knepley 
224d433e6cbSMatthew G. Knepley   Logically Collective
225d433e6cbSMatthew G. Knepley 
2264165533cSJose E. Roman   Input Parameters:
227d433e6cbSMatthew G. Knepley + dm   - DM to associate callback with
228d433e6cbSMatthew G. Knepley . func - local function evaluation
229d433e6cbSMatthew G. Knepley - ctx  - context for function evaluation
230d433e6cbSMatthew G. Knepley 
231d433e6cbSMatthew G. Knepley   Level: beginner
232d433e6cbSMatthew G. Knepley 
233d433e6cbSMatthew G. Knepley .seealso: DMTSSetIFunction(), DMTSSetIJacobianLocal()
234d433e6cbSMatthew G. Knepley @*/
235d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetIFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, Vec, void *), void *ctx)
236d433e6cbSMatthew G. Knepley {
237d433e6cbSMatthew G. Knepley   DMTS           tdm;
238d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts;
239d433e6cbSMatthew G. Knepley 
240d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
241d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
242*9566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
243*9566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
244d433e6cbSMatthew G. Knepley 
245d433e6cbSMatthew G. Knepley   dmlocalts->ifunctionlocal    = func;
246d433e6cbSMatthew G. Knepley   dmlocalts->ifunctionlocalctx = ctx;
247d433e6cbSMatthew G. Knepley 
248*9566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMLocal, dmlocalts));
249d433e6cbSMatthew G. Knepley   if (!tdm->ops->ijacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
250*9566063dSJacob Faibussowitsch     PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts));
251d433e6cbSMatthew G. Knepley   }
252d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
253d433e6cbSMatthew G. Knepley }
254d433e6cbSMatthew G. Knepley 
255d433e6cbSMatthew G. Knepley /*@C
256d433e6cbSMatthew G. Knepley   DMTSSetIJacobianLocal - set a local Jacobian evaluation function
257d433e6cbSMatthew G. Knepley 
258d433e6cbSMatthew G. Knepley   Logically Collective
259d433e6cbSMatthew G. Knepley 
2604165533cSJose E. Roman   Input Parameters:
261d433e6cbSMatthew G. Knepley + dm - DM to associate callback with
262d433e6cbSMatthew G. Knepley . func - local Jacobian evaluation
263d433e6cbSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation
264d433e6cbSMatthew G. Knepley 
265d433e6cbSMatthew G. Knepley   Level: beginner
266d433e6cbSMatthew G. Knepley 
267d433e6cbSMatthew G. Knepley .seealso: DMTSSetIFunctionLocal(), DMTSSetIJacobian(), DMTSSetIFunction()
268d433e6cbSMatthew G. Knepley @*/
269d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetIJacobianLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *ctx)
270d433e6cbSMatthew G. Knepley {
271d433e6cbSMatthew G. Knepley   DMTS           tdm;
272d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts;
273d433e6cbSMatthew G. Knepley 
274d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
275d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
276*9566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
277*9566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
278d433e6cbSMatthew G. Knepley 
279d433e6cbSMatthew G. Knepley   dmlocalts->ijacobianlocal    = func;
280d433e6cbSMatthew G. Knepley   dmlocalts->ijacobianlocalctx = ctx;
281d433e6cbSMatthew G. Knepley 
282*9566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts));
283d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
284d433e6cbSMatthew G. Knepley }
285d433e6cbSMatthew G. Knepley 
286d433e6cbSMatthew G. Knepley /*@C
287d433e6cbSMatthew G. Knepley   DMTSSetRHSFunctionLocal - set a local rhs function evaluation function. This function is called with local vector
288d433e6cbSMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
289d433e6cbSMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
290d433e6cbSMatthew G. Knepley 
291d433e6cbSMatthew G. Knepley   Logically Collective
292d433e6cbSMatthew G. Knepley 
2934165533cSJose E. Roman   Input Parameters:
294d433e6cbSMatthew G. Knepley + dm   - DM to associate callback with
295d433e6cbSMatthew G. Knepley . func - local function evaluation
296d433e6cbSMatthew G. Knepley - ctx  - context for function evaluation
297d433e6cbSMatthew G. Knepley 
298d433e6cbSMatthew G. Knepley   Level: beginner
299d433e6cbSMatthew G. Knepley 
300d433e6cbSMatthew G. Knepley .seealso: DMTSSetRHSFunction(), DMTSSetIFunction(), DMTSSetIJacobianLocal()
301d433e6cbSMatthew G. Knepley @*/
302d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx)
303d433e6cbSMatthew G. Knepley {
304d433e6cbSMatthew G. Knepley   DMTS           tdm;
305d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts;
306d433e6cbSMatthew G. Knepley 
307d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
308d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
309*9566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
310*9566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
311d433e6cbSMatthew G. Knepley 
312d433e6cbSMatthew G. Knepley   dmlocalts->rhsfunctionlocal    = func;
313d433e6cbSMatthew G. Knepley   dmlocalts->rhsfunctionlocalctx = ctx;
314d433e6cbSMatthew G. Knepley 
315*9566063dSJacob Faibussowitsch   PetscCall(DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMLocal, dmlocalts));
316d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
317d433e6cbSMatthew G. Knepley }
318b4937a87SMatthew G. Knepley 
319b4937a87SMatthew G. Knepley /*@C
320b4937a87SMatthew 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.
321b4937a87SMatthew G. Knepley 
322b4937a87SMatthew G. Knepley   Collective on dm
323b4937a87SMatthew G. Knepley 
324b4937a87SMatthew G. Knepley   Input Parameters:
325b4937a87SMatthew G. Knepley . dm   - DM providing the mass matrix
326b4937a87SMatthew G. Knepley 
327b4937a87SMatthew 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.
328b4937a87SMatthew G. Knepley 
329b4937a87SMatthew G. Knepley   Level: developer
330b4937a87SMatthew G. Knepley 
331b4937a87SMatthew G. Knepley .seealso: DMTSCreateRHSMassMatrixLumped(), DMTSDestroyRHSMassMatrix(), DMCreateMassMatrix(), DMTS
332b4937a87SMatthew G. Knepley @*/
333b4937a87SMatthew G. Knepley PetscErrorCode DMTSCreateRHSMassMatrix(DM dm)
334b4937a87SMatthew G. Knepley {
335b4937a87SMatthew G. Knepley   DMTS           tdm;
336b4937a87SMatthew G. Knepley   DMTS_Local    *dmlocalts;
337b4937a87SMatthew G. Knepley   const char    *prefix;
338b4937a87SMatthew G. Knepley 
339b4937a87SMatthew G. Knepley   PetscFunctionBegin;
340b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
341*9566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
342*9566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
343*9566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dm, dm, &dmlocalts->mass));
344*9566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PetscObjectComm((PetscObject) dm), &dmlocalts->kspmass));
345*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix));
346*9566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(dmlocalts->kspmass, prefix));
347*9566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(dmlocalts->kspmass, "mass_"));
348*9566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(dmlocalts->kspmass));
349*9566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(dmlocalts->kspmass, dmlocalts->mass, dmlocalts->mass));
350b4937a87SMatthew G. Knepley   PetscFunctionReturn(0);
351b4937a87SMatthew G. Knepley }
352b4937a87SMatthew G. Knepley 
353b4937a87SMatthew G. Knepley /*@C
354b4937a87SMatthew 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.
355b4937a87SMatthew G. Knepley 
356b4937a87SMatthew G. Knepley   Collective on dm
357b4937a87SMatthew G. Knepley 
358b4937a87SMatthew G. Knepley   Input Parameters:
359b4937a87SMatthew G. Knepley . dm   - DM providing the mass matrix
360b4937a87SMatthew G. Knepley 
361b4937a87SMatthew 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.
362b4937a87SMatthew G. Knepley   Since the matrix is lumped, inversion is trivial.
363b4937a87SMatthew G. Knepley 
364b4937a87SMatthew G. Knepley   Level: developer
365b4937a87SMatthew G. Knepley 
366b4937a87SMatthew G. Knepley .seealso: DMTSCreateRHSMassMatrix(), DMTSDestroyRHSMassMatrix(), DMCreateMassMatrix(), DMTS
367b4937a87SMatthew G. Knepley @*/
368b4937a87SMatthew G. Knepley PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM dm)
369b4937a87SMatthew G. Knepley {
370b4937a87SMatthew G. Knepley   DMTS           tdm;
371b4937a87SMatthew G. Knepley   DMTS_Local    *dmlocalts;
372b4937a87SMatthew G. Knepley 
373b4937a87SMatthew G. Knepley   PetscFunctionBegin;
374b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
375*9566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
376*9566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
377*9566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrixLumped(dm, &dmlocalts->lumpedmassinv));
378*9566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(dmlocalts->lumpedmassinv));
379*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(dmlocalts->lumpedmassinv, NULL, "-lumped_mass_inv_view"));
380b4937a87SMatthew G. Knepley   PetscFunctionReturn(0);
381b4937a87SMatthew G. Knepley }
382b4937a87SMatthew G. Knepley 
383b4937a87SMatthew G. Knepley /*@C
384b4937a87SMatthew G. Knepley   DMTSDestroyRHSMassMatrix - Destroys the mass matrix and solver stored in the DMTS context, if they exist.
385b4937a87SMatthew G. Knepley 
386b4937a87SMatthew G. Knepley   Logically Collective
387b4937a87SMatthew G. Knepley 
388b4937a87SMatthew G. Knepley   Input Parameters:
389b4937a87SMatthew G. Knepley . dm   - DM providing the mass matrix
390b4937a87SMatthew G. Knepley 
391b4937a87SMatthew G. Knepley   Level: developer
392b4937a87SMatthew G. Knepley 
393b4937a87SMatthew G. Knepley .seealso: DMTSCreateRHSMassMatrixLumped(), DMCreateMassMatrix(), DMCreateMassMatrix(), DMTS
394b4937a87SMatthew G. Knepley @*/
395b4937a87SMatthew G. Knepley PetscErrorCode DMTSDestroyRHSMassMatrix(DM dm)
396b4937a87SMatthew G. Knepley {
397b4937a87SMatthew G. Knepley   DMTS           tdm;
398b4937a87SMatthew G. Knepley   DMTS_Local    *dmlocalts;
399b4937a87SMatthew G. Knepley 
400b4937a87SMatthew G. Knepley   PetscFunctionBegin;
401b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
402*9566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
403*9566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
404*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&dmlocalts->lumpedmassinv));
405*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&dmlocalts->mass));
406*9566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&dmlocalts->kspmass));
407b4937a87SMatthew G. Knepley   PetscFunctionReturn(0);
408b4937a87SMatthew G. Knepley }
409