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