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