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