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