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