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