xref: /petsc/src/ts/utils/dmlocalts.c (revision a40652d484231dadb820fc3a9703f41e31810d43)
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 {
5*a40652d4SToby Isaac   PetscErrorCode (*iboundarylocal)(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*);
8*a40652d4SToby Isaac   PetscErrorCode (*rhsboundarylocal)(DM,PetscReal,Vec,void*);
9d433e6cbSMatthew G. Knepley   PetscErrorCode (*rhsfunctionlocal)(DM,PetscReal,Vec,Vec,void*);
10*a40652d4SToby Isaac   void *iboundarylocalctx;
11d433e6cbSMatthew G. Knepley   void *ifunctionlocalctx;
12d433e6cbSMatthew G. Knepley   void *ijacobianlocalctx;
13*a40652d4SToby Isaac   void *rhsboundarylocalctx;
14d433e6cbSMatthew G. Knepley   void *rhsfunctionlocalctx;
15d433e6cbSMatthew G. Knepley } DMTS_Local;
16d433e6cbSMatthew G. Knepley 
17d433e6cbSMatthew G. Knepley #undef __FUNCT__
18d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMTSDestroy_DMLocal"
19d433e6cbSMatthew G. Knepley static PetscErrorCode DMTSDestroy_DMLocal(DMTS tdm)
20d433e6cbSMatthew G. Knepley {
21d433e6cbSMatthew G. Knepley   PetscErrorCode ierr;
22d433e6cbSMatthew G. Knepley 
23d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
24d433e6cbSMatthew G. Knepley   ierr = PetscFree(tdm->data);CHKERRQ(ierr);
25d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
26d433e6cbSMatthew G. Knepley }
27d433e6cbSMatthew G. Knepley 
28d433e6cbSMatthew G. Knepley #undef __FUNCT__
29d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMTSDuplicate_DMLocal"
30d433e6cbSMatthew G. Knepley static PetscErrorCode DMTSDuplicate_DMLocal(DMTS oldtdm, DMTS tdm)
31d433e6cbSMatthew G. Knepley {
32d433e6cbSMatthew G. Knepley   PetscErrorCode ierr;
33d433e6cbSMatthew G. Knepley 
34d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
35d433e6cbSMatthew G. Knepley   ierr = PetscNewLog(tdm, (DMTS_Local **) &tdm->data);CHKERRQ(ierr);
36d433e6cbSMatthew G. Knepley   if (oldtdm->data) {ierr = PetscMemcpy(tdm->data, oldtdm->data, sizeof(DMTS_Local));CHKERRQ(ierr);}
37d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
38d433e6cbSMatthew G. Knepley }
39d433e6cbSMatthew G. Knepley 
40d433e6cbSMatthew G. Knepley #undef __FUNCT__
41d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMLocalTSGetContext"
42d433e6cbSMatthew G. Knepley static PetscErrorCode DMLocalTSGetContext(DM dm, DMTS tdm, DMTS_Local **dmlocalts)
43d433e6cbSMatthew G. Knepley {
44d433e6cbSMatthew G. Knepley   PetscErrorCode ierr;
45d433e6cbSMatthew G. Knepley 
46d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
47d433e6cbSMatthew G. Knepley   *dmlocalts = NULL;
48d433e6cbSMatthew G. Knepley   if (!tdm->data) {
49d433e6cbSMatthew G. Knepley     ierr = PetscNewLog(dm, (DMTS_Local **) &tdm->data);CHKERRQ(ierr);
50d433e6cbSMatthew G. Knepley 
51d433e6cbSMatthew G. Knepley     tdm->ops->destroy   = DMTSDestroy_DMLocal;
52d433e6cbSMatthew G. Knepley     tdm->ops->duplicate = DMTSDuplicate_DMLocal;
53d433e6cbSMatthew G. Knepley   }
54d433e6cbSMatthew G. Knepley   *dmlocalts = (DMTS_Local *) tdm->data;
55d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
56d433e6cbSMatthew G. Knepley }
57d433e6cbSMatthew G. Knepley 
58d433e6cbSMatthew G. Knepley #undef __FUNCT__
59d433e6cbSMatthew G. Knepley #define __FUNCT__ "TSComputeIFunction_DMLocal"
60d433e6cbSMatthew G. Knepley static PetscErrorCode TSComputeIFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, Vec F, void *ctx)
61d433e6cbSMatthew G. Knepley {
62d433e6cbSMatthew G. Knepley   DM             dm;
63d433e6cbSMatthew G. Knepley   Vec            locX, locX_t, locF;
64d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts = (DMTS_Local *) ctx;
65d433e6cbSMatthew G. Knepley   PetscErrorCode ierr;
66d433e6cbSMatthew G. Knepley 
67d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
68d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
69d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
70d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(X_t,VEC_CLASSID,4);
71d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(F,VEC_CLASSID,5);
72d433e6cbSMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
73d433e6cbSMatthew G. Knepley   ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr);
74d433e6cbSMatthew G. Knepley   ierr = DMGetLocalVector(dm, &locX_t);CHKERRQ(ierr);
75d433e6cbSMatthew G. Knepley   ierr = DMGetLocalVector(dm, &locF);CHKERRQ(ierr);
76d433e6cbSMatthew G. Knepley   ierr = VecZeroEntries(locX);CHKERRQ(ierr);
77d433e6cbSMatthew G. Knepley   ierr = VecZeroEntries(locX_t);CHKERRQ(ierr);
78*a40652d4SToby Isaac   if (dmlocalts->iboundarylocal) {ierr = (*dmlocalts->iboundarylocal)(dm, time, locX, locX_t,dmlocalts->iboundarylocalctx);CHKERRQ(ierr);}
79d433e6cbSMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
80d433e6cbSMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
81d433e6cbSMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t);CHKERRQ(ierr);
82d433e6cbSMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t);CHKERRQ(ierr);
83d433e6cbSMatthew G. Knepley   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
84d433e6cbSMatthew G. Knepley   CHKMEMQ;
85d433e6cbSMatthew G. Knepley   ierr = (*dmlocalts->ifunctionlocal)(dm, time, locX, locX_t, locF, dmlocalts->ifunctionlocalctx);CHKERRQ(ierr);
86d433e6cbSMatthew G. Knepley   CHKMEMQ;
87d433e6cbSMatthew G. Knepley   ierr = VecZeroEntries(F);CHKERRQ(ierr);
88d433e6cbSMatthew G. Knepley   ierr = DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F);CHKERRQ(ierr);
89d433e6cbSMatthew G. Knepley   ierr = DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F);CHKERRQ(ierr);
90d433e6cbSMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr);
91d433e6cbSMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &locX_t);CHKERRQ(ierr);
92d433e6cbSMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &locF);CHKERRQ(ierr);
93d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
94d433e6cbSMatthew G. Knepley }
95d433e6cbSMatthew G. Knepley 
96d433e6cbSMatthew G. Knepley #undef __FUNCT__
97d433e6cbSMatthew G. Knepley #define __FUNCT__ "TSComputeRHSFunction_DMLocal"
98d433e6cbSMatthew G. Knepley static PetscErrorCode TSComputeRHSFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec F, void *ctx)
99d433e6cbSMatthew G. Knepley {
100d433e6cbSMatthew G. Knepley   DM             dm;
101d433e6cbSMatthew G. Knepley   Vec            locX;
102d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts = (DMTS_Local *) ctx;
103d433e6cbSMatthew G. Knepley   PetscErrorCode ierr;
104d433e6cbSMatthew G. Knepley 
105d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
106d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
107d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
108d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(F,VEC_CLASSID,5);
109d433e6cbSMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
110d433e6cbSMatthew G. Knepley   ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr);
111d433e6cbSMatthew G. Knepley   ierr = VecZeroEntries(locX);CHKERRQ(ierr);
112*a40652d4SToby Isaac   if (dmlocalts->rhsboundarylocal) {ierr = (*dmlocalts->rhsboundarylocal)(dm, time, locX, dmlocalts->rhsboundarylocalctx);CHKERRQ(ierr);}
113d433e6cbSMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
114d433e6cbSMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
115d433e6cbSMatthew G. Knepley   ierr = VecZeroEntries(F);CHKERRQ(ierr);
116d433e6cbSMatthew G. Knepley   CHKMEMQ;
117d433e6cbSMatthew G. Knepley   ierr = (*dmlocalts->rhsfunctionlocal)(dm, time, locX, F, dmlocalts->rhsfunctionlocalctx);CHKERRQ(ierr);
118d433e6cbSMatthew G. Knepley   CHKMEMQ;
119d433e6cbSMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr);
120d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
121d433e6cbSMatthew G. Knepley }
122d433e6cbSMatthew G. Knepley 
123d433e6cbSMatthew G. Knepley #undef __FUNCT__
124d433e6cbSMatthew G. Knepley #define __FUNCT__ "TSComputeIJacobian_DMLocal"
125d433e6cbSMatthew G. Knepley static PetscErrorCode TSComputeIJacobian_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, PetscReal a, Mat A, Mat B, void *ctx)
126d433e6cbSMatthew G. Knepley {
127d433e6cbSMatthew G. Knepley   DM             dm;
128d433e6cbSMatthew G. Knepley   Vec            locX, locX_t;
129d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts = (DMTS_Local *) ctx;
130d433e6cbSMatthew G. Knepley   PetscErrorCode ierr;
131d433e6cbSMatthew G. Knepley 
132d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
133d433e6cbSMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
134d433e6cbSMatthew G. Knepley   if (dmlocalts->ijacobianlocal) {
135d433e6cbSMatthew G. Knepley     ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr);
136d433e6cbSMatthew G. Knepley     ierr = DMGetLocalVector(dm, &locX_t);CHKERRQ(ierr);
137d433e6cbSMatthew G. Knepley     ierr = VecZeroEntries(locX);CHKERRQ(ierr);
138d433e6cbSMatthew G. Knepley     ierr = VecZeroEntries(locX_t);CHKERRQ(ierr);
139d433e6cbSMatthew G. Knepley     ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
140d433e6cbSMatthew G. Knepley     ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
141d433e6cbSMatthew G. Knepley     ierr = DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t);CHKERRQ(ierr);
142d433e6cbSMatthew G. Knepley     ierr = DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t);CHKERRQ(ierr);
143d433e6cbSMatthew G. Knepley     CHKMEMQ;
144d433e6cbSMatthew G. Knepley     ierr = (*dmlocalts->ijacobianlocal)(dm, time, locX, locX_t, a, A, B, dmlocalts->ijacobianlocalctx);CHKERRQ(ierr);
145d433e6cbSMatthew G. Knepley     CHKMEMQ;
146d433e6cbSMatthew G. Knepley     ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr);
147d433e6cbSMatthew G. Knepley     ierr = DMRestoreLocalVector(dm, &locX_t);CHKERRQ(ierr);
148d433e6cbSMatthew G. Knepley   } else {
149d433e6cbSMatthew G. Knepley     MatFDColoring fdcoloring;
150d433e6cbSMatthew G. Knepley     ierr = PetscObjectQuery((PetscObject) dm, "DMDASNES_FDCOLORING", (PetscObject *) &fdcoloring);CHKERRQ(ierr);
151d433e6cbSMatthew G. Knepley     if (!fdcoloring) {
152d433e6cbSMatthew G. Knepley       ISColoring coloring;
153d433e6cbSMatthew G. Knepley 
154d433e6cbSMatthew G. Knepley       ierr = DMCreateColoring(dm, dm->coloringtype, &coloring);CHKERRQ(ierr);
155d433e6cbSMatthew G. Knepley       ierr = MatFDColoringCreate(B, coloring, &fdcoloring);CHKERRQ(ierr);
156d433e6cbSMatthew G. Knepley       ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
157d433e6cbSMatthew G. Knepley       switch (dm->coloringtype) {
158d433e6cbSMatthew G. Knepley       case IS_COLORING_GLOBAL:
159d433e6cbSMatthew G. Knepley         ierr = MatFDColoringSetFunction(fdcoloring, (PetscErrorCode (*)(void)) TSComputeIFunction_DMLocal, dmlocalts);CHKERRQ(ierr);
160d433e6cbSMatthew G. Knepley         break;
161d433e6cbSMatthew G. Knepley       default: SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]);
162d433e6cbSMatthew G. Knepley       }
163d433e6cbSMatthew G. Knepley       ierr = PetscObjectSetOptionsPrefix((PetscObject) fdcoloring, ((PetscObject) dm)->prefix);CHKERRQ(ierr);
164d433e6cbSMatthew G. Knepley       ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr);
165d433e6cbSMatthew G. Knepley       ierr = MatFDColoringSetUp(B, coloring, fdcoloring);CHKERRQ(ierr);
166d433e6cbSMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) dm, "DMDASNES_FDCOLORING", (PetscObject) fdcoloring);CHKERRQ(ierr);
167d433e6cbSMatthew G. Knepley       ierr = PetscObjectDereference((PetscObject) fdcoloring);CHKERRQ(ierr);
168d433e6cbSMatthew G. Knepley 
169d433e6cbSMatthew G. Knepley       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
170d433e6cbSMatthew G. Knepley        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
171d433e6cbSMatthew G. Knepley        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
172d433e6cbSMatthew G. Knepley        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
173d433e6cbSMatthew G. Knepley        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
174d433e6cbSMatthew G. Knepley        */
175d433e6cbSMatthew G. Knepley       ierr = PetscObjectDereference((PetscObject) dm);CHKERRQ(ierr);
176d433e6cbSMatthew G. Knepley     }
177d433e6cbSMatthew G. Knepley     ierr = MatFDColoringApply(B, fdcoloring, X, ts);CHKERRQ(ierr);
178d433e6cbSMatthew G. Knepley   }
179d433e6cbSMatthew G. Knepley   /* This will be redundant if the user called both, but it's too common to forget. */
180d433e6cbSMatthew G. Knepley   if (A != B) {
181d433e6cbSMatthew G. Knepley     ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
182d433e6cbSMatthew G. Knepley     ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
183d433e6cbSMatthew G. Knepley   }
184d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
185d433e6cbSMatthew G. Knepley }
186d433e6cbSMatthew G. Knepley 
187d433e6cbSMatthew G. Knepley #undef __FUNCT__
188*a40652d4SToby Isaac #define __FUNCT__ "DMTSSetIBoundaryLocal"
189*a40652d4SToby Isaac /*@C
190*a40652d4SToby Isaac   DMTSSetIBoundaryLocal - set the function for essential boundary data for a local implicit function evaluation.
191*a40652d4SToby 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).
192*a40652d4SToby Isaac     Vectors are initialized to zero before this function, so it is only needed for non homogeneous data.
193*a40652d4SToby Isaac 
194*a40652d4SToby Isaac   Note that this function is somewhat optional: boundary data could potentially be inserted by a function passed to
195*a40652d4SToby Isaac   DMTSSetIFunctionLocal().  The use case for this function is for discretizations with constraints (see
196*a40652d4SToby Isaac   DMGetDefaultConstraints()): this function inserts boundary values before constraint interpolation.
197*a40652d4SToby Isaac 
198*a40652d4SToby Isaac   Logically Collective
199*a40652d4SToby Isaac 
200*a40652d4SToby Isaac   Input Arguments:
201*a40652d4SToby Isaac + dm   - DM to associate callback with
202*a40652d4SToby Isaac . func - local function evaluation
203*a40652d4SToby Isaac - ctx  - context for function evaluation
204*a40652d4SToby Isaac 
205*a40652d4SToby Isaac   Level: intermediate
206*a40652d4SToby Isaac 
207*a40652d4SToby Isaac .seealso: DMTSSetIFunction(), DMTSSetIJacobianLocal()
208*a40652d4SToby Isaac @*/
209*a40652d4SToby Isaac PetscErrorCode DMTSSetIBoundaryLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx)
210*a40652d4SToby Isaac {
211*a40652d4SToby Isaac   DMTS           tdm;
212*a40652d4SToby Isaac   DMTS_Local    *dmlocalts;
213*a40652d4SToby Isaac   PetscErrorCode ierr;
214*a40652d4SToby Isaac 
215*a40652d4SToby Isaac   PetscFunctionBegin;
216*a40652d4SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
217*a40652d4SToby Isaac   ierr = DMGetDMTSWrite(dm, &tdm);CHKERRQ(ierr);
218*a40652d4SToby Isaac   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr);
219*a40652d4SToby Isaac 
220*a40652d4SToby Isaac   dmlocalts->iboundarylocal    = func;
221*a40652d4SToby Isaac   dmlocalts->iboundarylocalctx = ctx;
222*a40652d4SToby Isaac 
223*a40652d4SToby Isaac   PetscFunctionReturn(0);
224*a40652d4SToby Isaac }
225*a40652d4SToby Isaac 
226*a40652d4SToby Isaac #undef __FUNCT__
227d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMTSSetIFunctionLocal"
228d433e6cbSMatthew G. Knepley /*@C
229d433e6cbSMatthew G. Knepley   DMTSSetIFunctionLocal - set a local implicit function evaluation function. This function is called with local vector
230d433e6cbSMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
231d433e6cbSMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
232d433e6cbSMatthew G. Knepley 
233d433e6cbSMatthew G. Knepley   Logically Collective
234d433e6cbSMatthew G. Knepley 
235d433e6cbSMatthew G. Knepley   Input Arguments:
236d433e6cbSMatthew G. Knepley + dm   - DM to associate callback with
237d433e6cbSMatthew G. Knepley . func - local function evaluation
238d433e6cbSMatthew G. Knepley - ctx  - context for function evaluation
239d433e6cbSMatthew G. Knepley 
240d433e6cbSMatthew G. Knepley   Level: beginner
241d433e6cbSMatthew G. Knepley 
242d433e6cbSMatthew G. Knepley .seealso: DMTSSetIFunction(), DMTSSetIJacobianLocal()
243d433e6cbSMatthew G. Knepley @*/
244d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetIFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, Vec, void *), void *ctx)
245d433e6cbSMatthew G. Knepley {
246d433e6cbSMatthew G. Knepley   DMTS           tdm;
247d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts;
248d433e6cbSMatthew G. Knepley   PetscErrorCode ierr;
249d433e6cbSMatthew G. Knepley 
250d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
251d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
252d433e6cbSMatthew G. Knepley   ierr = DMGetDMTSWrite(dm, &tdm);CHKERRQ(ierr);
253d433e6cbSMatthew G. Knepley   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr);
254d433e6cbSMatthew G. Knepley 
255d433e6cbSMatthew G. Knepley   dmlocalts->ifunctionlocal    = func;
256d433e6cbSMatthew G. Knepley   dmlocalts->ifunctionlocalctx = ctx;
257d433e6cbSMatthew G. Knepley 
258d433e6cbSMatthew G. Knepley   ierr = DMTSSetIFunction(dm, TSComputeIFunction_DMLocal, dmlocalts);CHKERRQ(ierr);
259d433e6cbSMatthew G. Knepley   if (!tdm->ops->ijacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
260d433e6cbSMatthew G. Knepley     ierr = DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts);CHKERRQ(ierr);
261d433e6cbSMatthew G. Knepley   }
262d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
263d433e6cbSMatthew G. Knepley }
264d433e6cbSMatthew G. Knepley 
265d433e6cbSMatthew G. Knepley #undef __FUNCT__
266d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMTSSetIJacobianLocal"
267d433e6cbSMatthew G. Knepley /*@C
268d433e6cbSMatthew G. Knepley   DMTSSetIJacobianLocal - set a local Jacobian evaluation function
269d433e6cbSMatthew G. Knepley 
270d433e6cbSMatthew G. Knepley   Logically Collective
271d433e6cbSMatthew G. Knepley 
272d433e6cbSMatthew G. Knepley   Input Arguments:
273d433e6cbSMatthew G. Knepley + dm - DM to associate callback with
274d433e6cbSMatthew G. Knepley . func - local Jacobian evaluation
275d433e6cbSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation
276d433e6cbSMatthew G. Knepley 
277d433e6cbSMatthew G. Knepley   Level: beginner
278d433e6cbSMatthew G. Knepley 
279d433e6cbSMatthew G. Knepley .seealso: DMTSSetIFunctionLocal(), DMTSSetIJacobian(), DMTSSetIFunction()
280d433e6cbSMatthew G. Knepley @*/
281d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetIJacobianLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *ctx)
282d433e6cbSMatthew G. Knepley {
283d433e6cbSMatthew G. Knepley   DMTS           tdm;
284d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts;
285d433e6cbSMatthew G. Knepley   PetscErrorCode ierr;
286d433e6cbSMatthew G. Knepley 
287d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
288d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
289d433e6cbSMatthew G. Knepley   ierr = DMGetDMTSWrite(dm, &tdm);CHKERRQ(ierr);
290d433e6cbSMatthew G. Knepley   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr);
291d433e6cbSMatthew G. Knepley 
292d433e6cbSMatthew G. Knepley   dmlocalts->ijacobianlocal    = func;
293d433e6cbSMatthew G. Knepley   dmlocalts->ijacobianlocalctx = ctx;
294d433e6cbSMatthew G. Knepley 
295d433e6cbSMatthew G. Knepley   ierr = DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts);CHKERRQ(ierr);
296d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
297d433e6cbSMatthew G. Knepley }
298d433e6cbSMatthew G. Knepley 
299d433e6cbSMatthew G. Knepley #undef __FUNCT__
300*a40652d4SToby Isaac #define __FUNCT__ "DMTSSetRHSBoundaryLocal"
301*a40652d4SToby Isaac /*@C
302*a40652d4SToby Isaac   DMTSSetRHSBoundaryLocal - set the function for essential boundary data for a local righthand side function evaluation.
303*a40652d4SToby Isaac     It should set the essential boundary data for the local portion of the solution X.
304*a40652d4SToby Isaac     Vectors are initialized to zero before this function, so it is only needed for non homogeneous data.
305*a40652d4SToby Isaac 
306*a40652d4SToby Isaac   Note that this function is somewhat optional: boundary data could potentially be inserted by a function passed to
307*a40652d4SToby Isaac   DMTSSetRHSFunctionLocal().  The use case for this function is for discretizations with constraints (see
308*a40652d4SToby Isaac   DMGetDefaultConstraints()): this function inserts boundary values before constraint interpolation.
309*a40652d4SToby Isaac 
310*a40652d4SToby Isaac   Logically Collective
311*a40652d4SToby Isaac 
312*a40652d4SToby Isaac   Input Arguments:
313*a40652d4SToby Isaac + dm   - DM to associate callback with
314*a40652d4SToby Isaac . func - local function evaluation
315*a40652d4SToby Isaac - ctx  - context for function evaluation
316*a40652d4SToby Isaac 
317*a40652d4SToby Isaac   Level: intermediate
318*a40652d4SToby Isaac 
319*a40652d4SToby Isaac .seealso: DMTSSetRHSFunction()
320*a40652d4SToby Isaac @*/
321*a40652d4SToby Isaac PetscErrorCode DMTSSetRHSBoundaryLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, void *), void *ctx)
322*a40652d4SToby Isaac {
323*a40652d4SToby Isaac   DMTS           tdm;
324*a40652d4SToby Isaac   DMTS_Local    *dmlocalts;
325*a40652d4SToby Isaac   PetscErrorCode ierr;
326*a40652d4SToby Isaac 
327*a40652d4SToby Isaac   PetscFunctionBegin;
328*a40652d4SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
329*a40652d4SToby Isaac   ierr = DMGetDMTSWrite(dm, &tdm);CHKERRQ(ierr);
330*a40652d4SToby Isaac   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr);
331*a40652d4SToby Isaac 
332*a40652d4SToby Isaac   dmlocalts->rhsboundarylocal    = func;
333*a40652d4SToby Isaac   dmlocalts->rhsboundarylocalctx = ctx;
334*a40652d4SToby Isaac 
335*a40652d4SToby Isaac   PetscFunctionReturn(0);
336*a40652d4SToby Isaac }
337*a40652d4SToby Isaac 
338*a40652d4SToby Isaac #undef __FUNCT__
339d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMTSSetRHSFunctionLocal"
340d433e6cbSMatthew G. Knepley /*@C
341d433e6cbSMatthew G. Knepley   DMTSSetRHSFunctionLocal - set a local rhs function evaluation function. This function is called with local vector
342d433e6cbSMatthew G. Knepley       containing the local vector information PLUS ghost point information. It should compute a result for all local
343d433e6cbSMatthew G. Knepley       elements and DMTS will automatically accumulate the overlapping values.
344d433e6cbSMatthew G. Knepley 
345d433e6cbSMatthew G. Knepley   Logically Collective
346d433e6cbSMatthew G. Knepley 
347d433e6cbSMatthew G. Knepley   Input Arguments:
348d433e6cbSMatthew G. Knepley + dm   - DM to associate callback with
349d433e6cbSMatthew G. Knepley . func - local function evaluation
350d433e6cbSMatthew G. Knepley - ctx  - context for function evaluation
351d433e6cbSMatthew G. Knepley 
352d433e6cbSMatthew G. Knepley   Level: beginner
353d433e6cbSMatthew G. Knepley 
354d433e6cbSMatthew G. Knepley .seealso: DMTSSetRHSFunction(), DMTSSetIFunction(), DMTSSetIJacobianLocal()
355d433e6cbSMatthew G. Knepley @*/
356d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx)
357d433e6cbSMatthew G. Knepley {
358d433e6cbSMatthew G. Knepley   DMTS           tdm;
359d433e6cbSMatthew G. Knepley   DMTS_Local    *dmlocalts;
360d433e6cbSMatthew G. Knepley   PetscErrorCode ierr;
361d433e6cbSMatthew G. Knepley 
362d433e6cbSMatthew G. Knepley   PetscFunctionBegin;
363d433e6cbSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
364d433e6cbSMatthew G. Knepley   ierr = DMGetDMTSWrite(dm, &tdm);CHKERRQ(ierr);
365d433e6cbSMatthew G. Knepley   ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr);
366d433e6cbSMatthew G. Knepley 
367d433e6cbSMatthew G. Knepley   dmlocalts->rhsfunctionlocal    = func;
368d433e6cbSMatthew G. Knepley   dmlocalts->rhsfunctionlocalctx = ctx;
369d433e6cbSMatthew G. Knepley 
370d433e6cbSMatthew G. Knepley   ierr = DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMLocal, dmlocalts);CHKERRQ(ierr);
371d433e6cbSMatthew G. Knepley   PetscFunctionReturn(0);
372d433e6cbSMatthew G. Knepley }
373