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