xref: /petsc/src/snes/utils/dmdasnes.c (revision fe998a80077c9ee0917a39496df43fc256e1b478)
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   void       *residuallocalctx;
11   void       *jacobianlocalctx;
12   void       *objectivelocalctx;
13   InsertMode residuallocalimode;
14 
15   /*   For Picard iteration defined locally */
16   PetscErrorCode (*rhsplocal)(DMDALocalInfo*,void*,void*,void*);
17   PetscErrorCode (*jacobianplocal)(DMDALocalInfo*,void*,Mat,Mat,void*);
18   void *picardlocalctx;
19 } DMSNES_DA;
20 
21 #undef __FUNCT__
22 #define __FUNCT__ "DMSNESDestroy_DMDA"
23 static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm)
24 {
25   PetscErrorCode ierr;
26 
27   PetscFunctionBegin;
28   ierr = PetscFree(sdm->data);CHKERRQ(ierr);
29   PetscFunctionReturn(0);
30 }
31 
32 #undef __FUNCT__
33 #define __FUNCT__ "DMSNESDuplicate_DMDA"
34 static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm,DMSNES sdm)
35 {
36   PetscErrorCode ierr;
37 
38   PetscFunctionBegin;
39   ierr = PetscNewLog(sdm,(DMSNES_DA**)&sdm->data);CHKERRQ(ierr);
40   if (oldsdm->data) {
41     ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_DA));CHKERRQ(ierr);
42   }
43   PetscFunctionReturn(0);
44 }
45 
46 
47 #undef __FUNCT__
48 #define __FUNCT__ "DMDASNESGetContext"
49 static PetscErrorCode DMDASNESGetContext(DM dm,DMSNES sdm,DMSNES_DA  **dmdasnes)
50 {
51   PetscErrorCode ierr;
52 
53   PetscFunctionBegin;
54   *dmdasnes = NULL;
55   if (!sdm->data) {
56     ierr = PetscNewLog(dm,(DMSNES_DA**)&sdm->data);CHKERRQ(ierr);
57     sdm->ops->destroy   = DMSNESDestroy_DMDA;
58     sdm->ops->duplicate = DMSNESDuplicate_DMDA;
59   }
60   *dmdasnes = (DMSNES_DA*)sdm->data;
61   PetscFunctionReturn(0);
62 }
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "SNESComputeFunction_DMDA"
66 static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx)
67 {
68   PetscErrorCode ierr;
69   DM             dm;
70   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
71   DMDALocalInfo  info;
72   Vec            Xloc;
73   void           *x,*f;
74 
75   PetscFunctionBegin;
76   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
77   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
78   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
79   if (!dmdasnes->residuallocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
80   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
81   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
82   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
83   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
84   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
85   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
86   switch (dmdasnes->residuallocalimode) {
87   case INSERT_VALUES: {
88     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
89     CHKMEMQ;
90     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
91     CHKMEMQ;
92     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
93   } break;
94   case ADD_VALUES: {
95     Vec Floc;
96     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
97     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
98     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
99     CHKMEMQ;
100     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
101     CHKMEMQ;
102     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
103     ierr = VecZeroEntries(F);CHKERRQ(ierr);
104     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
105     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
106     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
107   } break;
108   default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
109   }
110   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
111   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
112   PetscFunctionReturn(0);
113 }
114 
115 #undef __FUNCT__
116 #define __FUNCT__ "SNESComputeObjective_DMDA"
117 static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx)
118 {
119   PetscErrorCode ierr;
120   DM             dm;
121   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
122   DMDALocalInfo  info;
123   Vec            Xloc;
124   void           *x;
125 
126   PetscFunctionBegin;
127   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
128   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
129   PetscValidPointer(ob,3);
130   if (!dmdasnes->objectivelocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
131   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
132   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
133   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
134   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
135   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
136   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
137   CHKMEMQ;
138   ierr = (*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx);CHKERRQ(ierr);
139   CHKMEMQ;
140   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
141   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
142   PetscFunctionReturn(0);
143 }
144 
145 
146 #undef __FUNCT__
147 #define __FUNCT__ "SNESComputeJacobian_DMDA"
148 static PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
149 {
150   PetscErrorCode ierr;
151   DM             dm;
152   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
153   DMDALocalInfo  info;
154   Vec            Xloc;
155   void           *x;
156 
157   PetscFunctionBegin;
158   if (!dmdasnes->residuallocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
159   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
160 
161   if (dmdasnes->jacobianlocal) {
162     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
163     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
164     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
165     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
166     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
167     CHKMEMQ;
168     ierr = (*dmdasnes->jacobianlocal)(&info,x,A,B,dmdasnes->jacobianlocalctx);CHKERRQ(ierr);
169     CHKMEMQ;
170     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
171     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
172   } else {
173     MatFDColoring fdcoloring;
174     ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr);
175     if (!fdcoloring) {
176       ISColoring coloring;
177 
178       ierr = DMCreateColoring(dm,dm->coloringtype,&coloring);CHKERRQ(ierr);
179       ierr = MatFDColoringCreate(B,coloring,&fdcoloring);CHKERRQ(ierr);
180       switch (dm->coloringtype) {
181       case IS_COLORING_GLOBAL:
182         ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
183         break;
184       default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
185       }
186       ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr);
187       ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr);
188       ierr = MatFDColoringSetUp(B,coloring,fdcoloring);CHKERRQ(ierr);
189       ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
190       ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr);
191       ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr);
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       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
200     }
201     ierr  = MatFDColoringApply(B,fdcoloring,X,snes);CHKERRQ(ierr);
202   }
203   /* This will be redundant if the user called both, but it's too common to forget. */
204   if (A != B) {
205     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
206     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
207   }
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "DMDASNESSetFunctionLocal"
213 /*@C
214    DMDASNESSetFunctionLocal - set a local residual evaluation function
215 
216    Logically Collective
217 
218    Input Arguments:
219 +  dm - DM to associate callback with
220 .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
221 .  func - local residual evaluation
222 -  ctx - optional context for local residual evaluation
223 
224    Calling sequence for func:
225 +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
226 .  x - dimensional pointer to state at which to evaluate residual
227 .  f - dimensional pointer to residual, write the residual here
228 -  ctx - optional context passed above
229 
230    Level: beginner
231 
232 .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
233 @*/
234 PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx)
235 {
236   PetscErrorCode ierr;
237   DMSNES         sdm;
238   DMSNES_DA      *dmdasnes;
239 
240   PetscFunctionBegin;
241   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
242   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
243   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
244 
245   dmdasnes->residuallocalimode = imode;
246   dmdasnes->residuallocal      = func;
247   dmdasnes->residuallocalctx   = ctx;
248 
249   ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
250   if (!sdm->ops->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
251     ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
252   }
253   PetscFunctionReturn(0);
254 }
255 
256 #undef __FUNCT__
257 #define __FUNCT__ "DMDASNESSetJacobianLocal"
258 /*@C
259    DMDASNESSetJacobianLocal - set a local Jacobian evaluation function
260 
261    Logically Collective
262 
263    Input Arguments:
264 +  dm - DM to associate callback with
265 .  func - local residual evaluation
266 -  ctx - optional context for local residual evaluation
267 
268    Calling sequence for func:
269 +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
270 .  x - dimensional pointer to state at which to evaluate residual
271 .  f - dimensional pointer to residual, write the residual here
272 -  ctx - optional context passed above
273 
274    Level: beginner
275 
276 .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
277 @*/
278 PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
279 {
280   PetscErrorCode ierr;
281   DMSNES         sdm;
282   DMSNES_DA      *dmdasnes;
283 
284   PetscFunctionBegin;
285   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
286   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
287   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
288 
289   dmdasnes->jacobianlocal    = func;
290   dmdasnes->jacobianlocalctx = ctx;
291 
292   ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
293   PetscFunctionReturn(0);
294 }
295 
296 
297 #undef __FUNCT__
298 #define __FUNCT__ "DMDASNESSetObjectiveLocal"
299 /*@C
300    DMDASNESSetObjectiveLocal - set a local residual evaluation function
301 
302    Logically Collective
303 
304    Input Arguments:
305 +  dm - DM to associate callback with
306 .  func - local objective evaluation
307 -  ctx - optional context for local residual evaluation
308 
309    Calling sequence for func:
310 +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
311 .  x - dimensional pointer to state at which to evaluate residual
312 .  ob - eventual objective value
313 -  ctx - optional context passed above
314 
315    Level: beginner
316 
317 .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
318 @*/
319 PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx)
320 {
321   PetscErrorCode ierr;
322   DMSNES         sdm;
323   DMSNES_DA      *dmdasnes;
324 
325   PetscFunctionBegin;
326   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
327   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
328   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
329 
330   dmdasnes->objectivelocal    = func;
331   dmdasnes->objectivelocalctx = ctx;
332 
333   ierr = DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes);CHKERRQ(ierr);
334   PetscFunctionReturn(0);
335 }
336 
337 #undef __FUNCT__
338 #define __FUNCT__ "SNESComputePicard_DMDA"
339 static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx)
340 {
341   PetscErrorCode ierr;
342   DM             dm;
343   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
344   DMDALocalInfo  info;
345   Vec            Xloc;
346   void           *x,*f;
347 
348   PetscFunctionBegin;
349   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
350   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
351   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
352   if (!dmdasnes->rhsplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
353   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
354   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
355   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
356   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
357   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
358   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
359   switch (dmdasnes->residuallocalimode) {
360   case INSERT_VALUES: {
361     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
362     CHKMEMQ;
363     ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr);
364     CHKMEMQ;
365     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
366   } break;
367   case ADD_VALUES: {
368     Vec Floc;
369     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
370     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
371     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
372     CHKMEMQ;
373     ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr);
374     CHKMEMQ;
375     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
376     ierr = VecZeroEntries(F);CHKERRQ(ierr);
377     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
378     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
379     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
380   } break;
381   default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
382   }
383   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
384   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
385   PetscFunctionReturn(0);
386 }
387 
388 #undef __FUNCT__
389 #define __FUNCT__ "SNESComputePicardJacobian_DMDA"
390 static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
391 {
392   PetscErrorCode ierr;
393   DM             dm;
394   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
395   DMDALocalInfo  info;
396   Vec            Xloc;
397   void           *x;
398 
399   PetscFunctionBegin;
400   if (!dmdasnes->jacobianplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
401   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
402 
403   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
404   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
405   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
406   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
407   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
408   CHKMEMQ;
409   ierr = (*dmdasnes->jacobianplocal)(&info,x,A,B,dmdasnes->picardlocalctx);CHKERRQ(ierr);
410   CHKMEMQ;
411   ierr  = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
412   ierr  = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
413   if (A != B) {
414     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
415     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
416   }
417   PetscFunctionReturn(0);
418 }
419 
420 #undef __FUNCT__
421 #define __FUNCT__ "DMDASNESSetPicardLocal"
422 /*@C
423    DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration
424 
425    Logically Collective
426 
427    Input Arguments:
428 +  dm - DM to associate callback with
429 .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
430 .  func - local residual evaluation
431 -  ctx - optional context for local residual evaluation
432 
433    Calling sequence for func:
434 +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
435 .  x - dimensional pointer to state at which to evaluate residual
436 .  f - dimensional pointer to residual, write the residual here
437 -  ctx - optional context passed above
438 
439    Notes:  The user must use
440     extern PetscErrorCode  SNESPicardComputeFunction(SNES,Vec,Vec,void*);
441     extern PetscErrorCode  SNESPicardComputeJacobian(SNES,Vec,Mat,Mat,MatStructure*,void*);
442     ierr = SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user);CHKERRQ(ierr);
443     in their code before calling this routine.
444 
445 
446    Level: beginner
447 
448 .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
449 @*/
450 PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),
451                                       PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
452 {
453   PetscErrorCode ierr;
454   DMSNES         sdm;
455   DMSNES_DA      *dmdasnes;
456 
457   PetscFunctionBegin;
458   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
459   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
460   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
461 
462   dmdasnes->residuallocalimode = imode;
463   dmdasnes->rhsplocal          = func;
464   dmdasnes->jacobianplocal     = jac;
465   dmdasnes->picardlocalctx     = ctx;
466 
467   ierr = DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
468   PetscFunctionReturn(0);
469 }
470