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