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