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