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