xref: /petsc/src/ts/utils/dmlocalts.c (revision 0c9409e4145098594ef6b2e930d3138131b6944a)
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;
219566063dSJacob 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;
289566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tdm, (DMTS_Local **) &tdm->data));
299566063dSJacob 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) {
389566063dSJacob 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);
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));
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);
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));
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;
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) {
154d433e6cbSMatthew G. Knepley       case IS_COLORING_GLOBAL:
1559566063dSJacob 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       }
1599566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject) fdcoloring, ((PetscObject) dm)->prefix));
1609566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1619566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring));
1629566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject) dm, "DMDASNES_FDCOLORING", (PetscObject) fdcoloring));
1639566063dSJacob 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        */
1719566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject) dm));
172d433e6cbSMatthew G. Knepley     }
1739566063dSJacob 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) {
1779566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1789566063dSJacob 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);
2109566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
2119566063dSJacob 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
220*0c9409e4SMatthew G. Knepley   DMTSGetIFunctionLocal - get the local implicit function evaluation function. This function is called with local vector
221*0c9409e4SMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
222*0c9409e4SMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
223*0c9409e4SMatthew G. Knepley 
224*0c9409e4SMatthew G. Knepley   Logically Collective
225*0c9409e4SMatthew G. Knepley 
226*0c9409e4SMatthew G. Knepley   Input Parameter:
227*0c9409e4SMatthew G. Knepley . dm   - DM to associate callback with
228*0c9409e4SMatthew G. Knepley 
229*0c9409e4SMatthew G. Knepley   Output Parameters:
230*0c9409e4SMatthew G. Knepley + func - local function evaluation
231*0c9409e4SMatthew G. Knepley - ctx  - context for function evaluation
232*0c9409e4SMatthew G. Knepley 
233*0c9409e4SMatthew G. Knepley   Level: beginner
234*0c9409e4SMatthew G. Knepley 
235*0c9409e4SMatthew G. Knepley .seealso: DMTSSetIFunctionLocal((), DMTSSetIFunction(), DMTSSetIJacobianLocal()
236*0c9409e4SMatthew G. Knepley @*/
237*0c9409e4SMatthew G. Knepley PetscErrorCode DMTSGetIFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, Vec, void *), void **ctx)
238*0c9409e4SMatthew G. Knepley {
239*0c9409e4SMatthew G. Knepley   DMTS           tdm;
240*0c9409e4SMatthew G. Knepley   DMTS_Local    *dmlocalts;
241*0c9409e4SMatthew G. Knepley   PetscErrorCode ierr;
242*0c9409e4SMatthew G. Knepley 
243*0c9409e4SMatthew G. Knepley   PetscFunctionBegin;
244*0c9409e4SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
245*0c9409e4SMatthew G. Knepley   ierr = DMGetDMTS(dm, &tdm);CHKERRQ(ierr);
246*0c9409e4SMatthew G. Knepley   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr);
247*0c9409e4SMatthew G. Knepley   if (func) {PetscValidPointer(func, 2); *func = dmlocalts->ifunctionlocal;}
248*0c9409e4SMatthew G. Knepley   if (ctx)  {PetscValidPointer(ctx, 3);  *ctx  = dmlocalts->ifunctionlocalctx;}
249*0c9409e4SMatthew G. Knepley   PetscFunctionReturn(0);
250*0c9409e4SMatthew G. Knepley }
251*0c9409e4SMatthew G. Knepley 
252*0c9409e4SMatthew G. Knepley /*@C
253d433e6cbSMatthew G. Knepley   DMTSSetIFunctionLocal - set a local implicit function evaluation function. This function is called with local vector
254d433e6cbSMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
255d433e6cbSMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
256d433e6cbSMatthew G. Knepley 
257d433e6cbSMatthew G. Knepley   Logically Collective
258d433e6cbSMatthew G. Knepley 
2594165533cSJose E. Roman   Input Parameters:
260d433e6cbSMatthew G. Knepley + dm   - DM to associate callback with
261d433e6cbSMatthew G. Knepley . func - local function evaluation
262d433e6cbSMatthew G. Knepley - ctx  - context for function evaluation
263d433e6cbSMatthew G. Knepley 
264d433e6cbSMatthew G. Knepley   Level: beginner
265d433e6cbSMatthew G. Knepley 
266*0c9409e4SMatthew G. Knepley .seealso: DMTSGetIFunctionLocal(), DMTSSetIFunction(), DMTSSetIJacobianLocal()
267d433e6cbSMatthew G. Knepley @*/
268d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetIFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, Vec, void *), void *ctx)
269d433e6cbSMatthew G. Knepley {
270d433e6cbSMatthew G. Knepley   DMTS           tdm;
271d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts;
272d433e6cbSMatthew G. Knepley 
273d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
274d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2759566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
2769566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
277d433e6cbSMatthew G. Knepley 
278d433e6cbSMatthew G. Knepley   dmlocalts->ifunctionlocal    = func;
279d433e6cbSMatthew G. Knepley   dmlocalts->ifunctionlocalctx = ctx;
280d433e6cbSMatthew G. Knepley 
2819566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMLocal, dmlocalts));
282d433e6cbSMatthew G. Knepley   if (!tdm->ops->ijacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
2839566063dSJacob Faibussowitsch     PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts));
284d433e6cbSMatthew G. Knepley   }
285d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
286d433e6cbSMatthew G. Knepley }
287d433e6cbSMatthew G. Knepley 
288d433e6cbSMatthew G. Knepley /*@C
289*0c9409e4SMatthew G. Knepley   DMTSGetIJacobianLocal - get a local Jacobian evaluation function
290*0c9409e4SMatthew G. Knepley 
291*0c9409e4SMatthew G. Knepley   Logically Collective
292*0c9409e4SMatthew G. Knepley 
293*0c9409e4SMatthew G. Knepley   Input Parameter:
294*0c9409e4SMatthew G. Knepley . dm - DM to associate callback with
295*0c9409e4SMatthew G. Knepley 
296*0c9409e4SMatthew G. Knepley   Output Parameters:
297*0c9409e4SMatthew G. Knepley + func - local Jacobian evaluation
298*0c9409e4SMatthew G. Knepley - ctx - optional context for local Jacobian evaluation
299*0c9409e4SMatthew G. Knepley 
300*0c9409e4SMatthew G. Knepley   Level: beginner
301*0c9409e4SMatthew G. Knepley 
302*0c9409e4SMatthew G. Knepley .seealso: DMTSSetIJacobianLocal(), DMTSSetIFunctionLocal(), DMTSSetIJacobian(), DMTSSetIFunction()
303*0c9409e4SMatthew G. Knepley @*/
304*0c9409e4SMatthew G. Knepley PetscErrorCode DMTSGetIJacobianLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **ctx)
305*0c9409e4SMatthew G. Knepley {
306*0c9409e4SMatthew G. Knepley   DMTS           tdm;
307*0c9409e4SMatthew G. Knepley   DMTS_Local    *dmlocalts;
308*0c9409e4SMatthew G. Knepley   PetscErrorCode ierr;
309*0c9409e4SMatthew G. Knepley 
310*0c9409e4SMatthew G. Knepley   PetscFunctionBegin;
311*0c9409e4SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
312*0c9409e4SMatthew G. Knepley   ierr = DMGetDMTS(dm, &tdm);CHKERRQ(ierr);
313*0c9409e4SMatthew G. Knepley   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr);
314*0c9409e4SMatthew G. Knepley   if (func) {PetscValidPointer(func, 2); *func = dmlocalts->ijacobianlocal;}
315*0c9409e4SMatthew G. Knepley   if (ctx)  {PetscValidPointer(ctx, 3);  *ctx  = dmlocalts->ijacobianlocalctx;}
316*0c9409e4SMatthew G. Knepley   PetscFunctionReturn(0);
317*0c9409e4SMatthew G. Knepley }
318*0c9409e4SMatthew G. Knepley 
319*0c9409e4SMatthew G. Knepley /*@C
320d433e6cbSMatthew G. Knepley   DMTSSetIJacobianLocal - set a local Jacobian evaluation function
321d433e6cbSMatthew G. Knepley 
322d433e6cbSMatthew G. Knepley   Logically Collective
323d433e6cbSMatthew G. Knepley 
3244165533cSJose E. Roman   Input Parameters:
325d433e6cbSMatthew G. Knepley + dm - DM to associate callback with
326d433e6cbSMatthew G. Knepley . func - local Jacobian evaluation
327d433e6cbSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation
328d433e6cbSMatthew G. Knepley 
329d433e6cbSMatthew G. Knepley   Level: beginner
330d433e6cbSMatthew G. Knepley 
331*0c9409e4SMatthew G. Knepley .seealso: DMTSGetIJacobianLocal(), DMTSSetIFunctionLocal(), DMTSSetIJacobian(), DMTSSetIFunction()
332d433e6cbSMatthew G. Knepley @*/
333d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetIJacobianLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *ctx)
334d433e6cbSMatthew G. Knepley {
335d433e6cbSMatthew G. Knepley   DMTS           tdm;
336d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts;
337d433e6cbSMatthew G. Knepley 
338d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
339d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3409566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
3419566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
342d433e6cbSMatthew G. Knepley 
343d433e6cbSMatthew G. Knepley   dmlocalts->ijacobianlocal    = func;
344d433e6cbSMatthew G. Knepley   dmlocalts->ijacobianlocalctx = ctx;
345d433e6cbSMatthew G. Knepley 
3469566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts));
347d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
348d433e6cbSMatthew G. Knepley }
349d433e6cbSMatthew G. Knepley 
350d433e6cbSMatthew G. Knepley /*@C
351*0c9409e4SMatthew G. Knepley   DMTSGetRHSFunctionLocal - get a local rhs function evaluation function. This function is called with local vector
352*0c9409e4SMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
353*0c9409e4SMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
354*0c9409e4SMatthew G. Knepley 
355*0c9409e4SMatthew G. Knepley   Logically Collective
356*0c9409e4SMatthew G. Knepley 
357*0c9409e4SMatthew G. Knepley   Input Parameter:
358*0c9409e4SMatthew G. Knepley . dm   - DM to associate callback with
359*0c9409e4SMatthew G. Knepley 
360*0c9409e4SMatthew G. Knepley   Output Parameters:
361*0c9409e4SMatthew G. Knepley + func - local function evaluation
362*0c9409e4SMatthew G. Knepley - ctx  - context for function evaluation
363*0c9409e4SMatthew G. Knepley 
364*0c9409e4SMatthew G. Knepley   Level: beginner
365*0c9409e4SMatthew G. Knepley 
366*0c9409e4SMatthew G. Knepley .seealso: DMTSSetRHSFunctionLocal(), DMTSSetRHSFunction(), DMTSSetIFunction(), DMTSSetIJacobianLocal()
367*0c9409e4SMatthew G. Knepley @*/
368*0c9409e4SMatthew G. Knepley PetscErrorCode DMTSGetRHSFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, void *), void **ctx)
369*0c9409e4SMatthew G. Knepley {
370*0c9409e4SMatthew G. Knepley   DMTS           tdm;
371*0c9409e4SMatthew G. Knepley   DMTS_Local    *dmlocalts;
372*0c9409e4SMatthew G. Knepley   PetscErrorCode ierr;
373*0c9409e4SMatthew G. Knepley 
374*0c9409e4SMatthew G. Knepley   PetscFunctionBegin;
375*0c9409e4SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
376*0c9409e4SMatthew G. Knepley   ierr = DMGetDMTS(dm, &tdm);CHKERRQ(ierr);
377*0c9409e4SMatthew G. Knepley   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr);
378*0c9409e4SMatthew G. Knepley   if (func) {PetscValidPointer(func, 2); *func = dmlocalts->rhsfunctionlocal;}
379*0c9409e4SMatthew G. Knepley   if (ctx)  {PetscValidPointer(ctx, 3);  *ctx  = dmlocalts->rhsfunctionlocalctx;}
380*0c9409e4SMatthew G. Knepley   PetscFunctionReturn(0);
381*0c9409e4SMatthew G. Knepley }
382*0c9409e4SMatthew G. Knepley 
383*0c9409e4SMatthew G. Knepley /*@C
384d433e6cbSMatthew G. Knepley   DMTSSetRHSFunctionLocal - set a local rhs function evaluation function. This function is called with local vector
385d433e6cbSMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
386d433e6cbSMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
387d433e6cbSMatthew G. Knepley 
388d433e6cbSMatthew G. Knepley   Logically Collective
389d433e6cbSMatthew G. Knepley 
3904165533cSJose E. Roman   Input Parameters:
391d433e6cbSMatthew G. Knepley + dm   - DM to associate callback with
392d433e6cbSMatthew G. Knepley . func - local function evaluation
393d433e6cbSMatthew G. Knepley - ctx  - context for function evaluation
394d433e6cbSMatthew G. Knepley 
395d433e6cbSMatthew G. Knepley   Level: beginner
396d433e6cbSMatthew G. Knepley 
397*0c9409e4SMatthew G. Knepley .seealso: DMTSGetRHSFunctionLocal(), DMTSSetRHSFunction(), DMTSSetIFunction(), DMTSSetIJacobianLocal()
398d433e6cbSMatthew G. Knepley @*/
399d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx)
400d433e6cbSMatthew G. Knepley {
401d433e6cbSMatthew G. Knepley   DMTS           tdm;
402d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts;
403d433e6cbSMatthew G. Knepley 
404d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
405d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4069566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
4079566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
408d433e6cbSMatthew G. Knepley 
409d433e6cbSMatthew G. Knepley   dmlocalts->rhsfunctionlocal    = func;
410d433e6cbSMatthew G. Knepley   dmlocalts->rhsfunctionlocalctx = ctx;
411d433e6cbSMatthew G. Knepley 
4129566063dSJacob Faibussowitsch   PetscCall(DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMLocal, dmlocalts));
413d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
414d433e6cbSMatthew G. Knepley }
415b4937a87SMatthew G. Knepley 
416b4937a87SMatthew G. Knepley /*@C
417b4937a87SMatthew 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.
418b4937a87SMatthew G. Knepley 
419b4937a87SMatthew G. Knepley   Collective on dm
420b4937a87SMatthew G. Knepley 
421b4937a87SMatthew G. Knepley   Input Parameters:
422b4937a87SMatthew G. Knepley . dm   - DM providing the mass matrix
423b4937a87SMatthew G. Knepley 
424b4937a87SMatthew 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.
425b4937a87SMatthew G. Knepley 
426b4937a87SMatthew G. Knepley   Level: developer
427b4937a87SMatthew G. Knepley 
428b4937a87SMatthew G. Knepley .seealso: DMTSCreateRHSMassMatrixLumped(), DMTSDestroyRHSMassMatrix(), DMCreateMassMatrix(), DMTS
429b4937a87SMatthew G. Knepley @*/
430b4937a87SMatthew G. Knepley PetscErrorCode DMTSCreateRHSMassMatrix(DM dm)
431b4937a87SMatthew G. Knepley {
432b4937a87SMatthew G. Knepley   DMTS           tdm;
433b4937a87SMatthew G. Knepley   DMTS_Local    *dmlocalts;
434b4937a87SMatthew G. Knepley   const char    *prefix;
435b4937a87SMatthew G. Knepley 
436b4937a87SMatthew G. Knepley   PetscFunctionBegin;
437b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4389566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
4399566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
4409566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dm, dm, &dmlocalts->mass));
4419566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PetscObjectComm((PetscObject) dm), &dmlocalts->kspmass));
4429566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix));
4439566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(dmlocalts->kspmass, prefix));
4449566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(dmlocalts->kspmass, "mass_"));
4459566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(dmlocalts->kspmass));
4469566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(dmlocalts->kspmass, dmlocalts->mass, dmlocalts->mass));
447b4937a87SMatthew G. Knepley   PetscFunctionReturn(0);
448b4937a87SMatthew G. Knepley }
449b4937a87SMatthew G. Knepley 
450b4937a87SMatthew G. Knepley /*@C
451b4937a87SMatthew 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.
452b4937a87SMatthew G. Knepley 
453b4937a87SMatthew G. Knepley   Collective on dm
454b4937a87SMatthew G. Knepley 
455b4937a87SMatthew G. Knepley   Input Parameters:
456b4937a87SMatthew G. Knepley . dm   - DM providing the mass matrix
457b4937a87SMatthew G. Knepley 
458b4937a87SMatthew 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.
459b4937a87SMatthew G. Knepley   Since the matrix is lumped, inversion is trivial.
460b4937a87SMatthew G. Knepley 
461b4937a87SMatthew G. Knepley   Level: developer
462b4937a87SMatthew G. Knepley 
463b4937a87SMatthew G. Knepley .seealso: DMTSCreateRHSMassMatrix(), DMTSDestroyRHSMassMatrix(), DMCreateMassMatrix(), DMTS
464b4937a87SMatthew G. Knepley @*/
465b4937a87SMatthew G. Knepley PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM dm)
466b4937a87SMatthew G. Knepley {
467b4937a87SMatthew G. Knepley   DMTS           tdm;
468b4937a87SMatthew G. Knepley   DMTS_Local    *dmlocalts;
469b4937a87SMatthew G. Knepley 
470b4937a87SMatthew G. Knepley   PetscFunctionBegin;
471b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4729566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
4739566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
4749566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrixLumped(dm, &dmlocalts->lumpedmassinv));
4759566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(dmlocalts->lumpedmassinv));
4769566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(dmlocalts->lumpedmassinv, NULL, "-lumped_mass_inv_view"));
477b4937a87SMatthew G. Knepley   PetscFunctionReturn(0);
478b4937a87SMatthew G. Knepley }
479b4937a87SMatthew G. Knepley 
480b4937a87SMatthew G. Knepley /*@C
481b4937a87SMatthew G. Knepley   DMTSDestroyRHSMassMatrix - Destroys the mass matrix and solver stored in the DMTS context, if they exist.
482b4937a87SMatthew G. Knepley 
483b4937a87SMatthew G. Knepley   Logically Collective
484b4937a87SMatthew G. Knepley 
485b4937a87SMatthew G. Knepley   Input Parameters:
486b4937a87SMatthew G. Knepley . dm   - DM providing the mass matrix
487b4937a87SMatthew G. Knepley 
488b4937a87SMatthew G. Knepley   Level: developer
489b4937a87SMatthew G. Knepley 
490b4937a87SMatthew G. Knepley .seealso: DMTSCreateRHSMassMatrixLumped(), DMCreateMassMatrix(), DMCreateMassMatrix(), DMTS
491b4937a87SMatthew G. Knepley @*/
492b4937a87SMatthew G. Knepley PetscErrorCode DMTSDestroyRHSMassMatrix(DM dm)
493b4937a87SMatthew G. Knepley {
494b4937a87SMatthew G. Knepley   DMTS           tdm;
495b4937a87SMatthew G. Knepley   DMTS_Local    *dmlocalts;
496b4937a87SMatthew G. Knepley 
497b4937a87SMatthew G. Knepley   PetscFunctionBegin;
498b4937a87SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4999566063dSJacob Faibussowitsch   PetscCall(DMGetDMTSWrite(dm, &tdm));
5009566063dSJacob Faibussowitsch   PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts));
5019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&dmlocalts->lumpedmassinv));
5029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&dmlocalts->mass));
5039566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&dmlocalts->kspmass));
504b4937a87SMatthew G. Knepley   PetscFunctionReturn(0);
505b4937a87SMatthew G. Knepley }
506