xref: /petsc/src/dm/impls/da/da2.c (revision 0e9bae810fdaeb60e2713eaa8ddb89f42e079fd1)
1 
2 #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMView_DA_2d"
6 PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer)
7 {
8   PetscErrorCode ierr;
9   PetscMPIInt    rank;
10   PetscBool      iascii,isdraw,isbinary;
11   DM_DA          *dd = (DM_DA*)da->data;
12 #if defined(PETSC_HAVE_MATLAB_ENGINE)
13   PetscBool      ismatlab;
14 #endif
15 
16   PetscFunctionBegin;
17   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
18 
19   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
20   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
21   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
22 #if defined(PETSC_HAVE_MATLAB_ENGINE)
23   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
24 #endif
25   if (iascii) {
26     PetscViewerFormat format;
27 
28     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
29     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
30       DMDALocalInfo info;
31       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
32       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
33       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D m %D n %D w %D s %D\n",rank,dd->M,dd->N,dd->m,dd->n,dd->w,dd->s);CHKERRQ(ierr);
34       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D\n",info.xs,info.xs+info.xm,info.ys,info.ys+info.ym);CHKERRQ(ierr);
35       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
36       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
37     } else {
38       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
39     }
40   } else if (isdraw) {
41     PetscDraw  draw;
42     double     ymin = -1*dd->s-1,ymax = dd->N+dd->s;
43     double     xmin = -1*dd->s-1,xmax = dd->M+dd->s;
44     double     x,y;
45     PetscInt   base,*idx;
46     char       node[10];
47     PetscBool  isnull;
48 
49     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
50     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
51     if (!dd->coordinates) {
52       ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
53     }
54     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
55 
56     /* first processor draw all node lines */
57     if (!rank) {
58       ymin = 0.0; ymax = dd->N - 1;
59       for (xmin=0; xmin<dd->M; xmin++) {
60         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
61       }
62       xmin = 0.0; xmax = dd->M - 1;
63       for (ymin=0; ymin<dd->N; ymin++) {
64         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
65       }
66     }
67     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
68     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
69 
70     /* draw my box */
71     ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w;
72     xmax =(dd->xe-1)/dd->w;
73     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
74     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
75     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
76     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
77 
78     /* put in numbers */
79     base = (dd->base)/dd->w;
80     for (y=ymin; y<=ymax; y++) {
81       for (x=xmin; x<=xmax; x++) {
82         sprintf(node,"%d",(int)base++);
83         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
84       }
85     }
86 
87     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
88     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
89     /* overlay ghost numbers, useful for error checking */
90     /* put in numbers */
91 
92     base = 0; idx = dd->idx;
93     ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe;
94     for (y=ymin; y<ymax; y++) {
95       for (x=xmin; x<xmax; x++) {
96         if ((base % dd->w) == 0) {
97           sprintf(node,"%d",(int)(idx[base]/dd->w));
98           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
99         }
100         base++;
101       }
102     }
103     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
104     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
105   } else if (isbinary){
106     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
107 #if defined(PETSC_HAVE_MATLAB_ENGINE)
108   } else if (ismatlab) {
109     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
110 #endif
111   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name);
112   PetscFunctionReturn(0);
113 }
114 
115 /*
116       M is number of grid points
117       m is number of processors
118 
119 */
120 #undef __FUNCT__
121 #define __FUNCT__ "DMDASplitComm2d"
122 PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
123 {
124   PetscErrorCode ierr;
125   PetscInt       m,n = 0,x = 0,y = 0;
126   PetscMPIInt    size,csize,rank;
127 
128   PetscFunctionBegin;
129   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
130   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
131 
132   csize = 4*size;
133   do {
134     if (csize % 4) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot split communicator of size %d tried %d %D %D",size,csize,x,y);
135     csize   = csize/4;
136 
137     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
138     if (!m) m = 1;
139     while (m > 0) {
140       n = csize/m;
141       if (m*n == csize) break;
142       m--;
143     }
144     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
145 
146     x = M/m + ((M % m) > ((csize-1) % m));
147     y = (N + (csize-1)/m)/n;
148   } while ((x < 4 || y < 4) && csize > 1);
149   if (size != csize) {
150     MPI_Group    entire_group,sub_group;
151     PetscMPIInt  i,*groupies;
152 
153     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
154     ierr = PetscMalloc(csize*sizeof(PetscInt),&groupies);CHKERRQ(ierr);
155     for (i=0; i<csize; i++) {
156       groupies[i] = (rank/csize)*csize + i;
157     }
158     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
159     ierr = PetscFree(groupies);CHKERRQ(ierr);
160     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
161     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
162     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
163     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
164   } else {
165     *outcomm = comm;
166   }
167   PetscFunctionReturn(0);
168 }
169 
170 #undef __FUNCT__
171 #define __FUNCT__ "DMDAFunction"
172 static PetscErrorCode DMDAFunction(DM dm,Vec x,Vec F)
173 {
174   PetscErrorCode ierr;
175   Vec            localX;
176 
177   PetscFunctionBegin;
178   ierr = DMGetLocalVector(dm,&localX);CHKERRQ(ierr);
179   ierr = DMGlobalToLocalBegin(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
180   ierr = DMGlobalToLocalEnd(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
181   ierr = DMDAComputeFunction1(dm,localX,F,dm->ctx);CHKERRQ(ierr);
182   ierr = DMRestoreLocalVector(dm,&localX);CHKERRQ(ierr);
183   PetscFunctionReturn(0);
184 }
185 
186 #undef __FUNCT__
187 #define __FUNCT__ "DMDASetLocalFunction"
188 /*@C
189        DMDASetLocalFunction - Caches in a DM a local function.
190 
191    Logically Collective on DMDA
192 
193    Input Parameter:
194 +  da - initial distributed array
195 -  lf - the local function
196 
197    Level: intermediate
198 
199    Notes:
200       If you use SNESSetDM() and did not use SNESSetFunction() then the context passed to your function is the context set with DMSetApplicationContext()
201 
202    Developer Notes: It is possibly confusing which context is passed to the user function, it would be nice to unify them somehow.
203 
204 .keywords:  distributed array, refine
205 
206 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunctioni()
207 @*/
208 PetscErrorCode  DMDASetLocalFunction(DM da,DMDALocalFunction1 lf)
209 {
210   PetscErrorCode ierr;
211   DM_DA          *dd = (DM_DA*)da->data;
212 
213   PetscFunctionBegin;
214   PetscValidHeaderSpecific(da,DM_CLASSID,1);
215   ierr = DMSetFunction(da,DMDAFunction);CHKERRQ(ierr);
216   dd->lf       = lf;
217   PetscFunctionReturn(0);
218 }
219 
220 #undef __FUNCT__
221 #define __FUNCT__ "DMDASetLocalFunctioni"
222 /*@C
223        DMDASetLocalFunctioni - Caches in a DM a local function that evaluates a single component
224 
225    Logically Collective on DMDA
226 
227    Input Parameter:
228 +  da - initial distributed array
229 -  lfi - the local function
230 
231    Level: intermediate
232 
233 .keywords:  distributed array, refine
234 
235 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
236 @*/
237 PetscErrorCode  DMDASetLocalFunctioni(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
238 {
239   DM_DA          *dd = (DM_DA*)da->data;
240   PetscFunctionBegin;
241   PetscValidHeaderSpecific(da,DM_CLASSID,1);
242   dd->lfi = lfi;
243   PetscFunctionReturn(0);
244 }
245 
246 #undef __FUNCT__
247 #define __FUNCT__ "DMDASetLocalFunctionib"
248 /*@C
249        DMDASetLocalFunctionib - Caches in a DM a block local function that evaluates a single component
250 
251    Logically Collective on DMDA
252 
253    Input Parameter:
254 +  da - initial distributed array
255 -  lfi - the local function
256 
257    Level: intermediate
258 
259 .keywords:  distributed array, refine
260 
261 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
262 @*/
263 PetscErrorCode  DMDASetLocalFunctionib(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
264 {
265   DM_DA          *dd = (DM_DA*)da->data;
266   PetscFunctionBegin;
267   PetscValidHeaderSpecific(da,DM_CLASSID,1);
268   dd->lfib = lfi;
269   PetscFunctionReturn(0);
270 }
271 
272 #undef __FUNCT__
273 #define __FUNCT__ "DMDASetLocalAdicFunction_Private"
274 PetscErrorCode DMDASetLocalAdicFunction_Private(DM da,DMDALocalFunction1 ad_lf)
275 {
276   DM_DA          *dd = (DM_DA*)da->data;
277   PetscFunctionBegin;
278   PetscValidHeaderSpecific(da,DM_CLASSID,1);
279   dd->adic_lf = ad_lf;
280   PetscFunctionReturn(0);
281 }
282 
283 /*MC
284        DMDASetLocalAdicFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
285 
286    Synopsis:
287    PetscErrorCode DMDASetLocalAdicFunctioni(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
288 
289    Logically Collective on DMDA
290 
291    Input Parameter:
292 +  da - initial distributed array
293 -  ad_lfi - the local function as computed by ADIC/ADIFOR
294 
295    Level: intermediate
296 
297 .keywords:  distributed array, refine
298 
299 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
300           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
301 M*/
302 
303 #undef __FUNCT__
304 #define __FUNCT__ "DMDASetLocalAdicFunctioni_Private"
305 PetscErrorCode DMDASetLocalAdicFunctioni_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
306 {
307   DM_DA          *dd = (DM_DA*)da->data;
308   PetscFunctionBegin;
309   PetscValidHeaderSpecific(da,DM_CLASSID,1);
310   dd->adic_lfi = ad_lfi;
311   PetscFunctionReturn(0);
312 }
313 
314 /*MC
315        DMDASetLocalAdicMFFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
316 
317    Synopsis:
318    PetscErrorCode  DMDASetLocalAdicFunctioni(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
319 
320    Logically Collective on DMDA
321 
322    Input Parameter:
323 +  da - initial distributed array
324 -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
325 
326    Level: intermediate
327 
328 .keywords:  distributed array, refine
329 
330 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
331           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
332 M*/
333 
334 #undef __FUNCT__
335 #define __FUNCT__ "DMDASetLocalAdicMFFunctioni_Private"
336 PetscErrorCode DMDASetLocalAdicMFFunctioni_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
337 {
338   DM_DA          *dd = (DM_DA*)da->data;
339   PetscFunctionBegin;
340   PetscValidHeaderSpecific(da,DM_CLASSID,1);
341   dd->adicmf_lfi = admf_lfi;
342   PetscFunctionReturn(0);
343 }
344 
345 /*MC
346        DMDASetLocalAdicFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
347 
348    Synopsis:
349    PetscErrorCode DMDASetLocalAdicFunctionib(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
350 
351    Logically Collective on DMDA
352 
353    Input Parameter:
354 +  da - initial distributed array
355 -  ad_lfi - the local function as computed by ADIC/ADIFOR
356 
357    Level: intermediate
358 
359 .keywords:  distributed array, refine
360 
361 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
362           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
363 M*/
364 
365 #undef __FUNCT__
366 #define __FUNCT__ "DMDASetLocalAdicFunctionib_Private"
367 PetscErrorCode DMDASetLocalAdicFunctionib_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
368 {
369   DM_DA          *dd = (DM_DA*)da->data;
370   PetscFunctionBegin;
371   PetscValidHeaderSpecific(da,DM_CLASSID,1);
372   dd->adic_lfib = ad_lfi;
373   PetscFunctionReturn(0);
374 }
375 
376 /*MC
377        DMDASetLocalAdicMFFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
378 
379    Synopsis:
380    PetscErrorCode  DMDASetLocalAdicFunctionib(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
381 
382    Logically Collective on DMDA
383 
384    Input Parameter:
385 +  da - initial distributed array
386 -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
387 
388    Level: intermediate
389 
390 .keywords:  distributed array, refine
391 
392 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
393           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
394 M*/
395 
396 #undef __FUNCT__
397 #define __FUNCT__ "DMDASetLocalAdicMFFunctionib_Private"
398 PetscErrorCode DMDASetLocalAdicMFFunctionib_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
399 {
400   DM_DA          *dd = (DM_DA*)da->data;
401   PetscFunctionBegin;
402   PetscValidHeaderSpecific(da,DM_CLASSID,1);
403   dd->adicmf_lfib = admf_lfi;
404   PetscFunctionReturn(0);
405 }
406 
407 /*MC
408        DMDASetLocalAdicMFFunction - Caches in a DM a local function computed by ADIC/ADIFOR
409 
410    Synopsis:
411    PetscErrorCode DMDASetLocalAdicMFFunction(DM da,DMDALocalFunction1 ad_lf)
412 
413    Logically Collective on DMDA
414 
415    Input Parameter:
416 +  da - initial distributed array
417 -  ad_lf - the local function as computed by ADIC/ADIFOR
418 
419    Level: intermediate
420 
421 .keywords:  distributed array, refine
422 
423 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
424           DMDASetLocalJacobian()
425 M*/
426 
427 #undef __FUNCT__
428 #define __FUNCT__ "DMDASetLocalAdicMFFunction_Private"
429 PetscErrorCode DMDASetLocalAdicMFFunction_Private(DM da,DMDALocalFunction1 ad_lf)
430 {
431   DM_DA          *dd = (DM_DA*)da->data;
432   PetscFunctionBegin;
433   PetscValidHeaderSpecific(da,DM_CLASSID,1);
434   dd->adicmf_lf = ad_lf;
435   PetscFunctionReturn(0);
436 }
437 
438 #undef __FUNCT__
439 #define __FUNCT__ "DMDAJacobianDefaultLocal"
440 PetscErrorCode DMDAJacobianLocal(DM dm,Vec x,Mat A,Mat B, MatStructure *str)
441 {
442   PetscErrorCode ierr;
443   Vec            localX;
444 
445   PetscFunctionBegin;
446   ierr = DMGetLocalVector(dm,&localX);CHKERRQ(ierr);
447   ierr = DMGlobalToLocalBegin(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
448   ierr = DMGlobalToLocalEnd(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
449   ierr = MatFDColoringApply(B,dm->fd,localX,str,dm);CHKERRQ(ierr);
450   ierr = DMRestoreLocalVector(dm,&localX);CHKERRQ(ierr);
451   /* Assemble true Jacobian; if it is different */
452   if (A != B) {
453     ierr  = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
454     ierr  = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
455   }
456   ierr  = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
457   *str = SAME_NONZERO_PATTERN;
458   PetscFunctionReturn(0);
459 }
460 
461 
462 #undef __FUNCT__
463 #define __FUNCT__ "DMDAJacobian"
464 static PetscErrorCode DMDAJacobian(DM dm,Vec x,Mat A,Mat B, MatStructure *str)
465 {
466   PetscErrorCode ierr;
467   Vec            localX;
468 
469   PetscFunctionBegin;
470   ierr = DMGetLocalVector(dm,&localX);CHKERRQ(ierr);
471   ierr = DMGlobalToLocalBegin(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
472   ierr = DMGlobalToLocalEnd(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
473   ierr = DMDAComputeJacobian1(dm,localX,B,dm->ctx);CHKERRQ(ierr);
474   ierr = DMRestoreLocalVector(dm,&localX);CHKERRQ(ierr);
475   /* Assemble true Jacobian; if it is different */
476   if (A != B) {
477     ierr  = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
478     ierr  = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
479   }
480   ierr  = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
481   *str = SAME_NONZERO_PATTERN;
482   PetscFunctionReturn(0);
483 }
484 
485 /*@C
486        DMDASetLocalJacobian - Caches in a DM a local Jacobian computation function
487 
488    Logically Collective on DMDA
489 
490 
491    Input Parameter:
492 +  da - initial distributed array
493 -  lj - the local Jacobian
494 
495    Level: intermediate
496 
497 .keywords:  distributed array, refine
498 
499 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
500 @*/
501 #undef __FUNCT__
502 #define __FUNCT__ "DMDASetLocalJacobian"
503 PetscErrorCode  DMDASetLocalJacobian(DM da,DMDALocalFunction1 lj)
504 {
505   PetscErrorCode ierr;
506   DM_DA          *dd = (DM_DA*)da->data;
507 
508   PetscFunctionBegin;
509   PetscValidHeaderSpecific(da,DM_CLASSID,1);
510   ierr = DMSetJacobian(da,DMDAJacobian);CHKERRQ(ierr);
511   dd->lj    = lj;
512   PetscFunctionReturn(0);
513 }
514 
515 #undef __FUNCT__
516 #define __FUNCT__ "DMDAGetLocalFunction"
517 /*@C
518        DMDAGetLocalFunction - Gets from a DM a local function and its ADIC/ADIFOR Jacobian
519 
520    Note Collective
521 
522    Input Parameter:
523 .  da - initial distributed array
524 
525    Output Parameter:
526 .  lf - the local function
527 
528    Level: intermediate
529 
530 .keywords:  distributed array, refine
531 
532 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalJacobian(), DMDASetLocalFunction()
533 @*/
534 PetscErrorCode  DMDAGetLocalFunction(DM da,DMDALocalFunction1 *lf)
535 {
536   DM_DA *dd = (DM_DA*)da->data;
537   PetscFunctionBegin;
538   PetscValidHeaderSpecific(da,DM_CLASSID,1);
539   if (lf) *lf = dd->lf;
540   PetscFunctionReturn(0);
541 }
542 
543 #undef __FUNCT__
544 #define __FUNCT__ "DMDAGetLocalJacobian"
545 /*@C
546        DMDAGetLocalJacobian - Gets from a DM a local jacobian
547 
548    Not Collective
549 
550    Input Parameter:
551 .  da - initial distributed array
552 
553    Output Parameter:
554 .  lj - the local jacobian
555 
556    Level: intermediate
557 
558 .keywords:  distributed array, refine
559 
560 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalJacobian()
561 @*/
562 PetscErrorCode  DMDAGetLocalJacobian(DM da,DMDALocalFunction1 *lj)
563 {
564   DM_DA *dd = (DM_DA*)da->data;
565   PetscFunctionBegin;
566   PetscValidHeaderSpecific(da,DM_CLASSID,1);
567   if (lj) *lj = dd->lj;
568   PetscFunctionReturn(0);
569 }
570 
571 #undef __FUNCT__
572 #define __FUNCT__ "DMDAComputeFunction"
573 /*@
574     DMDAComputeFunction - Evaluates a user provided function on each processor that
575         share a DMDA
576 
577    Input Parameters:
578 +    da - the DM that defines the grid
579 .    vu - input vector
580 .    vfu - output vector
581 -    w - any user data
582 
583     Notes: Does NOT do ghost updates on vu upon entry
584 
585            This should eventually replace DMDAComputeFunction1
586 
587     Level: advanced
588 
589 .seealso: DMDAComputeJacobian1WithAdic()
590 
591 @*/
592 PetscErrorCode  DMDAComputeFunction(DM da,PetscErrorCode (*lf)(void),Vec vu,Vec vfu,void *w)
593 {
594   PetscErrorCode ierr;
595   void           *u,*fu;
596   DMDALocalInfo  info;
597   PetscErrorCode (*f)(DMDALocalInfo*,void*,void*,void*) = (PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))lf;
598 
599   PetscFunctionBegin;
600   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
601   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
602   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
603 
604   ierr = (*f)(&info,u,fu,w);CHKERRQ(ierr);
605 
606   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
607   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
608   PetscFunctionReturn(0);
609 }
610 
611 #undef __FUNCT__
612 #define __FUNCT__ "DMDAComputeFunctionLocal"
613 /*@C
614    DMDAComputeFunctionLocal - This is a universal function evaluation routine for
615    a local DM function.
616 
617    Collective on DMDA
618 
619    Input Parameters:
620 +  da - the DM context
621 .  func - The local function
622 .  X - input vector
623 .  F - function vector
624 -  ctx - A user context
625 
626    Level: intermediate
627 
628 .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
629           SNESSetFunction(), SNESSetJacobian()
630 
631 @*/
632 PetscErrorCode  DMDAComputeFunctionLocal(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
633 {
634   Vec            localX;
635   DMDALocalInfo  info;
636   void           *u;
637   void           *fu;
638   PetscErrorCode ierr;
639 
640   PetscFunctionBegin;
641   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
642   /*
643      Scatter ghost points to local vector, using the 2-step process
644         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
645   */
646   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
647   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
648   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
649   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
650   ierr = DMDAVecGetArray(da,F,&fu);CHKERRQ(ierr);
651   ierr = (*func)(&info,u,fu,ctx);CHKERRQ(ierr);
652   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
653   ierr = DMDAVecRestoreArray(da,F,&fu);CHKERRQ(ierr);
654   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
655   PetscFunctionReturn(0);
656 }
657 
658 #undef __FUNCT__
659 #define __FUNCT__ "DMDAComputeFunctionLocalGhost"
660 /*@C
661    DMDAComputeFunctionLocalGhost - This is a universal function evaluation routine for
662    a local DM function, but the ghost values of the output are communicated and added.
663 
664    Collective on DMDA
665 
666    Input Parameters:
667 +  da - the DM context
668 .  func - The local function
669 .  X - input vector
670 .  F - function vector
671 -  ctx - A user context
672 
673    Level: intermediate
674 
675 .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
676           SNESSetFunction(), SNESSetJacobian()
677 
678 @*/
679 PetscErrorCode  DMDAComputeFunctionLocalGhost(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
680 {
681   Vec            localX, localF;
682   DMDALocalInfo  info;
683   void           *u;
684   void           *fu;
685   PetscErrorCode ierr;
686 
687   PetscFunctionBegin;
688   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
689   ierr = DMGetLocalVector(da,&localF);CHKERRQ(ierr);
690   /*
691      Scatter ghost points to local vector, using the 2-step process
692         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
693   */
694   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
695   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
696   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
697   ierr = VecSet(localF, 0.0);CHKERRQ(ierr);
698   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
699   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
700   ierr = DMDAVecGetArray(da,localF,&fu);CHKERRQ(ierr);
701   ierr = (*func)(&info,u,fu,ctx);CHKERRQ(ierr);
702   ierr = DMLocalToGlobalBegin(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
703   ierr = DMLocalToGlobalEnd(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
704   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
705   ierr = DMDAVecRestoreArray(da,localF,&fu);CHKERRQ(ierr);
706   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
707   ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr);
708   PetscFunctionReturn(0);
709 }
710 
711 #undef __FUNCT__
712 #define __FUNCT__ "DMDAComputeFunction1"
713 /*@
714     DMDAComputeFunction1 - Evaluates a user provided function on each processor that
715         share a DMDA
716 
717    Input Parameters:
718 +    da - the DM that defines the grid
719 .    vu - input vector
720 .    vfu - output vector
721 -    w - any user data
722 
723     Notes: Does NOT do ghost updates on vu upon entry
724 
725     Level: advanced
726 
727 .seealso: DMDAComputeJacobian1WithAdic()
728 
729 @*/
730 PetscErrorCode  DMDAComputeFunction1(DM da,Vec vu,Vec vfu,void *w)
731 {
732   PetscErrorCode ierr;
733   void           *u,*fu;
734   DMDALocalInfo  info;
735   DM_DA          *dd = (DM_DA*)da->data;
736 
737   PetscFunctionBegin;
738   if (!dd->lf) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_NULL,"DMDASetLocalFunction() never called");
739   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
740   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
741   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
742 
743   CHKMEMQ;
744   ierr = (*dd->lf)(&info,u,fu,w);CHKERRQ(ierr);
745   CHKMEMQ;
746 
747   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
748   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
749   PetscFunctionReturn(0);
750 }
751 
752 #undef __FUNCT__
753 #define __FUNCT__ "DMDAComputeFunctioniTest1"
754 PetscErrorCode  DMDAComputeFunctioniTest1(DM da,void *w)
755 {
756   Vec            vu,fu,fui;
757   PetscErrorCode ierr;
758   PetscInt       i,n;
759   PetscScalar    *ui;
760   PetscRandom    rnd;
761   PetscReal      norm;
762 
763   PetscFunctionBegin;
764   ierr = DMGetLocalVector(da,&vu);CHKERRQ(ierr);
765   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rnd);CHKERRQ(ierr);
766   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
767   ierr = VecSetRandom(vu,rnd);CHKERRQ(ierr);
768   ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr);
769 
770   ierr = DMGetGlobalVector(da,&fu);CHKERRQ(ierr);
771   ierr = DMGetGlobalVector(da,&fui);CHKERRQ(ierr);
772 
773   ierr = DMDAComputeFunction1(da,vu,fu,w);CHKERRQ(ierr);
774 
775   ierr = VecGetArray(fui,&ui);CHKERRQ(ierr);
776   ierr = VecGetLocalSize(fui,&n);CHKERRQ(ierr);
777   for (i=0; i<n; i++) {
778     ierr = DMDAComputeFunctioni1(da,i,vu,ui+i,w);CHKERRQ(ierr);
779   }
780   ierr = VecRestoreArray(fui,&ui);CHKERRQ(ierr);
781 
782   ierr = VecAXPY(fui,-1.0,fu);CHKERRQ(ierr);
783   ierr = VecNorm(fui,NORM_2,&norm);CHKERRQ(ierr);
784   ierr = PetscPrintf(((PetscObject)da)->comm,"Norm of difference in vectors %G\n",norm);CHKERRQ(ierr);
785   ierr = VecView(fu,0);CHKERRQ(ierr);
786   ierr = VecView(fui,0);CHKERRQ(ierr);
787 
788   ierr = DMRestoreLocalVector(da,&vu);CHKERRQ(ierr);
789   ierr = DMRestoreGlobalVector(da,&fu);CHKERRQ(ierr);
790   ierr = DMRestoreGlobalVector(da,&fui);CHKERRQ(ierr);
791   PetscFunctionReturn(0);
792 }
793 
794 #undef __FUNCT__
795 #define __FUNCT__ "DMDAComputeFunctioni1"
796 /*@
797     DMDAComputeFunctioni1 - Evaluates a user provided point-wise function
798 
799    Input Parameters:
800 +    da - the DM that defines the grid
801 .    i - the component of the function we wish to compute (must be local)
802 .    vu - input vector
803 .    vfu - output value
804 -    w - any user data
805 
806     Notes: Does NOT do ghost updates on vu upon entry
807 
808     Level: advanced
809 
810 .seealso: DMDAComputeJacobian1WithAdic()
811 
812 @*/
813 PetscErrorCode  DMDAComputeFunctioni1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
814 {
815   PetscErrorCode ierr;
816   void           *u;
817   DMDALocalInfo  info;
818   MatStencil     stencil;
819   DM_DA          *dd = (DM_DA*)da->data;
820 
821   PetscFunctionBegin;
822 
823   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
824   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
825 
826   /* figure out stencil value from i */
827   stencil.c = i % info.dof;
828   stencil.i = (i % (info.xm*info.dof))/info.dof;
829   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
830   stencil.k = i/(info.xm*info.ym*info.dof);
831 
832   ierr = (*dd->lfi)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
833 
834   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
835   PetscFunctionReturn(0);
836 }
837 
838 #undef __FUNCT__
839 #define __FUNCT__ "DMDAComputeFunctionib1"
840 /*@
841     DMDAComputeFunctionib1 - Evaluates a user provided point-block function
842 
843    Input Parameters:
844 +    da - the DM that defines the grid
845 .    i - the component of the function we wish to compute (must be local)
846 .    vu - input vector
847 .    vfu - output value
848 -    w - any user data
849 
850     Notes: Does NOT do ghost updates on vu upon entry
851 
852     Level: advanced
853 
854 .seealso: DMDAComputeJacobian1WithAdic()
855 
856 @*/
857 PetscErrorCode  DMDAComputeFunctionib1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
858 {
859   PetscErrorCode ierr;
860   void           *u;
861   DMDALocalInfo  info;
862   MatStencil     stencil;
863   DM_DA          *dd = (DM_DA*)da->data;
864 
865   PetscFunctionBegin;
866   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
867   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
868 
869   /* figure out stencil value from i */
870   stencil.c = i % info.dof;
871   if (stencil.c) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Point-block functions can only be called for the entire block");
872   stencil.i = (i % (info.xm*info.dof))/info.dof;
873   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
874   stencil.k = i/(info.xm*info.ym*info.dof);
875 
876   ierr = (*dd->lfib)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
877 
878   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
879   PetscFunctionReturn(0);
880 }
881 
882 #if defined(new)
883 #undef __FUNCT__
884 #define __FUNCT__ "DMDAGetDiagonal_MFFD"
885 /*
886   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
887     function lives on a DMDA
888 
889         y ~= (F(u + ha) - F(u))/h,
890   where F = nonlinear function, as set by SNESSetFunction()
891         u = current iterate
892         h = difference interval
893 */
894 PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
895 {
896   PetscScalar    h,*aa,*ww,v;
897   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
898   PetscErrorCode ierr;
899   PetscInt       gI,nI;
900   MatStencil     stencil;
901   DMDALocalInfo  info;
902 
903   PetscFunctionBegin;
904   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
905   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
906 
907   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
908   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
909 
910   nI = 0;
911     h  = ww[gI];
912     if (h == 0.0) h = 1.0;
913 #if !defined(PETSC_USE_COMPLEX)
914     if (h < umin && h >= 0.0)      h = umin;
915     else if (h < 0.0 && h > -umin) h = -umin;
916 #else
917     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
918     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
919 #endif
920     h     *= epsilon;
921 
922     ww[gI] += h;
923     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
924     aa[nI]  = (v - aa[nI])/h;
925     ww[gI] -= h;
926     nI++;
927   }
928   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
929   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
930   PetscFunctionReturn(0);
931 }
932 #endif
933 
934 #if defined(PETSC_HAVE_ADIC)
935 EXTERN_C_BEGIN
936 #include <adic/ad_utils.h>
937 EXTERN_C_END
938 
939 #undef __FUNCT__
940 #define __FUNCT__ "DMDAComputeJacobian1WithAdic"
941 /*@C
942     DMDAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that
943         share a DMDA
944 
945    Input Parameters:
946 +    da - the DM that defines the grid
947 .    vu - input vector (ghosted)
948 .    J - output matrix
949 -    w - any user data
950 
951    Level: advanced
952 
953     Notes: Does NOT do ghost updates on vu upon entry
954 
955 .seealso: DMDAComputeFunction1()
956 
957 @*/
958 PetscErrorCode  DMDAComputeJacobian1WithAdic(DM da,Vec vu,Mat J,void *w)
959 {
960   PetscErrorCode ierr;
961   PetscInt       gtdof,tdof;
962   PetscScalar    *ustart;
963   DMDALocalInfo  info;
964   void           *ad_u,*ad_f,*ad_ustart,*ad_fstart;
965   ISColoring     iscoloring;
966 
967   PetscFunctionBegin;
968   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
969 
970   PetscADResetIndep();
971 
972   /* get space for derivative objects.  */
973   ierr = DMDAGetAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
974   ierr = DMDAGetAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
975   ierr = VecGetArray(vu,&ustart);CHKERRQ(ierr);
976   ierr = DMCreateColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
977 
978   PetscADSetValueAndColor(ad_ustart,gtdof,iscoloring->colors,ustart);
979 
980   ierr = VecRestoreArray(vu,&ustart);CHKERRQ(ierr);
981   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
982   ierr = PetscADIncrementTotalGradSize(iscoloring->n);CHKERRQ(ierr);
983   PetscADSetIndepDone();
984 
985   ierr = PetscLogEventBegin(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
986   ierr = (*dd->adic_lf)(&info,ad_u,ad_f,w);CHKERRQ(ierr);
987   ierr = PetscLogEventEnd(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
988 
989   /* stick the values into the matrix */
990   ierr = MatSetValuesAdic(J,(PetscScalar**)ad_fstart);CHKERRQ(ierr);
991 
992   /* return space for derivative objects.  */
993   ierr = DMDARestoreAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
994   ierr = DMDARestoreAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
995   PetscFunctionReturn(0);
996 }
997 
998 #undef __FUNCT__
999 #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdic"
1000 /*@C
1001     DMDAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on
1002     each processor that shares a DMDA.
1003 
1004     Input Parameters:
1005 +   da - the DM that defines the grid
1006 .   vu - Jacobian is computed at this point (ghosted)
1007 .   v - product is done on this vector (ghosted)
1008 .   fu - output vector = J(vu)*v (not ghosted)
1009 -   w - any user data
1010 
1011     Notes:
1012     This routine does NOT do ghost updates on vu upon entry.
1013 
1014    Level: advanced
1015 
1016 .seealso: DMDAComputeFunction1()
1017 
1018 @*/
1019 PetscErrorCode  DMDAMultiplyByJacobian1WithAdic(DM da,Vec vu,Vec v,Vec f,void *w)
1020 {
1021   PetscErrorCode ierr;
1022   PetscInt       i,gtdof,tdof;
1023   PetscScalar    *avu,*av,*af,*ad_vustart,*ad_fstart;
1024   DMDALocalInfo  info;
1025   void           *ad_vu,*ad_f;
1026 
1027   PetscFunctionBegin;
1028   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1029 
1030   /* get space for derivative objects.  */
1031   ierr = DMDAGetAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
1032   ierr = DMDAGetAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
1033 
1034   /* copy input vector into derivative object */
1035   ierr = VecGetArray(vu,&avu);CHKERRQ(ierr);
1036   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
1037   for (i=0; i<gtdof; i++) {
1038     ad_vustart[2*i]   = avu[i];
1039     ad_vustart[2*i+1] = av[i];
1040   }
1041   ierr = VecRestoreArray(vu,&avu);CHKERRQ(ierr);
1042   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
1043 
1044   PetscADResetIndep();
1045   ierr = PetscADIncrementTotalGradSize(1);CHKERRQ(ierr);
1046   PetscADSetIndepDone();
1047 
1048   ierr = (*dd->adicmf_lf)(&info,ad_vu,ad_f,w);CHKERRQ(ierr);
1049 
1050   /* stick the values into the vector */
1051   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
1052   for (i=0; i<tdof; i++) {
1053     af[i] = ad_fstart[2*i+1];
1054   }
1055   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
1056 
1057   /* return space for derivative objects.  */
1058   ierr = DMDARestoreAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
1059   ierr = DMDARestoreAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
1060   PetscFunctionReturn(0);
1061 }
1062 #endif
1063 
1064 #undef __FUNCT__
1065 #define __FUNCT__ "DMDAComputeJacobian1"
1066 /*@
1067     DMDAComputeJacobian1 - Evaluates a local Jacobian function on each processor that
1068         share a DMDA
1069 
1070    Input Parameters:
1071 +    da - the DM that defines the grid
1072 .    vu - input vector (ghosted)
1073 .    J - output matrix
1074 -    w - any user data
1075 
1076     Notes: Does NOT do ghost updates on vu upon entry
1077 
1078     Level: advanced
1079 
1080 .seealso: DMDAComputeFunction1()
1081 
1082 @*/
1083 PetscErrorCode  DMDAComputeJacobian1(DM da,Vec vu,Mat J,void *w)
1084 {
1085   PetscErrorCode ierr;
1086   void           *u;
1087   DMDALocalInfo  info;
1088   DM_DA          *dd = (DM_DA*)da->data;
1089 
1090   PetscFunctionBegin;
1091   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1092   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
1093   ierr = (*dd->lj)(&info,u,J,w);CHKERRQ(ierr);
1094   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 
1099 #undef __FUNCT__
1100 #define __FUNCT__ "DMDAComputeJacobian1WithAdifor"
1101 /*
1102     DMDAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that
1103         share a DMDA
1104 
1105    Input Parameters:
1106 +    da - the DM that defines the grid
1107 .    vu - input vector (ghosted)
1108 .    J - output matrix
1109 -    w - any user data
1110 
1111     Notes: Does NOT do ghost updates on vu upon entry
1112 
1113 .seealso: DMDAComputeFunction1()
1114 
1115 */
1116 PetscErrorCode  DMDAComputeJacobian1WithAdifor(DM da,Vec vu,Mat J,void *w)
1117 {
1118   PetscErrorCode  ierr;
1119   PetscInt        i,Nc,N;
1120   ISColoringValue *color;
1121   DMDALocalInfo   info;
1122   PetscScalar     *u,*g_u,*g_f,*f = 0,*p_u;
1123   ISColoring      iscoloring;
1124   DM_DA          *dd = (DM_DA*)da->data;
1125   void            (*lf)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*) =
1126                   (void (*)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*))*dd->adifor_lf;
1127 
1128   PetscFunctionBegin;
1129   ierr = DMCreateColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
1130   Nc   = iscoloring->n;
1131   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1132   N    = info.gxm*info.gym*info.gzm*info.dof;
1133 
1134   /* get space for derivative objects.  */
1135   ierr  = PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);CHKERRQ(ierr);
1136   ierr  = PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));CHKERRQ(ierr);
1137   p_u   = g_u;
1138   color = iscoloring->colors;
1139   for (i=0; i<N; i++) {
1140     p_u[*color++] = 1.0;
1141     p_u          += Nc;
1142   }
1143   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
1144   ierr = PetscMalloc2(Nc*info.xm*info.ym*info.zm*info.dof,PetscScalar,&g_f,info.xm*info.ym*info.zm*info.dof,PetscScalar,&f);CHKERRQ(ierr);
1145 
1146   /* Seed the input array g_u with coloring information */
1147 
1148   ierr = VecGetArray(vu,&u);CHKERRQ(ierr);
1149   (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);CHKERRQ(ierr);
1150   ierr = VecRestoreArray(vu,&u);CHKERRQ(ierr);
1151 
1152   /* stick the values into the matrix */
1153   /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */
1154   ierr = MatSetValuesAdifor(J,Nc,g_f);CHKERRQ(ierr);
1155 
1156   /* return space for derivative objects.  */
1157   ierr = PetscFree(g_u);CHKERRQ(ierr);
1158   ierr = PetscFree2(g_f,f);CHKERRQ(ierr);
1159   PetscFunctionReturn(0);
1160 }
1161 
1162 #undef __FUNCT__
1163 #define __FUNCT__ "DMDAFormJacobianLocal"
1164 /*@C
1165    DMDAFormjacobianLocal - This is a universal Jacobian evaluation routine for
1166    a local DM function.
1167 
1168    Collective on DMDA
1169 
1170    Input Parameters:
1171 +  da - the DM context
1172 .  func - The local function
1173 .  X - input vector
1174 .  J - Jacobian matrix
1175 -  ctx - A user context
1176 
1177    Level: intermediate
1178 
1179 .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
1180           SNESSetFunction(), SNESSetJacobian()
1181 
1182 @*/
1183 PetscErrorCode  DMDAFormJacobianLocal(DM da, DMDALocalFunction1 func, Vec X, Mat J, void *ctx)
1184 {
1185   Vec            localX;
1186   DMDALocalInfo  info;
1187   void           *u;
1188   PetscErrorCode ierr;
1189 
1190   PetscFunctionBegin;
1191   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
1192   /*
1193      Scatter ghost points to local vector, using the 2-step process
1194         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
1195   */
1196   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1197   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1198   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1199   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
1200   ierr = (*func)(&info,u,J,ctx);CHKERRQ(ierr);
1201   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
1202   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
1203   PetscFunctionReturn(0);
1204 }
1205 
1206 #undef __FUNCT__
1207 #define __FUNCT__ "DMDAMultiplyByJacobian1WithAD"
1208 /*@C
1209     DMDAMultiplyByJacobian1WithAD - Applies a Jacobian function supplied by ADIFOR or ADIC
1210     to a vector on each processor that shares a DMDA.
1211 
1212    Input Parameters:
1213 +    da - the DM that defines the grid
1214 .    vu - Jacobian is computed at this point (ghosted)
1215 .    v - product is done on this vector (ghosted)
1216 .    fu - output vector = J(vu)*v (not ghosted)
1217 -    w - any user data
1218 
1219     Notes:
1220     This routine does NOT do ghost updates on vu and v upon entry.
1221 
1222     Automatically calls DMDAMultiplyByJacobian1WithAdifor() or DMDAMultiplyByJacobian1WithAdic()
1223     depending on whether DMDASetLocalAdicMFFunction() or DMDASetLocalAdiforMFFunction() was called.
1224 
1225    Level: advanced
1226 
1227 .seealso: DMDAComputeFunction1(), DMDAMultiplyByJacobian1WithAdifor(), DMDAMultiplyByJacobian1WithAdic()
1228 
1229 @*/
1230 PetscErrorCode  DMDAMultiplyByJacobian1WithAD(DM da,Vec u,Vec v,Vec f,void *w)
1231 {
1232   PetscErrorCode ierr;
1233   DM_DA          *dd = (DM_DA*)da->data;
1234 
1235   PetscFunctionBegin;
1236   if (dd->adicmf_lf) {
1237 #if defined(PETSC_HAVE_ADIC)
1238     ierr = DMDAMultiplyByJacobian1WithAdic(da,u,v,f,w);CHKERRQ(ierr);
1239 #else
1240     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP_SYS,"Requires ADIC to be installed and cannot use complex numbers");
1241 #endif
1242   } else if (dd->adiformf_lf) {
1243     ierr = DMDAMultiplyByJacobian1WithAdifor(da,u,v,f,w);CHKERRQ(ierr);
1244   } else SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ORDER,"Must call DMDASetLocalAdiforMFFunction() or DMDASetLocalAdicMFFunction() before using");
1245   PetscFunctionReturn(0);
1246 }
1247 
1248 
1249 #undef __FUNCT__
1250 #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdifor"
1251 /*@C
1252     DMDAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that
1253         share a DM to a vector
1254 
1255    Input Parameters:
1256 +    da - the DM that defines the grid
1257 .    vu - Jacobian is computed at this point (ghosted)
1258 .    v - product is done on this vector (ghosted)
1259 .    fu - output vector = J(vu)*v (not ghosted)
1260 -    w - any user data
1261 
1262     Notes: Does NOT do ghost updates on vu and v upon entry
1263 
1264    Level: advanced
1265 
1266 .seealso: DMDAComputeFunction1()
1267 
1268 @*/
1269 PetscErrorCode  DMDAMultiplyByJacobian1WithAdifor(DM da,Vec u,Vec v,Vec f,void *w)
1270 {
1271   PetscErrorCode ierr;
1272   PetscScalar    *au,*av,*af,*awork;
1273   Vec            work;
1274   DMDALocalInfo  info;
1275   DM_DA          *dd = (DM_DA*)da->data;
1276   void           (*lf)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*) =
1277                  (void (*)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*))*dd->adiformf_lf;
1278 
1279   PetscFunctionBegin;
1280   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1281 
1282   ierr = DMGetGlobalVector(da,&work);CHKERRQ(ierr);
1283   ierr = VecGetArray(u,&au);CHKERRQ(ierr);
1284   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
1285   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
1286   ierr = VecGetArray(work,&awork);CHKERRQ(ierr);
1287   (lf)(&info,au,av,awork,af,w,&ierr);CHKERRQ(ierr);
1288   ierr = VecRestoreArray(u,&au);CHKERRQ(ierr);
1289   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
1290   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
1291   ierr = VecRestoreArray(work,&awork);CHKERRQ(ierr);
1292   ierr = DMRestoreGlobalVector(da,&work);CHKERRQ(ierr);
1293 
1294   PetscFunctionReturn(0);
1295 }
1296 
1297 #undef __FUNCT__
1298 #define __FUNCT__ "DMSetUp_DA_2D"
1299 PetscErrorCode  DMSetUp_DA_2D(DM da)
1300 {
1301   DM_DA                  *dd = (DM_DA*)da->data;
1302   const PetscInt         M            = dd->M;
1303   const PetscInt         N            = dd->N;
1304   PetscInt               m            = dd->m;
1305   PetscInt               n            = dd->n;
1306   PetscInt               o            = dd->overlap;
1307   const PetscInt         dof          = dd->w;
1308   const PetscInt         s            = dd->s;
1309   const DMDABoundaryType bx         = dd->bx;
1310   const DMDABoundaryType by         = dd->by;
1311   const DMDAStencilType  stencil_type = dd->stencil_type;
1312   PetscInt               *lx           = dd->lx;
1313   PetscInt               *ly           = dd->ly;
1314   MPI_Comm               comm;
1315   PetscMPIInt            rank,size;
1316   PetscInt               xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end,IXs,IXe,IYs,IYe;
1317   PetscInt               up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn,*idx_cpy;
1318   const PetscInt         *idx_full;
1319   PetscInt               xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
1320   PetscInt               s_x,s_y; /* s proportionalized to w */
1321   PetscInt               sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
1322   Vec                    local,global;
1323   VecScatter             ltog,gtol;
1324   IS                     to,from,ltogis;
1325   PetscErrorCode         ierr;
1326 
1327   PetscFunctionBegin;
1328   if (dof < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
1329   if (s < 0) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
1330   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
1331 #if !defined(PETSC_USE_64BIT_INDICES)
1332   if (((Petsc64bitInt) M)*((Petsc64bitInt) N)*((Petsc64bitInt) dof) > (Petsc64bitInt) PETSC_MPI_INT_MAX) SETERRQ3(comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %D by %D by %D (dof) is too large for 32 bit indices",M,N,dof);
1333 #endif
1334 
1335   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
1336   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
1337 
1338   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1339   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1340 
1341   dd->dim         = 2;
1342   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
1343   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
1344 
1345   if (m != PETSC_DECIDE) {
1346     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
1347     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
1348   }
1349   if (n != PETSC_DECIDE) {
1350     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
1351     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
1352   }
1353 
1354   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
1355     if (n != PETSC_DECIDE) {
1356       m = size/n;
1357     } else if (m != PETSC_DECIDE) {
1358       n = size/m;
1359     } else {
1360       /* try for squarish distribution */
1361       m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
1362       if (!m) m = 1;
1363       while (m > 0) {
1364 	n = size/m;
1365 	if (m*n == size) break;
1366 	m--;
1367       }
1368       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
1369     }
1370     if (m*n != size) SETERRQ(comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n ");
1371   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
1372 
1373   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
1374   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
1375 
1376   /*
1377      Determine locally owned region
1378      xs is the first local node number, x is the number of local nodes
1379   */
1380   if (!lx) {
1381     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
1382     lx = dd->lx;
1383     for (i=0; i<m; i++) {
1384       lx[i] = M/m + ((M % m) > i);
1385     }
1386   }
1387   x  = lx[rank % m];
1388   xs = 0;
1389   for (i=0; i<(rank % m); i++) {
1390     xs += lx[i];
1391   }
1392 #if defined(PETSC_USE_DEBUG)
1393   left = xs;
1394   for (i=(rank % m); i<m; i++) {
1395     left += lx[i];
1396   }
1397   if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M);
1398 #endif
1399 
1400   /*
1401      Determine locally owned region
1402      ys is the first local node number, y is the number of local nodes
1403   */
1404   if (!ly) {
1405     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
1406     ly = dd->ly;
1407     for (i=0; i<n; i++) {
1408       ly[i] = N/n + ((N % n) > i);
1409     }
1410   }
1411   y  = ly[rank/m];
1412   ys = 0;
1413   for (i=0; i<(rank/m); i++) {
1414     ys += ly[i];
1415   }
1416 #if defined(PETSC_USE_DEBUG)
1417   left = ys;
1418   for (i=(rank/m); i<n; i++) {
1419     left += ly[i];
1420   }
1421   if (left != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N);
1422 #endif
1423 
1424   /*
1425    check if the scatter requires more than one process neighbor or wraps around
1426    the domain more than once
1427   */
1428   if ((x < s+o) && ((m > 1) || (bx == DMDA_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s+o);
1429   if ((y < s+o) && ((n > 1) || (by == DMDA_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s+o);
1430   xe = xs + x;
1431   ye = ys + y;
1432 
1433   /* determine ghost region (Xs) and region scattered into (IXs)  */
1434   if (xs-s-o > 0) {
1435     Xs = xs - s - o; IXs = xs - s - o;
1436   } else {
1437     if (bx) {
1438       Xs = xs - s;
1439     } else {
1440       Xs = 0;
1441     }
1442     IXs = 0;
1443   }
1444   if (xe+s+o <= M) {
1445     Xe = xe + s + o; IXe = xe + s + o;
1446   } else {
1447     if (bx) {
1448       Xs = xs - s - o; Xe = xe + s;
1449     } else {
1450       Xe = M;
1451     }
1452     IXe = M;
1453   }
1454 
1455   if (bx == DMDA_BOUNDARY_PERIODIC) {
1456     IXs = xs - s - o;
1457     IXe = xe + s + o;
1458     Xs = xs - s - o;
1459     Xe = xe + s + o;
1460   }
1461 
1462   if (ys-s-o > 0) {
1463     Ys = ys - s - o; IYs = ys - s - o;
1464   } else {
1465     if (by) {
1466       Ys = ys - s;
1467     } else {
1468       Ys = 0;
1469     }
1470     IYs = 0;
1471   }
1472   if (ye+s+o <= N) {
1473     Ye = ye + s + o; IYe = ye + s + o;
1474   } else {
1475     if (by) {
1476       Ye = ye + s;
1477     } else {
1478       Ye = N;
1479     }
1480     IYe = N;
1481   }
1482 
1483   if (by == DMDA_BOUNDARY_PERIODIC) {
1484     IYs = ys - s - o;
1485     IYe = ye + s + o;
1486     Ys = ys - s - o;
1487     Ye = ye + s + o;
1488   }
1489 
1490   /* stencil length in each direction */
1491   s_x = s + o;
1492   s_y = s + o;
1493 
1494   /* determine starting point of each processor */
1495   nn    = x*y;
1496   ierr  = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
1497   ierr  = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
1498   bases[0] = 0;
1499   for (i=1; i<=size; i++) {
1500     bases[i] = ldims[i-1];
1501   }
1502   for (i=1; i<=size; i++) {
1503     bases[i] += bases[i-1];
1504   }
1505   base = bases[rank]*dof;
1506 
1507   /* allocate the base parallel and sequential vectors */
1508   dd->Nlocal = x*y*dof;
1509   ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
1510   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
1511   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr);
1512 
1513   /* generate appropriate vector scatters */
1514   /* local to global inserts non-ghost point region into global */
1515   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
1516   ierr = ISCreateStride(comm,x*y*dof,start,1,&to);CHKERRQ(ierr);
1517 
1518   count = x*y;
1519   ierr = PetscMalloc(x*y*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1520   left = xs - Xs; right = left + x;
1521   down = ys - Ys; up = down + y;
1522   count = 0;
1523   for (i=down; i<up; i++) {
1524     for (j=left; j<right; j++) {
1525       idx[count++] = i*(Xe-Xs) + j;
1526     }
1527   }
1528 
1529   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
1530   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
1531   ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr);
1532   ierr = ISDestroy(&from);CHKERRQ(ierr);
1533   ierr = ISDestroy(&to);CHKERRQ(ierr);
1534 
1535   /* global to local must include ghost points within the domain,
1536      but not ghost points outside the domain that aren't periodic */
1537   if (stencil_type == DMDA_STENCIL_BOX || o > 0) {
1538     count = (IXe-IXs)*(IYe-IYs);
1539     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1540 
1541     left = IXs - Xs; right = left + (IXe-IXs);
1542     down = IYs - Ys; up = down + (IYe-IYs);
1543     count = 0;
1544     for (i=down; i<up; i++) {
1545       for (j=left; j<right; j++) {
1546         idx[count++] = j + i*(Xe-Xs);
1547       }
1548     }
1549     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
1550 
1551   } else {
1552     /* must drop into cross shape region */
1553     /*       ---------|
1554             |  top    |
1555          |---         ---| up
1556          |   middle      |
1557          |               |
1558          ----         ---- down
1559             | bottom  |
1560             -----------
1561          Xs xs        xe Xe */
1562     count = (ys-IYs)*x + y*(IXe-IXs) + (IYe-ye)*x;
1563     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1564 
1565     left = xs - Xs; right = left + x;
1566     down = ys - Ys; up = down + y;
1567     count = 0;
1568     /* bottom */
1569     for (i=(IYs-Ys); i<down; i++) {
1570       for (j=left; j<right; j++) {
1571         idx[count++] = j + i*(Xe-Xs);
1572       }
1573     }
1574     /* middle */
1575     for (i=down; i<up; i++) {
1576       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
1577         idx[count++] = j + i*(Xe-Xs);
1578       }
1579     }
1580     /* top */
1581     for (i=up; i<up+IYe-ye; i++) {
1582       for (j=left; j<right; j++) {
1583         idx[count++] = j + i*(Xe-Xs);
1584       }
1585     }
1586     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
1587   }
1588 
1589 
1590   /* determine who lies on each side of us stored in    n6 n7 n8
1591                                                         n3    n5
1592                                                         n0 n1 n2
1593   */
1594 
1595   /* Assume the Non-Periodic Case */
1596   n1 = rank - m;
1597   if (rank % m) {
1598     n0 = n1 - 1;
1599   } else {
1600     n0 = -1;
1601   }
1602   if ((rank+1) % m) {
1603     n2 = n1 + 1;
1604     n5 = rank + 1;
1605     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
1606   } else {
1607     n2 = -1; n5 = -1; n8 = -1;
1608   }
1609   if (rank % m) {
1610     n3 = rank - 1;
1611     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
1612   } else {
1613     n3 = -1; n6 = -1;
1614   }
1615   n7 = rank + m; if (n7 >= m*n) n7 = -1;
1616 
1617   if (bx == DMDA_BOUNDARY_PERIODIC && by == DMDA_BOUNDARY_PERIODIC) {
1618   /* Modify for Periodic Cases */
1619     /* Handle all four corners */
1620     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
1621     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
1622     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
1623     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
1624 
1625     /* Handle Top and Bottom Sides */
1626     if (n1 < 0) n1 = rank + m * (n-1);
1627     if (n7 < 0) n7 = rank - m * (n-1);
1628     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
1629     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
1630     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
1631     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
1632 
1633     /* Handle Left and Right Sides */
1634     if (n3 < 0) n3 = rank + (m-1);
1635     if (n5 < 0) n5 = rank - (m-1);
1636     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
1637     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
1638     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
1639     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
1640   } else if (by == DMDA_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
1641     if (n1 < 0) n1 = rank + m * (n-1);
1642     if (n7 < 0) n7 = rank - m * (n-1);
1643     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
1644     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
1645     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
1646     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
1647   } else if (bx == DMDA_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
1648     if (n3 < 0) n3 = rank + (m-1);
1649     if (n5 < 0) n5 = rank - (m-1);
1650     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
1651     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
1652     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
1653     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
1654   }
1655 
1656   ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
1657   dd->neighbors[0] = n0;
1658   dd->neighbors[1] = n1;
1659   dd->neighbors[2] = n2;
1660   dd->neighbors[3] = n3;
1661   dd->neighbors[4] = rank;
1662   dd->neighbors[5] = n5;
1663   dd->neighbors[6] = n6;
1664   dd->neighbors[7] = n7;
1665   dd->neighbors[8] = n8;
1666 
1667   if (stencil_type == DMDA_STENCIL_STAR && o == 0) {
1668     /* save corner processor numbers */
1669     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
1670     n0 = n2 = n6 = n8 = -1;
1671   }
1672 
1673   ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1674   ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*sizeof(PetscInt));CHKERRQ(ierr);
1675 
1676   nn = 0;
1677   xbase = bases[rank];
1678   for (i=1; i<=s_y; i++) {
1679     if (n0 >= 0) { /* left below */
1680       x_t = lx[n0 % m];
1681       y_t = ly[(n0/m)];
1682       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
1683       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1684     }
1685     if (n1 >= 0) { /* directly below */
1686       x_t = x;
1687       y_t = ly[(n1/m)];
1688       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
1689       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1690     }
1691     if (n2 >= 0) { /* right below */
1692       x_t = lx[n2 % m];
1693       y_t = ly[(n2/m)];
1694       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
1695       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1696     }
1697   }
1698 
1699   for (i=0; i<y; i++) {
1700     if (n3 >= 0) { /* directly left */
1701       x_t = lx[n3 % m];
1702       /* y_t = y; */
1703       s_t = bases[n3] + (i+1)*x_t - s_x;
1704       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1705     }
1706 
1707     for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
1708 
1709     if (n5 >= 0) { /* directly right */
1710       x_t = lx[n5 % m];
1711       /* y_t = y; */
1712       s_t = bases[n5] + (i)*x_t;
1713       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1714     }
1715   }
1716 
1717   for (i=1; i<=s_y; i++) {
1718     if (n6 >= 0) { /* left above */
1719       x_t = lx[n6 % m];
1720       /* y_t = ly[(n6/m)]; */
1721       s_t = bases[n6] + (i)*x_t - s_x;
1722       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1723     }
1724     if (n7 >= 0) { /* directly above */
1725       x_t = x;
1726       /* y_t = ly[(n7/m)]; */
1727       s_t = bases[n7] + (i-1)*x_t;
1728       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1729     }
1730     if (n8 >= 0) { /* right above */
1731       x_t = lx[n8 % m];
1732       /* y_t = ly[(n8/m)]; */
1733       s_t = bases[n8] + (i-1)*x_t;
1734       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1735     }
1736   }
1737 
1738   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
1739   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
1740   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
1741   ierr = ISDestroy(&to);CHKERRQ(ierr);
1742   ierr = ISDestroy(&from);CHKERRQ(ierr);
1743 
1744   if (stencil_type == DMDA_STENCIL_STAR && o == 0) {
1745     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
1746   }
1747 
1748   if (((stencil_type == DMDA_STENCIL_STAR) ||
1749        (bx && bx != DMDA_BOUNDARY_PERIODIC) ||
1750        (by && by != DMDA_BOUNDARY_PERIODIC)) && o == 0) {
1751     /*
1752         Recompute the local to global mappings, this time keeping the
1753       information about the cross corner processor numbers and any ghosted
1754       but not periodic indices.
1755     */
1756     nn = 0;
1757     xbase = bases[rank];
1758     for (i=1; i<=s_y; i++) {
1759       if (n0 >= 0) { /* left below */
1760         x_t = lx[n0 % m];
1761         y_t = ly[(n0/m)];
1762         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
1763         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1764       } else if (xs-Xs > 0 && ys-Ys > 0) {
1765         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
1766       }
1767       if (n1 >= 0) { /* directly below */
1768         x_t = x;
1769         y_t = ly[(n1/m)];
1770         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
1771         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1772       } else if (ys-Ys > 0) {
1773         for (j=0; j<x; j++) { idx[nn++] = -1;}
1774       }
1775       if (n2 >= 0) { /* right below */
1776         x_t = lx[n2 % m];
1777         y_t = ly[(n2/m)];
1778         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
1779         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1780       } else if (Xe-xe> 0 && ys-Ys > 0) {
1781         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
1782       }
1783     }
1784 
1785     for (i=0; i<y; i++) {
1786       if (n3 >= 0) { /* directly left */
1787         x_t = lx[n3 % m];
1788         /* y_t = y; */
1789         s_t = bases[n3] + (i+1)*x_t - s_x;
1790         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1791       } else if (xs-Xs > 0) {
1792         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
1793       }
1794 
1795       for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
1796 
1797       if (n5 >= 0) { /* directly right */
1798         x_t = lx[n5 % m];
1799         /* y_t = y; */
1800         s_t = bases[n5] + (i)*x_t;
1801         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1802       } else if (Xe-xe > 0) {
1803         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
1804       }
1805     }
1806 
1807     for (i=1; i<=s_y; i++) {
1808       if (n6 >= 0) { /* left above */
1809         x_t = lx[n6 % m];
1810         /* y_t = ly[(n6/m)]; */
1811         s_t = bases[n6] + (i)*x_t - s_x;
1812         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1813       } else if (xs-Xs > 0 && Ye-ye > 0) {
1814         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
1815       }
1816       if (n7 >= 0) { /* directly above */
1817         x_t = x;
1818         /* y_t = ly[(n7/m)]; */
1819         s_t = bases[n7] + (i-1)*x_t;
1820         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1821       } else if (Ye-ye > 0) {
1822         for (j=0; j<x; j++) { idx[nn++] = -1;}
1823       }
1824       if (n8 >= 0) { /* right above */
1825         x_t = lx[n8 % m];
1826         /* y_t = ly[(n8/m)]; */
1827         s_t = bases[n8] + (i-1)*x_t;
1828         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1829       } else if (Xe-xe > 0 && Ye-ye > 0) {
1830         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
1831       }
1832     }
1833   }
1834   /*
1835      Set the local to global ordering in the global vector, this allows use
1836      of VecSetValuesLocal().
1837   */
1838   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,&ltogis);CHKERRQ(ierr);
1839   ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr);
1840   ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1841   ierr = ISGetIndices(ltogis, &idx_full);
1842   ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1843   ierr = ISRestoreIndices(ltogis, &idx_full);
1844   ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr);
1845   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
1846   ierr = ISDestroy(&ltogis);CHKERRQ(ierr);
1847   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
1848   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
1849 
1850   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
1851   dd->m  = m;  dd->n  = n;
1852   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
1853   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
1854   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
1855 
1856   ierr = VecDestroy(&local);CHKERRQ(ierr);
1857   ierr = VecDestroy(&global);CHKERRQ(ierr);
1858 
1859   dd->gtol      = gtol;
1860   dd->ltog      = ltog;
1861   dd->idx       = idx_cpy;
1862   dd->Nl        = nn*dof;
1863   dd->base      = base;
1864   da->ops->view = DMView_DA_2d;
1865   dd->ltol = PETSC_NULL;
1866   dd->ao   = PETSC_NULL;
1867 
1868   PetscFunctionReturn(0);
1869 }
1870 
1871 #undef __FUNCT__
1872 #define __FUNCT__ "DMDACreate2d"
1873 /*@C
1874    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
1875    regular array data that is distributed across some processors.
1876 
1877    Collective on MPI_Comm
1878 
1879    Input Parameters:
1880 +  comm - MPI communicator
1881 .  bx,by - type of ghost nodes the array have.
1882          Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC.
1883 .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
1884 .  M,N - global dimension in each direction of the array (use -M and or -N to indicate that it may be set to a different value
1885             from the command line with -da_grid_x <M> -da_grid_y <N>)
1886 .  m,n - corresponding number of processors in each dimension
1887          (or PETSC_DECIDE to have calculated)
1888 .  dof - number of degrees of freedom per node
1889 .  s - stencil width
1890 -  lx, ly - arrays containing the number of nodes in each cell along
1891            the x and y coordinates, or PETSC_NULL. If non-null, these
1892            must be of length as m and n, and the corresponding
1893            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
1894            must be M, and the sum of the ly[] entries must be N.
1895 
1896    Output Parameter:
1897 .  da - the resulting distributed array object
1898 
1899    Options Database Key:
1900 +  -da_view - Calls DMView() at the conclusion of DMDACreate2d()
1901 .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
1902 .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
1903 .  -da_processors_x <nx> - number of processors in x direction
1904 .  -da_processors_y <ny> - number of processors in y direction
1905 .  -da_refine_x <rx> - refinement ratio in x direction
1906 .  -da_refine_y <ry> - refinement ratio in y direction
1907 -  -da_refine <n> - refine the DMDA n times before creating, if M or N < 0
1908 
1909 
1910    Level: beginner
1911 
1912    Notes:
1913    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1914    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
1915    the standard 9-pt stencil.
1916 
1917    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1918    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1919    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
1920 
1921 .keywords: distributed array, create, two-dimensional
1922 
1923 .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
1924           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(),
1925           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
1926 
1927 @*/
1928 
1929 PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDAStencilType stencil_type,
1930                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
1931 {
1932   PetscErrorCode ierr;
1933 
1934   PetscFunctionBegin;
1935   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1936   ierr = DMDASetDim(*da, 2);CHKERRQ(ierr);
1937   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
1938   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
1939   ierr = DMDASetBoundaryType(*da, bx, by, DMDA_BOUNDARY_NONE);CHKERRQ(ierr);
1940   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1941   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1942   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1943   ierr = DMDASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr);
1944   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
1945   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
1946   ierr = DMSetUp(*da);CHKERRQ(ierr);
1947   ierr = DMView_DA_Private(*da);CHKERRQ(ierr);
1948   PetscFunctionReturn(0);
1949 }
1950