xref: /petsc/src/dm/impls/da/da.c (revision 5fd668637986a8d8518383a9159eebc368e1d5b4)
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 
244   PetscFunctionReturn(0);
245 }
246 
247 #undef __FUNCT__
248 #define __FUNCT__ "DMDAGetOffset"
249 /*@
250   DMDAGetOffset - Gets the index offset of the DA.
251 
252   Not collective
253 
254   Input Parameter:
255 . da  - The DMDA
256 
257   Output Parameters:
258 + xo  - The offset in the x direction
259 . yo  - The offset in the y direction
260 . zo  - The offset in the z direction
261 . Mo  - The global size in the x direction
262 . No  - The global size in the y direction
263 - Po  - The global size in the z direction
264 
265   Level: intermediate
266 
267 .keywords:  distributed array, degrees of freedom
268 .seealso: DMDASetOffset(), DMDAVecGetArray()
269 @*/
270 PetscErrorCode  DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po)
271 {
272   DM_DA *dd = (DM_DA*)da->data;
273 
274   PetscFunctionBegin;
275   PetscValidHeaderSpecific(da,DM_CLASSID,1);
276   if (xo) *xo = dd->xo;
277   if (yo) *yo = dd->yo;
278   if (zo) *zo = dd->zo;
279   if (Mo) *Mo = dd->Mo;
280   if (No) *No = dd->No;
281   if (Po) *Po = dd->Po;
282   PetscFunctionReturn(0);
283 }
284 
285 
286 #undef __FUNCT__
287 #define __FUNCT__ "DMDASetStencilType"
288 /*@
289   DMDASetStencilType - Sets the type of the communication stencil
290 
291   Logically Collective on DMDA
292 
293   Input Parameter:
294 + da    - The DMDA
295 - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
296 
297   Level: intermediate
298 
299 .keywords:  distributed array, stencil
300 .seealso: DMDACreate(), DMDestroy(), DMDA
301 @*/
302 PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
303 {
304   DM_DA *dd = (DM_DA*)da->data;
305 
306   PetscFunctionBegin;
307   PetscValidHeaderSpecific(da,DM_CLASSID,1);
308   PetscValidLogicalCollectiveEnum(da,stype,2);
309   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
310   dd->stencil_type = stype;
311   PetscFunctionReturn(0);
312 }
313 
314 #undef __FUNCT__
315 #define __FUNCT__ "DMDASetStencilWidth"
316 /*@
317   DMDASetStencilWidth - Sets the width of the communication stencil
318 
319   Logically Collective on DMDA
320 
321   Input Parameter:
322 + da    - The DMDA
323 - width - The stencil width
324 
325   Level: intermediate
326 
327 .keywords:  distributed array, stencil
328 .seealso: DMDACreate(), DMDestroy(), DMDA
329 @*/
330 PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
331 {
332   DM_DA *dd = (DM_DA*)da->data;
333 
334   PetscFunctionBegin;
335   PetscValidHeaderSpecific(da,DM_CLASSID,1);
336   PetscValidLogicalCollectiveInt(da,width,2);
337   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
338   dd->s = width;
339   PetscFunctionReturn(0);
340 }
341 
342 #undef __FUNCT__
343 #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
344 static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
345 {
346   PetscInt i,sum;
347 
348   PetscFunctionBegin;
349   if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
350   for (i=sum=0; i<m; i++) sum += lx[i];
351   if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
352   PetscFunctionReturn(0);
353 }
354 
355 #undef __FUNCT__
356 #define __FUNCT__ "DMDASetOwnershipRanges"
357 /*@
358   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
359 
360   Logically Collective on DMDA
361 
362   Input Parameter:
363 + da - The DMDA
364 . 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
365 . 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
366 - 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.
367 
368   Level: intermediate
369 
370   Note: these numbers are NOT multiplied by the number of dof per node.
371 
372 .keywords:  distributed array
373 .seealso: DMDACreate(), DMDestroy(), DMDA
374 @*/
375 PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
376 {
377   PetscErrorCode ierr;
378   DM_DA          *dd = (DM_DA*)da->data;
379 
380   PetscFunctionBegin;
381   PetscValidHeaderSpecific(da,DM_CLASSID,1);
382   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
383   if (lx) {
384     if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
385     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
386     if (!dd->lx) {
387       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
388     }
389     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
390   }
391   if (ly) {
392     if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
393     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
394     if (!dd->ly) {
395       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
396     }
397     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
398   }
399   if (lz) {
400     if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
401     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
402     if (!dd->lz) {
403       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
404     }
405     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
406   }
407   PetscFunctionReturn(0);
408 }
409 
410 #undef __FUNCT__
411 #define __FUNCT__ "DMDASetInterpolationType"
412 /*@
413        DMDASetInterpolationType - Sets the type of interpolation that will be
414           returned by DMCreateInterpolation()
415 
416    Logically Collective on DMDA
417 
418    Input Parameter:
419 +  da - initial distributed array
420 .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
421 
422    Level: intermediate
423 
424    Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation()
425 
426 .keywords:  distributed array, interpolation
427 
428 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
429 @*/
430 PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
431 {
432   DM_DA *dd = (DM_DA*)da->data;
433 
434   PetscFunctionBegin;
435   PetscValidHeaderSpecific(da,DM_CLASSID,1);
436   PetscValidLogicalCollectiveEnum(da,ctype,2);
437   dd->interptype = ctype;
438   PetscFunctionReturn(0);
439 }
440 
441 #undef __FUNCT__
442 #define __FUNCT__ "DMDAGetInterpolationType"
443 /*@
444        DMDAGetInterpolationType - Gets the type of interpolation that will be
445           used by DMCreateInterpolation()
446 
447    Not Collective
448 
449    Input Parameter:
450 .  da - distributed array
451 
452    Output Parameter:
453 .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
454 
455    Level: intermediate
456 
457 .keywords:  distributed array, interpolation
458 
459 .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation()
460 @*/
461 PetscErrorCode  DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype)
462 {
463   DM_DA *dd = (DM_DA*)da->data;
464 
465   PetscFunctionBegin;
466   PetscValidHeaderSpecific(da,DM_CLASSID,1);
467   PetscValidPointer(ctype,2);
468   *ctype = dd->interptype;
469   PetscFunctionReturn(0);
470 }
471 
472 #undef __FUNCT__
473 #define __FUNCT__ "DMDAGetNeighbors"
474 /*@C
475       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
476         processes neighbors.
477 
478     Not Collective
479 
480    Input Parameter:
481 .     da - the DMDA object
482 
483    Output Parameters:
484 .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
485               this process itself is in the list
486 
487    Notes: In 2d the array is of length 9, in 3d of length 27
488           Not supported in 1d
489           Do not free the array, it is freed when the DMDA is destroyed.
490 
491    Fortran Notes: In fortran you must pass in an array of the appropriate length.
492 
493    Level: intermediate
494 
495 @*/
496 PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
497 {
498   DM_DA *dd = (DM_DA*)da->data;
499 
500   PetscFunctionBegin;
501   PetscValidHeaderSpecific(da,DM_CLASSID,1);
502   *ranks = dd->neighbors;
503   PetscFunctionReturn(0);
504 }
505 
506 #undef __FUNCT__
507 #define __FUNCT__ "DMDAGetOwnershipRanges"
508 /*@C
509       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
510 
511     Not Collective
512 
513    Input Parameter:
514 .     da - the DMDA object
515 
516    Output Parameter:
517 +     lx - ownership along x direction (optional)
518 .     ly - ownership along y direction (optional)
519 -     lz - ownership along z direction (optional)
520 
521    Level: intermediate
522 
523     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
524 
525     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
526     eighth arguments from DMDAGetInfo()
527 
528      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
529     DMDA they came from still exists (has not been destroyed).
530 
531     These numbers are NOT multiplied by the number of dof per node.
532 
533 .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
534 @*/
535 PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
536 {
537   DM_DA *dd = (DM_DA*)da->data;
538 
539   PetscFunctionBegin;
540   PetscValidHeaderSpecific(da,DM_CLASSID,1);
541   if (lx) *lx = dd->lx;
542   if (ly) *ly = dd->ly;
543   if (lz) *lz = dd->lz;
544   PetscFunctionReturn(0);
545 }
546 
547 #undef __FUNCT__
548 #define __FUNCT__ "DMDASetRefinementFactor"
549 /*@
550      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
551 
552     Logically Collective on DMDA
553 
554   Input Parameters:
555 +    da - the DMDA object
556 .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
557 .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
558 -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
559 
560   Options Database:
561 +  -da_refine_x - refinement ratio in x direction
562 .  -da_refine_y - refinement ratio in y direction
563 -  -da_refine_z - refinement ratio in z direction
564 
565   Level: intermediate
566 
567     Notes: Pass PETSC_IGNORE to leave a value unchanged
568 
569 .seealso: DMRefine(), DMDAGetRefinementFactor()
570 @*/
571 PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
572 {
573   DM_DA *dd = (DM_DA*)da->data;
574 
575   PetscFunctionBegin;
576   PetscValidHeaderSpecific(da,DM_CLASSID,1);
577   PetscValidLogicalCollectiveInt(da,refine_x,2);
578   PetscValidLogicalCollectiveInt(da,refine_y,3);
579   PetscValidLogicalCollectiveInt(da,refine_z,4);
580 
581   if (refine_x > 0) dd->refine_x = refine_x;
582   if (refine_y > 0) dd->refine_y = refine_y;
583   if (refine_z > 0) dd->refine_z = refine_z;
584   PetscFunctionReturn(0);
585 }
586 
587 #undef __FUNCT__
588 #define __FUNCT__ "DMDAGetRefinementFactor"
589 /*@C
590      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
591 
592     Not Collective
593 
594   Input Parameter:
595 .    da - the DMDA object
596 
597   Output Parameters:
598 +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
599 .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
600 -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
601 
602   Level: intermediate
603 
604     Notes: Pass PETSC_NULL for values you do not need
605 
606 .seealso: DMRefine(), DMDASetRefinementFactor()
607 @*/
608 PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
609 {
610   DM_DA *dd = (DM_DA*)da->data;
611 
612   PetscFunctionBegin;
613   PetscValidHeaderSpecific(da,DM_CLASSID,1);
614   if (refine_x) *refine_x = dd->refine_x;
615   if (refine_y) *refine_y = dd->refine_y;
616   if (refine_z) *refine_z = dd->refine_z;
617   PetscFunctionReturn(0);
618 }
619 
620 #undef __FUNCT__
621 #define __FUNCT__ "DMDASetGetMatrix"
622 /*@C
623      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
624 
625     Logically Collective on DMDA
626 
627   Input Parameters:
628 +    da - the DMDA object
629 -    f - the function that allocates the matrix for that specific DMDA
630 
631   Level: developer
632 
633    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
634        the diagonal and off-diagonal blocks of the matrix
635 
636 .seealso: DMCreateMatrix(), DMDASetBlockFills()
637 @*/
638 PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, MatType,Mat*))
639 {
640   PetscFunctionBegin;
641   PetscValidHeaderSpecific(da,DM_CLASSID,1);
642   da->ops->creatematrix = f;
643   PetscFunctionReturn(0);
644 }
645 
646 #undef __FUNCT__
647 #define __FUNCT__ "DMDARefineOwnershipRanges"
648 /*
649   Creates "balanced" ownership ranges after refinement, constrained by the need for the
650   fine grid boundaries to fall within one stencil width of the coarse partition.
651 
652   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
653 */
654 static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
655 {
656   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
657   PetscErrorCode ierr;
658 
659   PetscFunctionBegin;
660   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
661   if (ratio == 1) {
662     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
663     PetscFunctionReturn(0);
664   }
665   for (i=0; i<m; i++) totalc += lc[i];
666   remaining = (!periodic) + ratio * (totalc - (!periodic));
667   for (i=0; i<m; i++) {
668     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
669     if (i == m-1) lf[i] = want;
670     else {
671       const PetscInt nextc = startc + lc[i];
672       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
673        * coarse stencil width of the first coarse node in the next subdomain. */
674       while ((startf+want)/ratio < nextc - stencil_width) want++;
675       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
676        * coarse stencil width of the last coarse node in the current subdomain. */
677       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
678       /* Make sure all constraints are satisfied */
679       if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width)
680           || ((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");
681     }
682     lf[i] = want;
683     startc += lc[i];
684     startf += lf[i];
685     remaining -= lf[i];
686   }
687   PetscFunctionReturn(0);
688 }
689 
690 #undef __FUNCT__
691 #define __FUNCT__ "DMDACoarsenOwnershipRanges"
692 /*
693   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
694   fine grid boundaries to fall within one stencil width of the coarse partition.
695 
696   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
697 */
698 static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
699 {
700   PetscInt       i,totalf,remaining,startc,startf;
701   PetscErrorCode ierr;
702 
703   PetscFunctionBegin;
704   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
705   if (ratio == 1) {
706     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
707     PetscFunctionReturn(0);
708   }
709   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
710   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
711   for (i=0,startc=0,startf=0; i<m; i++) {
712     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
713     if (i == m-1) lc[i] = want;
714     else {
715       const PetscInt nextf = startf+lf[i];
716       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
717        * node is within one stencil width. */
718       while (nextf/ratio < startc+want-stencil_width) want--;
719       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
720        * fine node is within one stencil width. */
721       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
722       if (want < 0 || want > remaining
723           || (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");
724     }
725     lc[i] = want;
726     startc += lc[i];
727     startf += lf[i];
728     remaining -= lc[i];
729   }
730   PetscFunctionReturn(0);
731 }
732 
733 #undef __FUNCT__
734 #define __FUNCT__ "DMRefine_DA"
735 PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
736 {
737   PetscErrorCode ierr;
738   PetscInt       M,N,P,i;
739   DM             da2;
740   DM_DA          *dd = (DM_DA*)da->data,*dd2;
741 
742   PetscFunctionBegin;
743   PetscValidHeaderSpecific(da,DM_CLASSID,1);
744   PetscValidPointer(daref,3);
745 
746   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
747     M = dd->refine_x*dd->M;
748   } else {
749     M = 1 + dd->refine_x*(dd->M - 1);
750   }
751   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
752     if (dd->dim > 1) {
753       N = dd->refine_y*dd->N;
754     } else {
755       N = 1;
756     }
757   } else {
758     N = 1 + dd->refine_y*(dd->N - 1);
759   }
760   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
761     if (dd->dim > 2) {
762       P = dd->refine_z*dd->P;
763     } else {
764       P = 1;
765     }
766   } else {
767     P = 1 + dd->refine_z*(dd->P - 1);
768   }
769   ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr);
770   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
771   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
772   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
773   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
774   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
775   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
776   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
777   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
778   if (dd->dim == 3) {
779     PetscInt *lx,*ly,*lz;
780     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
781     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);
782     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);
783     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);
784     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
785     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
786   } else if (dd->dim == 2) {
787     PetscInt *lx,*ly;
788     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
789     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);
790     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);
791     ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr);
792     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
793   } else if (dd->dim == 1) {
794     PetscInt *lx;
795     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
796     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);
797     ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
798     ierr = PetscFree(lx);CHKERRQ(ierr);
799   }
800   dd2 = (DM_DA*)da2->data;
801 
802   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
803   da2->ops->creatematrix        = da->ops->creatematrix;
804   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
805   da2->ops->getcoloring      = da->ops->getcoloring;
806   dd2->interptype            = dd->interptype;
807 
808   /* copy fill information if given */
809   if (dd->dfill) {
810     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
811     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
812   }
813   if (dd->ofill) {
814     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
815     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
816   }
817   /* copy the refine information */
818   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
819   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
820   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
821 
822   /* copy vector type information */
823   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
824   ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr);
825 
826   dd2->lf = dd->lf;
827   dd2->lj = dd->lj;
828 
829   da2->leveldown = da->leveldown;
830   da2->levelup   = da->levelup + 1;
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   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
962   ierr = DMSetUp(da2);CHKERRQ(ierr);
963   ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr);
964 
965   /* inject coordinates if they are set on the fine grid */
966   if (da->coordinates) {
967     DM         cdaf,cdac;
968     Vec        coordsc,coordsf;
969     VecScatter inject;
970 
971     ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr);
972     ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr);
973     ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr);
974     /* force creation of the coordinate vector */
975     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
976     ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
977 
978     ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
979     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
980     ierr = VecScatterEnd(inject  ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
981     ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
982   }
983 
984   for (i=0; i<da->bs; i++) {
985     const char *fieldname;
986     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
987     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
988   }
989 
990   *daref = da2;
991   PetscFunctionReturn(0);
992 }
993 
994 #undef __FUNCT__
995 #define __FUNCT__ "DMRefineHierarchy_DA"
996 PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
997 {
998   PetscErrorCode ierr;
999   PetscInt       i,n,*refx,*refy,*refz;
1000 
1001   PetscFunctionBegin;
1002   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1003   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1004   if (nlevels == 0) PetscFunctionReturn(0);
1005   PetscValidPointer(daf,3);
1006 
1007   /* Get refinement factors, defaults taken from the coarse DMDA */
1008   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
1009   for (i=0; i<nlevels; i++) {
1010     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
1011   }
1012   n = nlevels;
1013   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
1014   n = nlevels;
1015   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
1016   n = nlevels;
1017   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
1018 
1019   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
1020   ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr);
1021   for (i=1; i<nlevels; i++) {
1022     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
1023     ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr);
1024   }
1025   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
1026   PetscFunctionReturn(0);
1027 }
1028 
1029 #undef __FUNCT__
1030 #define __FUNCT__ "DMCoarsenHierarchy_DA"
1031 PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
1032 {
1033   PetscErrorCode ierr;
1034   PetscInt       i;
1035 
1036   PetscFunctionBegin;
1037   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1038   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1039   if (nlevels == 0) PetscFunctionReturn(0);
1040   PetscValidPointer(dac,3);
1041   ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr);
1042   for (i=1; i<nlevels; i++) {
1043     ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr);
1044   }
1045   PetscFunctionReturn(0);
1046 }
1047