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