xref: /petsc/src/snes/utils/dmdasnes.c (revision 6a4a1270853b5b1995c94de98f44d28bec57470a)
1 #include <petscdmda.h> /*I "petscdmda.h" I*/
2 #include <petsc/private/dmimpl.h>
3 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
4 
5 /* This structure holds the user-provided DMDA callbacks */
6 typedef struct {
7   PetscErrorCode (*residuallocal)(DMDALocalInfo *, void *, void *, void *);
8   PetscErrorCode (*jacobianlocal)(DMDALocalInfo *, void *, Mat, Mat, void *);
9   PetscErrorCode (*objectivelocal)(DMDALocalInfo *, void *, PetscReal *, void *);
10 
11   PetscErrorCode (*residuallocalvec)(DMDALocalInfo *, Vec, Vec, void *);
12   PetscErrorCode (*jacobianlocalvec)(DMDALocalInfo *, Vec, Mat, Mat, void *);
13   PetscErrorCode (*objectivelocalvec)(DMDALocalInfo *, Vec, PetscReal *, void *);
14   void      *residuallocalctx;
15   void      *jacobianlocalctx;
16   void      *objectivelocalctx;
17   InsertMode residuallocalimode;
18 
19   /*   For Picard iteration defined locally */
20   PetscErrorCode (*rhsplocal)(DMDALocalInfo *, void *, void *, void *);
21   PetscErrorCode (*jacobianplocal)(DMDALocalInfo *, void *, Mat, Mat, void *);
22   void *picardlocalctx;
23 } DMSNES_DA;
24 
25 static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm) {
26   PetscFunctionBegin;
27   PetscCall(PetscFree(sdm->data));
28   PetscFunctionReturn(0);
29 }
30 
31 static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm, DMSNES sdm) {
32   PetscFunctionBegin;
33   PetscCall(PetscNewLog(sdm, (DMSNES_DA **)&sdm->data));
34   if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data, oldsdm->data, sizeof(DMSNES_DA)));
35   PetscFunctionReturn(0);
36 }
37 
38 static PetscErrorCode DMDASNESGetContext(DM dm, DMSNES sdm, DMSNES_DA **dmdasnes) {
39   PetscFunctionBegin;
40   *dmdasnes = NULL;
41   if (!sdm->data) {
42     PetscCall(PetscNewLog(dm, (DMSNES_DA **)&sdm->data));
43     sdm->ops->destroy   = DMSNESDestroy_DMDA;
44     sdm->ops->duplicate = DMSNESDuplicate_DMDA;
45   }
46   *dmdasnes = (DMSNES_DA *)sdm->data;
47   PetscFunctionReturn(0);
48 }
49 
50 static PetscErrorCode SNESComputeFunction_DMDA(SNES snes, Vec X, Vec F, void *ctx) {
51   DM            dm;
52   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
53   DMDALocalInfo info;
54   Vec           Xloc;
55   void         *x, *f;
56 
57   PetscFunctionBegin;
58   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
59   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
60   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
61   PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
62   PetscCall(SNESGetDM(snes, &dm));
63   PetscCall(DMGetLocalVector(dm, &Xloc));
64   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
65   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
66   PetscCall(DMDAGetLocalInfo(dm, &info));
67   switch (dmdasnes->residuallocalimode) {
68   case INSERT_VALUES: {
69     PetscCall(PetscLogEventBegin(SNES_FunctionEval, snes, X, F, 0));
70     if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocalvec)(&info, Xloc, F, dmdasnes->residuallocalctx));
71     else {
72       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
73       PetscCall(DMDAVecGetArray(dm, F, &f));
74       PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocal)(&info, x, f, dmdasnes->residuallocalctx));
75       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
76       PetscCall(DMDAVecRestoreArray(dm, F, &f));
77     }
78     PetscCall(PetscLogEventEnd(SNES_FunctionEval, snes, X, F, 0));
79   } break;
80   case ADD_VALUES: {
81     Vec Floc;
82     PetscCall(DMGetLocalVector(dm, &Floc));
83     PetscCall(VecZeroEntries(Floc));
84     PetscCall(PetscLogEventBegin(SNES_FunctionEval, snes, X, F, 0));
85     if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocalvec)(&info, Xloc, Floc, dmdasnes->residuallocalctx));
86     else {
87       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
88       PetscCall(DMDAVecGetArray(dm, Floc, &f));
89       PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocal)(&info, x, f, dmdasnes->residuallocalctx));
90       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
91       PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
92     }
93     PetscCall(PetscLogEventEnd(SNES_FunctionEval, snes, X, F, 0));
94     PetscCall(VecZeroEntries(F));
95     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
96     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
97     PetscCall(DMRestoreLocalVector(dm, &Floc));
98   } break;
99   default: SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdasnes->residuallocalimode);
100   }
101   PetscCall(DMRestoreLocalVector(dm, &Xloc));
102   if (snes->domainerror) PetscCall(VecSetInf(F));
103   PetscFunctionReturn(0);
104 }
105 
106 static PetscErrorCode SNESComputeObjective_DMDA(SNES snes, Vec X, PetscReal *ob, void *ctx) {
107   DM            dm;
108   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
109   DMDALocalInfo info;
110   Vec           Xloc;
111   void         *x;
112 
113   PetscFunctionBegin;
114   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
115   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
116   PetscValidRealPointer(ob, 3);
117   PetscCheck(dmdasnes->objectivelocal || dmdasnes->objectivelocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
118   PetscCall(SNESGetDM(snes, &dm));
119   PetscCall(DMGetLocalVector(dm, &Xloc));
120   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
121   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
122   PetscCall(DMDAGetLocalInfo(dm, &info));
123   if (dmdasnes->objectivelocalvec) PetscCallBack("SNES DMDA local callback objective", (*dmdasnes->objectivelocalvec)(&info, Xloc, ob, dmdasnes->objectivelocalctx));
124   else {
125     PetscCall(DMDAVecGetArray(dm, Xloc, &x));
126     PetscCallBack("SNES DMDA local callback objective", (*dmdasnes->objectivelocal)(&info, x, ob, dmdasnes->objectivelocalctx));
127     PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
128   }
129   PetscCall(DMRestoreLocalVector(dm, &Xloc));
130   PetscFunctionReturn(0);
131 }
132 
133 /* Routine is called by example, hence must be labeled PETSC_EXTERN */
134 PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes, Vec X, Mat A, Mat B, void *ctx) {
135   DM            dm;
136   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
137   DMDALocalInfo info;
138   Vec           Xloc;
139   void         *x;
140 
141   PetscFunctionBegin;
142   PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
143   PetscCall(SNESGetDM(snes, &dm));
144 
145   if (dmdasnes->jacobianlocal || dmdasnes->jacobianlocalctx) {
146     PetscCall(DMGetLocalVector(dm, &Xloc));
147     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
148     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
149     PetscCall(DMDAGetLocalInfo(dm, &info));
150     if (dmdasnes->jacobianlocalvec) PetscCallBack("SNES DMDA local callback Jacobian", (*dmdasnes->jacobianlocalvec)(&info, Xloc, A, B, dmdasnes->jacobianlocalctx));
151     else {
152       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
153       PetscCallBack("SNES DMDA local callback Jacobian", (*dmdasnes->jacobianlocal)(&info, x, A, B, dmdasnes->jacobianlocalctx));
154       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
155     }
156     PetscCall(DMRestoreLocalVector(dm, &Xloc));
157   } else {
158     MatFDColoring fdcoloring;
159     PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring));
160     if (!fdcoloring) {
161       ISColoring coloring;
162 
163       PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring));
164       PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring));
165       switch (dm->coloringtype) {
166       case IS_COLORING_GLOBAL: PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))SNESComputeFunction_DMDA, dmdasnes)); break;
167       default: SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]);
168       }
169       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix));
170       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
171       PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring));
172       PetscCall(ISColoringDestroy(&coloring));
173       PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring));
174       PetscCall(PetscObjectDereference((PetscObject)fdcoloring));
175 
176       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
177        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
178        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
179        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
180        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
181        */
182       PetscCall(PetscObjectDereference((PetscObject)dm));
183     }
184     PetscCall(MatFDColoringApply(B, fdcoloring, X, snes));
185   }
186   /* This will be redundant if the user called both, but it's too common to forget. */
187   if (A != B) {
188     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
189     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
190   }
191   PetscFunctionReturn(0);
192 }
193 
194 /*@C
195    DMDASNESSetFunctionLocal - set a local residual evaluation function for use with `DMDA`
196 
197    Logically Collective
198 
199    Input Parameters:
200 +  dm - `DM` to associate callback with
201 .  imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
202 .  func - local residual evaluation
203 -  ctx - optional context for local residual evaluation
204 
205    Calling sequence:
206    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x, void *f, void *ctx),
207 +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
208 .  x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x)
209 .  f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f)
210 -  ctx - optional context passed above
211 
212    Level: beginner
213 
214 .seealso: `DMDA`, `DMDASNESSetJacobianLocal()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
215 @*/
216 PetscErrorCode DMDASNESSetFunctionLocal(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *, void *, void *, void *), void *ctx) {
217   DMSNES     sdm;
218   DMSNES_DA *dmdasnes;
219 
220   PetscFunctionBegin;
221   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
222   PetscCall(DMGetDMSNESWrite(dm, &sdm));
223   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
224 
225   dmdasnes->residuallocalimode = imode;
226   dmdasnes->residuallocal      = func;
227   dmdasnes->residuallocalctx   = ctx;
228 
229   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
230   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
231     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
232   }
233   PetscFunctionReturn(0);
234 }
235 
236 /*@C
237    DMDASNESSetFunctionLocalVec - set a local residual evaluation function that operates on a local vector for `DMDA`
238 
239    Logically Collective
240 
241    Input Parameters:
242 +  dm - `DM` to associate callback with
243 .  imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
244 .  func - local residual evaluation
245 -  ctx - optional context for local residual evaluation
246 
247    Calling sequence:
248    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x, Vec f, void *ctx),
249 +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
250 .  x - state vector at which to evaluate residual
251 .  f - residual vector
252 -  ctx - optional context passed above
253 
254    Level: beginner
255 
256 .seealso: `DMDA`, `DMDASNESSetFunctionLocal()`, `DMDASNESSetJacobianLocalVec()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
257 @*/
258 PetscErrorCode DMDASNESSetFunctionLocalVec(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *, Vec, Vec, void *), void *ctx) {
259   DMSNES     sdm;
260   DMSNES_DA *dmdasnes;
261 
262   PetscFunctionBegin;
263   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
264   PetscCall(DMGetDMSNESWrite(dm, &sdm));
265   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
266 
267   dmdasnes->residuallocalimode = imode;
268   dmdasnes->residuallocalvec   = func;
269   dmdasnes->residuallocalctx   = ctx;
270 
271   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
272   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
273     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
274   }
275   PetscFunctionReturn(0);
276 }
277 
278 /*@C
279    DMDASNESSetJacobianLocal - set a local Jacobian evaluation function for use with `DMDA`
280 
281    Logically Collective
282 
283    Input Parameters:
284 +  dm - `DM` to associate callback with
285 .  func - local Jacobian evaluation
286 -  ctx - optional context for local Jacobian evaluation
287 
288    Calling sequence:
289    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x,Mat J,Mat M,void *ctx),
290 +  info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
291 .  x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x)
292 .  J - Mat object for the Jacobian
293 .  M - Mat object for the Jacobian preconditioner matrix
294 -  ctx - optional context passed above
295 
296    Level: beginner
297 
298 .seealso: `DMDA`, `DMDASNESSetFunctionLocal()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
299 @*/
300 PetscErrorCode DMDASNESSetJacobianLocal(DM dm, PetscErrorCode (*func)(DMDALocalInfo *, void *, Mat, Mat, void *), void *ctx) {
301   DMSNES     sdm;
302   DMSNES_DA *dmdasnes;
303 
304   PetscFunctionBegin;
305   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
306   PetscCall(DMGetDMSNESWrite(dm, &sdm));
307   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
308 
309   dmdasnes->jacobianlocal    = func;
310   dmdasnes->jacobianlocalctx = ctx;
311 
312   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
313   PetscFunctionReturn(0);
314 }
315 
316 /*@C
317    DMDASNESSetJacobianLocalVec - set a local Jacobian evaluation function that operates on a local vector with `DMDA`
318 
319    Logically Collective
320 
321    Input Parameters:
322 +  dm - `DM` to associate callback with
323 .  func - local Jacobian evaluation
324 -  ctx - optional context for local Jacobian evaluation
325 
326    Calling sequence:
327    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,Mat J,Mat M,void *ctx),
328 +  info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
329 .  x - state vector at which to evaluate Jacobian
330 .  J - Mat object for the Jacobian
331 .  M - Mat object for the Jacobian preconditioner matrix
332 -  ctx - optional context passed above
333 
334    Level: beginner
335 
336 .seealso: `DMDA`, `DMDASNESSetJacobianLocal()`, `DMDASNESSetFunctionLocalVec()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
337 @*/
338 PetscErrorCode DMDASNESSetJacobianLocalVec(DM dm, PetscErrorCode (*func)(DMDALocalInfo *, Vec, Mat, Mat, void *), void *ctx) {
339   DMSNES     sdm;
340   DMSNES_DA *dmdasnes;
341 
342   PetscFunctionBegin;
343   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
344   PetscCall(DMGetDMSNESWrite(dm, &sdm));
345   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
346 
347   dmdasnes->jacobianlocalvec = func;
348   dmdasnes->jacobianlocalctx = ctx;
349 
350   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
351   PetscFunctionReturn(0);
352 }
353 
354 /*@C
355    DMDASNESSetObjectiveLocal - set a local residual evaluation function to used with a `DMDA`
356 
357    Logically Collective
358 
359    Input Parameters:
360 +  dm - `DM` to associate callback with
361 .  func - local objective evaluation
362 -  ctx - optional context for local residual evaluation
363 
364    Calling sequence for func:
365 +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
366 .  x - dimensional pointer to state at which to evaluate residual
367 .  ob - eventual objective value
368 -  ctx - optional context passed above
369 
370    Level: beginner
371 
372 .seealso: `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
373 @*/
374 PetscErrorCode DMDASNESSetObjectiveLocal(DM dm, DMDASNESObjective func, void *ctx) {
375   DMSNES     sdm;
376   DMSNES_DA *dmdasnes;
377 
378   PetscFunctionBegin;
379   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
380   PetscCall(DMGetDMSNESWrite(dm, &sdm));
381   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
382 
383   dmdasnes->objectivelocal    = func;
384   dmdasnes->objectivelocalctx = ctx;
385 
386   PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes));
387   PetscFunctionReturn(0);
388 }
389 
390 /*@C
391    DMDASNESSetObjectiveLocal - set a local residual evaluation function that operates on a local vector with `DMDA`
392 
393    Logically Collective
394 
395    Input Parameters:
396 +  dm - `DM` to associate callback with
397 .  func - local objective evaluation
398 -  ctx - optional context for local residual evaluation
399 
400    Calling sequence
401    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,PetscReal *ob,void *ctx);
402 +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
403 .  x - state vector at which to evaluate residual
404 .  ob - eventual objective value
405 -  ctx - optional context passed above
406 
407    Level: beginner
408 
409 .seealso: `DMDA`, `DMDASNESSetObjectiveLocal()`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocalVec()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
410 @*/
411 PetscErrorCode DMDASNESSetObjectiveLocalVec(DM dm, DMDASNESObjectiveVec func, void *ctx) {
412   DMSNES     sdm;
413   DMSNES_DA *dmdasnes;
414 
415   PetscFunctionBegin;
416   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
417   PetscCall(DMGetDMSNESWrite(dm, &sdm));
418   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
419 
420   dmdasnes->objectivelocalvec = func;
421   dmdasnes->objectivelocalctx = ctx;
422 
423   PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes));
424   PetscFunctionReturn(0);
425 }
426 
427 static PetscErrorCode SNESComputePicard_DMDA(SNES snes, Vec X, Vec F, void *ctx) {
428   DM            dm;
429   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
430   DMDALocalInfo info;
431   Vec           Xloc;
432   void         *x, *f;
433 
434   PetscFunctionBegin;
435   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
436   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
437   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
438   PetscCheck(dmdasnes->rhsplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
439   PetscCall(SNESGetDM(snes, &dm));
440   PetscCall(DMGetLocalVector(dm, &Xloc));
441   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
442   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
443   PetscCall(DMDAGetLocalInfo(dm, &info));
444   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
445   switch (dmdasnes->residuallocalimode) {
446   case INSERT_VALUES: {
447     PetscCall(DMDAVecGetArray(dm, F, &f));
448     PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx));
449     PetscCall(DMDAVecRestoreArray(dm, F, &f));
450   } break;
451   case ADD_VALUES: {
452     Vec Floc;
453     PetscCall(DMGetLocalVector(dm, &Floc));
454     PetscCall(VecZeroEntries(Floc));
455     PetscCall(DMDAVecGetArray(dm, Floc, &f));
456     PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx));
457     PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
458     PetscCall(VecZeroEntries(F));
459     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
460     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
461     PetscCall(DMRestoreLocalVector(dm, &Floc));
462   } break;
463   default: SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdasnes->residuallocalimode);
464   }
465   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
466   PetscCall(DMRestoreLocalVector(dm, &Xloc));
467   PetscFunctionReturn(0);
468 }
469 
470 static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes, Vec X, Mat A, Mat B, void *ctx) {
471   DM            dm;
472   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
473   DMDALocalInfo info;
474   Vec           Xloc;
475   void         *x;
476 
477   PetscFunctionBegin;
478   PetscCheck(dmdasnes->jacobianplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
479   PetscCall(SNESGetDM(snes, &dm));
480 
481   PetscCall(DMGetLocalVector(dm, &Xloc));
482   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
483   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
484   PetscCall(DMDAGetLocalInfo(dm, &info));
485   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
486   PetscCallBack("SNES Picard DMDA local callback Jacobian", (*dmdasnes->jacobianplocal)(&info, x, A, B, dmdasnes->picardlocalctx));
487   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
488   PetscCall(DMRestoreLocalVector(dm, &Xloc));
489   if (A != B) {
490     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
491     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
492   }
493   PetscFunctionReturn(0);
494 }
495 
496 /*@C
497    DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration with `DMDA`
498 
499    Logically Collective
500 
501    Input Parameters:
502 +  dm - `DM` to associate callback with
503 .  imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
504 .  func - local residual evaluation
505 -  ctx - optional context for local residual evaluation
506 
507    Calling sequence for func:
508 +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
509 .  x - dimensional pointer to state at which to evaluate residual
510 .  f - dimensional pointer to residual, write the residual here
511 -  ctx - optional context passed above
512 
513    Note:
514     The user must use `SNESSetFunction`(snes,NULL,`SNESPicardComputeFunction`,&user));
515     in their code before calling this routine.
516 
517    Level: beginner
518 
519 .seealso: `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
520 @*/
521 PetscErrorCode DMDASNESSetPicardLocal(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *, void *, void *, void *), PetscErrorCode (*jac)(DMDALocalInfo *, void *, Mat, Mat, void *), void *ctx) {
522   DMSNES     sdm;
523   DMSNES_DA *dmdasnes;
524 
525   PetscFunctionBegin;
526   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
527   PetscCall(DMGetDMSNESWrite(dm, &sdm));
528   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
529 
530   dmdasnes->residuallocalimode = imode;
531   dmdasnes->rhsplocal          = func;
532   dmdasnes->jacobianplocal     = jac;
533   dmdasnes->picardlocalctx     = ctx;
534 
535   PetscCall(DMSNESSetPicard(dm, SNESComputePicard_DMDA, SNESComputePicardJacobian_DMDA, dmdasnes));
536   PetscCall(DMSNESSetMFFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
537   PetscFunctionReturn(0);
538 }
539