xref: /petsc/src/dm/impls/da/da.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1 #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
2 
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMDASetDim"
6 /*@
7   DMDASetDim - Sets the dimension
8 
9   Collective on DMDA
10 
11   Input Parameters:
12 + da - the DMDA
13 - dim - the dimension (or PETSC_DECIDE)
14 
15   Level: intermediate
16 
17 .seealso: DaGetDim(), DMDASetSizes()
18 @*/
19 PetscErrorCode  DMDASetDim(DM da, PetscInt dim)
20 {
21   DM_DA *dd = (DM_DA*)da->data;
22 
23   PetscFunctionBegin;
24   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
25   if (dd->dim > 0 && dim != dd->dim) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change DMDA dim from %D after it was set to %D",dd->dim,dim);
26   dd->dim = dim;
27   PetscFunctionReturn(0);
28 }
29 
30 #undef __FUNCT__
31 #define __FUNCT__ "DMDASetSizes"
32 /*@
33   DMDASetSizes - Sets the global sizes
34 
35   Logically Collective on DMDA
36 
37   Input Parameters:
38 + da - the DMDA
39 . M - the global X size (or PETSC_DECIDE)
40 . N - the global Y size (or PETSC_DECIDE)
41 - P - the global Z size (or PETSC_DECIDE)
42 
43   Level: intermediate
44 
45 .seealso: DMDAGetSize(), PetscSplitOwnership()
46 @*/
47 PetscErrorCode  DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P)
48 {
49   DM_DA *dd = (DM_DA*)da->data;
50 
51   PetscFunctionBegin;
52   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
53   PetscValidLogicalCollectiveInt(da,M,2);
54   PetscValidLogicalCollectiveInt(da,N,3);
55   PetscValidLogicalCollectiveInt(da,P,4);
56   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
57 
58   dd->M = M;
59   dd->N = N;
60   dd->P = P;
61   PetscFunctionReturn(0);
62 }
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "DMDASetNumProcs"
66 /*@
67   DMDASetNumProcs - Sets the number of processes in each dimension
68 
69   Logically Collective on DMDA
70 
71   Input Parameters:
72 + da - the DMDA
73 . m - the number of X procs (or PETSC_DECIDE)
74 . n - the number of Y procs (or PETSC_DECIDE)
75 - p - the number of Z procs (or PETSC_DECIDE)
76 
77   Level: intermediate
78 
79 .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership()
80 @*/
81 PetscErrorCode  DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p)
82 {
83   DM_DA          *dd = (DM_DA*)da->data;
84   PetscErrorCode ierr;
85 
86   PetscFunctionBegin;
87   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
88   PetscValidLogicalCollectiveInt(da,m,2);
89   PetscValidLogicalCollectiveInt(da,n,3);
90   PetscValidLogicalCollectiveInt(da,p,4);
91   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
92   dd->m = m;
93   dd->n = n;
94   dd->p = p;
95   if (dd->dim == 2) {
96     PetscMPIInt size;
97     ierr = MPI_Comm_size(((PetscObject)da)->comm,&size);CHKERRQ(ierr);
98     if ((dd->m > 0) && (dd->n < 0)) {
99       dd->n = size/dd->m;
100       if (dd->n*dd->m != size) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"%D processes in X direction not divisible into comm size %d",m,size);
101     }
102     if ((dd->n > 0) && (dd->m < 0)) {
103       dd->m = size/dd->n;
104       if (dd->n*dd->m != size) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"%D processes in Y direction not divisible into comm size %d",n,size);
105     }
106   }
107   PetscFunctionReturn(0);
108 }
109 
110 #undef __FUNCT__
111 #define __FUNCT__ "DMDASetBoundaryType"
112 /*@
113   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
114 
115   Not collective
116 
117   Input Parameter:
118 + da    - The DMDA
119 - bx,by,bz - One of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC
120 
121   Level: intermediate
122 
123 .keywords:  distributed array, periodicity
124 .seealso: DMDACreate(), DMDestroy(), DMDA, DMDABoundaryType
125 @*/
126 PetscErrorCode  DMDASetBoundaryType(DM da,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz)
127 {
128   DM_DA *dd = (DM_DA*)da->data;
129 
130   PetscFunctionBegin;
131   PetscValidHeaderSpecific(da,DM_CLASSID,1);
132   PetscValidLogicalCollectiveEnum(da,bx,2);
133   PetscValidLogicalCollectiveEnum(da,by,3);
134   PetscValidLogicalCollectiveEnum(da,bz,4);
135   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
136   dd->bx = bx;
137   dd->by = by;
138   dd->bz = bz;
139   PetscFunctionReturn(0);
140 }
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "DMDASetDof"
144 /*@
145   DMDASetDof - Sets the number of degrees of freedom per vertex
146 
147   Not collective
148 
149   Input Parameter:
150 + da  - The DMDA
151 - dof - Number of degrees of freedom
152 
153   Level: intermediate
154 
155 .keywords:  distributed array, degrees of freedom
156 .seealso: DMDACreate(), DMDestroy(), DMDA
157 @*/
158 PetscErrorCode  DMDASetDof(DM da, PetscInt dof)
159 {
160   DM_DA *dd = (DM_DA*)da->data;
161 
162   PetscFunctionBegin;
163   PetscValidHeaderSpecific(da,DM_CLASSID,1);
164   PetscValidLogicalCollectiveInt(da,dof,2);
165   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
166   dd->w  = dof;
167   da->bs = dof;
168   PetscFunctionReturn(0);
169 }
170 
171 #undef __FUNCT__
172 #define __FUNCT__ "DMDASetOverlap"
173 /*@
174   DMDASetOverlap - Sets the size of the per-processor overlap.
175 
176   Not collective
177 
178   Input Parameter:
179 + da  - The DMDA
180 - dof - Number of degrees of freedom
181 
182   Level: intermediate
183 
184 .keywords:  distributed array, degrees of freedom
185 .seealso: DMDACreate(), DMDestroy(), DMDA
186 @*/
187 PetscErrorCode  DMDASetOverlap(DM da, PetscInt overlap)
188 {
189   DM_DA *dd = (DM_DA*)da->data;
190 
191   PetscFunctionBegin;
192   PetscValidHeaderSpecific(da,DM_CLASSID,1);
193   PetscValidLogicalCollectiveInt(da,overlap,2);
194   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
195   dd->overlap = overlap;
196   PetscFunctionReturn(0);
197 }
198 
199 #undef __FUNCT__
200 #define __FUNCT__ "DMDASetOffset"
201 /*@
202   DMDASetOffset - Sets the index offset of the DA.
203 
204   Collective on DA
205 
206   Input Parameter:
207 + da  - The DMDA
208 . xo  - The offset in the x direction
209 . yo  - The offset in the y direction
210 - zo  - The offset in the z direction
211 
212   Level: intermediate
213 
214   Notes: This is used primarily to overlap a computation on a local DA with that on a global DA without
215   changing boundary conditions or subdomain features that depend upon the global offsets.
216 
217 .keywords:  distributed array, degrees of freedom
218 .seealso: DMDAGetOffset(), DMDAVecGetArray()
219 @*/
220 PetscErrorCode  DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po)
221 {
222   PetscErrorCode ierr;
223   DM_DA          *dd = (DM_DA*)da->data;
224 
225   PetscFunctionBegin;
226   PetscValidHeaderSpecific(da,DM_CLASSID,1);
227   PetscValidLogicalCollectiveInt(da,xo,2);
228   PetscValidLogicalCollectiveInt(da,yo,3);
229   PetscValidLogicalCollectiveInt(da,zo,4);
230   PetscValidLogicalCollectiveInt(da,Mo,5);
231   PetscValidLogicalCollectiveInt(da,No,6);
232   PetscValidLogicalCollectiveInt(da,Po,7);
233   dd->xo = xo;
234   dd->yo = yo;
235   dd->zo = zo;
236   dd->Mo = Mo;
237   dd->No = No;
238   dd->Po = Po;
239 
240   if (da->coordinateDM) {
241     ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr);
242   }
243   PetscFunctionReturn(0);
244 }
245 
246 #undef __FUNCT__
247 #define __FUNCT__ "DMDAGetOffset"
248 /*@
249   DMDAGetOffset - Gets the index offset of the DA.
250 
251   Not collective
252 
253   Input Parameter:
254 . da  - The DMDA
255 
256   Output Parameters:
257 + xo  - The offset in the x direction
258 . yo  - The offset in the y direction
259 . zo  - The offset in the z direction
260 . Mo  - The global size in the x direction
261 . No  - The global size in the y direction
262 - Po  - The global size in the z direction
263 
264   Level: intermediate
265 
266 .keywords:  distributed array, degrees of freedom
267 .seealso: DMDASetOffset(), DMDAVecGetArray()
268 @*/
269 PetscErrorCode  DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po)
270 {
271   DM_DA *dd = (DM_DA*)da->data;
272 
273   PetscFunctionBegin;
274   PetscValidHeaderSpecific(da,DM_CLASSID,1);
275   if (xo) *xo = dd->xo;
276   if (yo) *yo = dd->yo;
277   if (zo) *zo = dd->zo;
278   if (Mo) *Mo = dd->Mo;
279   if (No) *No = dd->No;
280   if (Po) *Po = dd->Po;
281   PetscFunctionReturn(0);
282 }
283 
284 
285 #undef __FUNCT__
286 #define __FUNCT__ "DMDASetStencilType"
287 /*@
288   DMDASetStencilType - Sets the type of the communication stencil
289 
290   Logically Collective on DMDA
291 
292   Input Parameter:
293 + da    - The DMDA
294 - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
295 
296   Level: intermediate
297 
298 .keywords:  distributed array, stencil
299 .seealso: DMDACreate(), DMDestroy(), DMDA
300 @*/
301 PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
302 {
303   DM_DA *dd = (DM_DA*)da->data;
304 
305   PetscFunctionBegin;
306   PetscValidHeaderSpecific(da,DM_CLASSID,1);
307   PetscValidLogicalCollectiveEnum(da,stype,2);
308   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
309   dd->stencil_type = stype;
310   PetscFunctionReturn(0);
311 }
312 
313 #undef __FUNCT__
314 #define __FUNCT__ "DMDASetStencilWidth"
315 /*@
316   DMDASetStencilWidth - Sets the width of the communication stencil
317 
318   Logically Collective on DMDA
319 
320   Input Parameter:
321 + da    - The DMDA
322 - width - The stencil width
323 
324   Level: intermediate
325 
326 .keywords:  distributed array, stencil
327 .seealso: DMDACreate(), DMDestroy(), DMDA
328 @*/
329 PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
330 {
331   DM_DA *dd = (DM_DA*)da->data;
332 
333   PetscFunctionBegin;
334   PetscValidHeaderSpecific(da,DM_CLASSID,1);
335   PetscValidLogicalCollectiveInt(da,width,2);
336   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
337   dd->s = width;
338   PetscFunctionReturn(0);
339 }
340 
341 #undef __FUNCT__
342 #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
343 static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
344 {
345   PetscInt i,sum;
346 
347   PetscFunctionBegin;
348   if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
349   for (i=sum=0; i<m; i++) sum += lx[i];
350   if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
351   PetscFunctionReturn(0);
352 }
353 
354 #undef __FUNCT__
355 #define __FUNCT__ "DMDASetOwnershipRanges"
356 /*@
357   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
358 
359   Logically Collective on DMDA
360 
361   Input Parameter:
362 + da - The DMDA
363 . lx - array containing number of nodes in the X direction on each process, or PETSC_NULL. If non-null, must be of length da->m
364 . ly - array containing number of nodes in the Y direction on each process, or PETSC_NULL. If non-null, must be of length da->n
365 - lz - array containing number of nodes in the Z direction on each process, or PETSC_NULL. If non-null, must be of length da->p.
366 
367   Level: intermediate
368 
369   Note: these numbers are NOT multiplied by the number of dof per node.
370 
371 .keywords:  distributed array
372 .seealso: DMDACreate(), DMDestroy(), DMDA
373 @*/
374 PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
375 {
376   PetscErrorCode ierr;
377   DM_DA          *dd = (DM_DA*)da->data;
378 
379   PetscFunctionBegin;
380   PetscValidHeaderSpecific(da,DM_CLASSID,1);
381   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
382   if (lx) {
383     if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
384     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
385     if (!dd->lx) {
386       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
387     }
388     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
389   }
390   if (ly) {
391     if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
392     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
393     if (!dd->ly) {
394       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
395     }
396     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
397   }
398   if (lz) {
399     if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
400     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
401     if (!dd->lz) {
402       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
403     }
404     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
405   }
406   PetscFunctionReturn(0);
407 }
408 
409 #undef __FUNCT__
410 #define __FUNCT__ "DMDASetInterpolationType"
411 /*@
412        DMDASetInterpolationType - Sets the type of interpolation that will be
413           returned by DMCreateInterpolation()
414 
415    Logically Collective on DMDA
416 
417    Input Parameter:
418 +  da - initial distributed array
419 .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
420 
421    Level: intermediate
422 
423    Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation()
424 
425 .keywords:  distributed array, interpolation
426 
427 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
428 @*/
429 PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
430 {
431   DM_DA *dd = (DM_DA*)da->data;
432 
433   PetscFunctionBegin;
434   PetscValidHeaderSpecific(da,DM_CLASSID,1);
435   PetscValidLogicalCollectiveEnum(da,ctype,2);
436   dd->interptype = ctype;
437   PetscFunctionReturn(0);
438 }
439 
440 #undef __FUNCT__
441 #define __FUNCT__ "DMDAGetInterpolationType"
442 /*@
443        DMDAGetInterpolationType - Gets the type of interpolation that will be
444           used by DMCreateInterpolation()
445 
446    Not Collective
447 
448    Input Parameter:
449 .  da - distributed array
450 
451    Output Parameter:
452 .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
453 
454    Level: intermediate
455 
456 .keywords:  distributed array, interpolation
457 
458 .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation()
459 @*/
460 PetscErrorCode  DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype)
461 {
462   DM_DA *dd = (DM_DA*)da->data;
463 
464   PetscFunctionBegin;
465   PetscValidHeaderSpecific(da,DM_CLASSID,1);
466   PetscValidPointer(ctype,2);
467   *ctype = dd->interptype;
468   PetscFunctionReturn(0);
469 }
470 
471 #undef __FUNCT__
472 #define __FUNCT__ "DMDAGetNeighbors"
473 /*@C
474       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
475         processes neighbors.
476 
477     Not Collective
478 
479    Input Parameter:
480 .     da - the DMDA object
481 
482    Output Parameters:
483 .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
484               this process itself is in the list
485 
486    Notes: In 2d the array is of length 9, in 3d of length 27
487           Not supported in 1d
488           Do not free the array, it is freed when the DMDA is destroyed.
489 
490    Fortran Notes: In fortran you must pass in an array of the appropriate length.
491 
492    Level: intermediate
493 
494 @*/
495 PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
496 {
497   DM_DA *dd = (DM_DA*)da->data;
498 
499   PetscFunctionBegin;
500   PetscValidHeaderSpecific(da,DM_CLASSID,1);
501   *ranks = dd->neighbors;
502   PetscFunctionReturn(0);
503 }
504 
505 #undef __FUNCT__
506 #define __FUNCT__ "DMDAGetOwnershipRanges"
507 /*@C
508       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
509 
510     Not Collective
511 
512    Input Parameter:
513 .     da - the DMDA object
514 
515    Output Parameter:
516 +     lx - ownership along x direction (optional)
517 .     ly - ownership along y direction (optional)
518 -     lz - ownership along z direction (optional)
519 
520    Level: intermediate
521 
522     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
523 
524     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
525     eighth arguments from DMDAGetInfo()
526 
527      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
528     DMDA they came from still exists (has not been destroyed).
529 
530     These numbers are NOT multiplied by the number of dof per node.
531 
532 .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
533 @*/
534 PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
535 {
536   DM_DA *dd = (DM_DA*)da->data;
537 
538   PetscFunctionBegin;
539   PetscValidHeaderSpecific(da,DM_CLASSID,1);
540   if (lx) *lx = dd->lx;
541   if (ly) *ly = dd->ly;
542   if (lz) *lz = dd->lz;
543   PetscFunctionReturn(0);
544 }
545 
546 #undef __FUNCT__
547 #define __FUNCT__ "DMDASetRefinementFactor"
548 /*@
549      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
550 
551     Logically Collective on DMDA
552 
553   Input Parameters:
554 +    da - the DMDA object
555 .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
556 .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
557 -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
558 
559   Options Database:
560 +  -da_refine_x - refinement ratio in x direction
561 .  -da_refine_y - refinement ratio in y direction
562 -  -da_refine_z - refinement ratio in z direction
563 
564   Level: intermediate
565 
566     Notes: Pass PETSC_IGNORE to leave a value unchanged
567 
568 .seealso: DMRefine(), DMDAGetRefinementFactor()
569 @*/
570 PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
571 {
572   DM_DA *dd = (DM_DA*)da->data;
573 
574   PetscFunctionBegin;
575   PetscValidHeaderSpecific(da,DM_CLASSID,1);
576   PetscValidLogicalCollectiveInt(da,refine_x,2);
577   PetscValidLogicalCollectiveInt(da,refine_y,3);
578   PetscValidLogicalCollectiveInt(da,refine_z,4);
579 
580   if (refine_x > 0) dd->refine_x = refine_x;
581   if (refine_y > 0) dd->refine_y = refine_y;
582   if (refine_z > 0) dd->refine_z = refine_z;
583   PetscFunctionReturn(0);
584 }
585 
586 #undef __FUNCT__
587 #define __FUNCT__ "DMDAGetRefinementFactor"
588 /*@C
589      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
590 
591     Not Collective
592 
593   Input Parameter:
594 .    da - the DMDA object
595 
596   Output Parameters:
597 +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
598 .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
599 -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
600 
601   Level: intermediate
602 
603     Notes: Pass PETSC_NULL for values you do not need
604 
605 .seealso: DMRefine(), DMDASetRefinementFactor()
606 @*/
607 PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
608 {
609   DM_DA *dd = (DM_DA*)da->data;
610 
611   PetscFunctionBegin;
612   PetscValidHeaderSpecific(da,DM_CLASSID,1);
613   if (refine_x) *refine_x = dd->refine_x;
614   if (refine_y) *refine_y = dd->refine_y;
615   if (refine_z) *refine_z = dd->refine_z;
616   PetscFunctionReturn(0);
617 }
618 
619 #undef __FUNCT__
620 #define __FUNCT__ "DMDASetGetMatrix"
621 /*@C
622      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
623 
624     Logically Collective on DMDA
625 
626   Input Parameters:
627 +    da - the DMDA object
628 -    f - the function that allocates the matrix for that specific DMDA
629 
630   Level: developer
631 
632    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
633        the diagonal and off-diagonal blocks of the matrix
634 
635 .seealso: DMCreateMatrix(), DMDASetBlockFills()
636 @*/
637 PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, MatType,Mat*))
638 {
639   PetscFunctionBegin;
640   PetscValidHeaderSpecific(da,DM_CLASSID,1);
641   da->ops->creatematrix = f;
642   PetscFunctionReturn(0);
643 }
644 
645 #undef __FUNCT__
646 #define __FUNCT__ "DMDARefineOwnershipRanges"
647 /*
648   Creates "balanced" ownership ranges after refinement, constrained by the need for the
649   fine grid boundaries to fall within one stencil width of the coarse partition.
650 
651   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
652 */
653 static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
654 {
655   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
656   PetscErrorCode ierr;
657 
658   PetscFunctionBegin;
659   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
660   if (ratio == 1) {
661     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
662     PetscFunctionReturn(0);
663   }
664   for (i=0; i<m; i++) totalc += lc[i];
665   remaining = (!periodic) + ratio * (totalc - (!periodic));
666   for (i=0; i<m; i++) {
667     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
668     if (i == m-1) lf[i] = want;
669     else {
670       const PetscInt nextc = startc + lc[i];
671       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
672        * coarse stencil width of the first coarse node in the next subdomain. */
673       while ((startf+want)/ratio < nextc - stencil_width) want++;
674       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
675        * coarse stencil width of the last coarse node in the current subdomain. */
676       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
677       /* Make sure all constraints are satisfied */
678       if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width)
679           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
680     }
681     lf[i]      = want;
682     startc    += lc[i];
683     startf    += lf[i];
684     remaining -= lf[i];
685   }
686   PetscFunctionReturn(0);
687 }
688 
689 #undef __FUNCT__
690 #define __FUNCT__ "DMDACoarsenOwnershipRanges"
691 /*
692   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
693   fine grid boundaries to fall within one stencil width of the coarse partition.
694 
695   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
696 */
697 static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
698 {
699   PetscInt       i,totalf,remaining,startc,startf;
700   PetscErrorCode ierr;
701 
702   PetscFunctionBegin;
703   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
704   if (ratio == 1) {
705     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
706     PetscFunctionReturn(0);
707   }
708   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
709   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
710   for (i=0,startc=0,startf=0; i<m; i++) {
711     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
712     if (i == m-1) lc[i] = want;
713     else {
714       const PetscInt nextf = startf+lf[i];
715       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
716        * node is within one stencil width. */
717       while (nextf/ratio < startc+want-stencil_width) want--;
718       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
719        * fine node is within one stencil width. */
720       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
721       if (want < 0 || want > remaining
722           || (nextf/ratio < startc+want-stencil_width) || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
723     }
724     lc[i]      = want;
725     startc    += lc[i];
726     startf    += lf[i];
727     remaining -= lc[i];
728   }
729   PetscFunctionReturn(0);
730 }
731 
732 #undef __FUNCT__
733 #define __FUNCT__ "DMRefine_DA"
734 PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
735 {
736   PetscErrorCode ierr;
737   PetscInt       M,N,P,i;
738   DM             da2;
739   DM_DA          *dd = (DM_DA*)da->data,*dd2;
740 
741   PetscFunctionBegin;
742   PetscValidHeaderSpecific(da,DM_CLASSID,1);
743   PetscValidPointer(daref,3);
744 
745   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
746     M = dd->refine_x*dd->M;
747   } else {
748     M = 1 + dd->refine_x*(dd->M - 1);
749   }
750   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
751     if (dd->dim > 1) {
752       N = dd->refine_y*dd->N;
753     } else {
754       N = 1;
755     }
756   } else {
757     N = 1 + dd->refine_y*(dd->N - 1);
758   }
759   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
760     if (dd->dim > 2) {
761       P = dd->refine_z*dd->P;
762     } else {
763       P = 1;
764     }
765   } else {
766     P = 1 + dd->refine_z*(dd->P - 1);
767   }
768   ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr);
769   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
770   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
771   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
772   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
773   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
774   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
775   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
776   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
777   if (dd->dim == 3) {
778     PetscInt *lx,*ly,*lz;
779     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
780     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
781     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
782     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
783     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
784     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
785   } else if (dd->dim == 2) {
786     PetscInt *lx,*ly;
787     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
788     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
789     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
790     ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr);
791     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
792   } else if (dd->dim == 1) {
793     PetscInt *lx;
794     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
795     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
796     ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
797     ierr = PetscFree(lx);CHKERRQ(ierr);
798   }
799   dd2 = (DM_DA*)da2->data;
800 
801   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
802   da2->ops->creatematrix = da->ops->creatematrix;
803   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
804   da2->ops->getcoloring = da->ops->getcoloring;
805   dd2->interptype       = dd->interptype;
806 
807   /* copy fill information if given */
808   if (dd->dfill) {
809     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
810     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
811   }
812   if (dd->ofill) {
813     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
814     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
815   }
816   /* copy the refine information */
817   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
818   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
819   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
820 
821   /* copy vector type information */
822   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
823   ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr);
824 
825   dd2->lf = dd->lf;
826   dd2->lj = dd->lj;
827 
828   da2->leveldown = da->leveldown;
829   da2->levelup   = da->levelup + 1;
830 
831   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
832   ierr = DMSetUp(da2);CHKERRQ(ierr);
833   ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr);
834 
835   /* interpolate coordinates if they are set on the coarse grid */
836   if (da->coordinates) {
837     DM  cdaf,cdac;
838     Vec coordsc,coordsf;
839     Mat II;
840 
841     ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr);
842     ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr);
843     ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr);
844     /* force creation of the coordinate vector */
845     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
846     ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
847     ierr = DMCreateInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr);
848     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
849     ierr = MatDestroy(&II);CHKERRQ(ierr);
850   }
851 
852   for (i=0; i<da->bs; i++) {
853     const char *fieldname;
854     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
855     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
856   }
857 
858   *daref = da2;
859   PetscFunctionReturn(0);
860 }
861 
862 
863 #undef __FUNCT__
864 #define __FUNCT__ "DMCoarsen_DA"
865 PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
866 {
867   PetscErrorCode ierr;
868   PetscInt       M,N,P,i;
869   DM             da2;
870   DM_DA          *dd = (DM_DA*)da->data,*dd2;
871 
872   PetscFunctionBegin;
873   PetscValidHeaderSpecific(da,DM_CLASSID,1);
874   PetscValidPointer(daref,3);
875 
876   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
877     M = dd->M / dd->coarsen_x;
878   } else {
879     M = 1 + (dd->M - 1) / dd->coarsen_x;
880   }
881   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
882     if (dd->dim > 1) {
883       N = dd->N / dd->coarsen_y;
884     } else {
885       N = 1;
886     }
887   } else {
888     N = 1 + (dd->N - 1) / dd->coarsen_y;
889   }
890   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
891     if (dd->dim > 2) {
892       P = dd->P / dd->coarsen_z;
893     } else {
894       P = 1;
895     }
896   } else {
897     P = 1 + (dd->P - 1) / dd->coarsen_z;
898   }
899   ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr);
900   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
901   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
902   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
903   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
904   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
905   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
906   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
907   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
908   if (dd->dim == 3) {
909     PetscInt *lx,*ly,*lz;
910     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
911     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
912     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
913     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
914     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
915     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
916   } else if (dd->dim == 2) {
917     PetscInt *lx,*ly;
918     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
919     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
920     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
921     ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr);
922     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
923   } else if (dd->dim == 1) {
924     PetscInt *lx;
925     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
926     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
927     ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
928     ierr = PetscFree(lx);CHKERRQ(ierr);
929   }
930   dd2 = (DM_DA*)da2->data;
931 
932   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
933   /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */
934   da2->ops->creatematrix = da->ops->creatematrix;
935   da2->ops->getcoloring  = da->ops->getcoloring;
936   dd2->interptype        = dd->interptype;
937 
938   /* copy fill information if given */
939   if (dd->dfill) {
940     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
941     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
942   }
943   if (dd->ofill) {
944     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
945     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
946   }
947   /* copy the refine information */
948   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
949   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
950   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
951 
952   /* copy vector type information */
953   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
954   ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr);
955 
956   dd2->lf = dd->lf;
957   dd2->lj = dd->lj;
958 
959   da2->leveldown = da->leveldown + 1;
960   da2->levelup   = da->levelup;
961 
962   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
963   ierr = DMSetUp(da2);CHKERRQ(ierr);
964   ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr);
965 
966   /* inject coordinates if they are set on the fine grid */
967   if (da->coordinates) {
968     DM         cdaf,cdac;
969     Vec        coordsc,coordsf;
970     VecScatter inject;
971 
972     ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr);
973     ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr);
974     ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr);
975     /* force creation of the coordinate vector */
976     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
977     ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
978 
979     ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
980     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
981     ierr = VecScatterEnd(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
982     ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
983   }
984 
985   for (i=0; i<da->bs; i++) {
986     const char *fieldname;
987     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
988     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
989   }
990 
991   *daref = da2;
992   PetscFunctionReturn(0);
993 }
994 
995 #undef __FUNCT__
996 #define __FUNCT__ "DMRefineHierarchy_DA"
997 PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
998 {
999   PetscErrorCode ierr;
1000   PetscInt       i,n,*refx,*refy,*refz;
1001 
1002   PetscFunctionBegin;
1003   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1004   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1005   if (nlevels == 0) PetscFunctionReturn(0);
1006   PetscValidPointer(daf,3);
1007 
1008   /* Get refinement factors, defaults taken from the coarse DMDA */
1009   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
1010   for (i=0; i<nlevels; i++) {
1011     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
1012   }
1013   n    = nlevels;
1014   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
1015   n    = nlevels;
1016   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
1017   n    = nlevels;
1018   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
1019 
1020   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
1021   ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr);
1022   for (i=1; i<nlevels; i++) {
1023     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
1024     ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr);
1025   }
1026   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
1027   PetscFunctionReturn(0);
1028 }
1029 
1030 #undef __FUNCT__
1031 #define __FUNCT__ "DMCoarsenHierarchy_DA"
1032 PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
1033 {
1034   PetscErrorCode ierr;
1035   PetscInt       i;
1036 
1037   PetscFunctionBegin;
1038   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1039   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1040   if (nlevels == 0) PetscFunctionReturn(0);
1041   PetscValidPointer(dac,3);
1042   ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr);
1043   for (i=1; i<nlevels; i++) {
1044     ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr);
1045   }
1046   PetscFunctionReturn(0);
1047 }
1048