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