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