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