xref: /petsc/src/dm/impls/da/da.c (revision 3e08d2bebf9958b2a7617e5d4703c5f9c1d28ec4)
1 #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
2 
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMDASetDim"
6 /*@
7   DMDASetDim - Sets the dimension
8 
9   Collective on DMDA
10 
11   Input Parameters:
12 + da - the DMDA
13 - dim - the dimension (or PETSC_DECIDE)
14 
15   Level: intermediate
16 
17 .seealso: DaGetDim(), DMDASetSizes()
18 @*/
19 PetscErrorCode  DMDASetDim(DM da, PetscInt dim)
20 {
21   DM_DA *dd = (DM_DA*)da->data;
22 
23   PetscFunctionBegin;
24   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
25   if (dd->dim > 0 && dim != dd->dim) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change DMDA dim from %D after it was set to %D",dd->dim,dim);
26   dd->dim = dim;
27   PetscFunctionReturn(0);
28 }
29 
30 #undef __FUNCT__
31 #define __FUNCT__ "DMDASetSizes"
32 /*@
33   DMDASetSizes - Sets the global sizes
34 
35   Logically Collective on DMDA
36 
37   Input Parameters:
38 + da - the DMDA
39 . M - the global X size (or PETSC_DECIDE)
40 . N - the global Y size (or PETSC_DECIDE)
41 - P - the global Z size (or PETSC_DECIDE)
42 
43   Level: intermediate
44 
45 .seealso: DMDAGetSize(), PetscSplitOwnership()
46 @*/
47 PetscErrorCode  DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P)
48 {
49   DM_DA *dd = (DM_DA*)da->data;
50 
51   PetscFunctionBegin;
52   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
53   PetscValidLogicalCollectiveInt(da,M,2);
54   PetscValidLogicalCollectiveInt(da,N,3);
55   PetscValidLogicalCollectiveInt(da,P,4);
56   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
57 
58   dd->M = M;
59   dd->N = N;
60   dd->P = P;
61   PetscFunctionReturn(0);
62 }
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "DMDASetNumProcs"
66 /*@
67   DMDASetNumProcs - Sets the number of processes in each dimension
68 
69   Logically Collective on DMDA
70 
71   Input Parameters:
72 + da - the DMDA
73 . m - the number of X procs (or PETSC_DECIDE)
74 . n - the number of Y procs (or PETSC_DECIDE)
75 - p - the number of Z procs (or PETSC_DECIDE)
76 
77   Level: intermediate
78 
79 .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership()
80 @*/
81 PetscErrorCode  DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p)
82 {
83   DM_DA *dd = (DM_DA*)da->data;
84 
85   PetscFunctionBegin;
86   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
87   PetscValidLogicalCollectiveInt(da,m,2);
88   PetscValidLogicalCollectiveInt(da,n,3);
89   PetscValidLogicalCollectiveInt(da,p,4);
90   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
91   dd->m = m;
92   dd->n = n;
93   dd->p = p;
94   PetscFunctionReturn(0);
95 }
96 
97 #undef __FUNCT__
98 #define __FUNCT__ "DMDASetBoundaryType"
99 /*@
100   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
101 
102   Not collective
103 
104   Input Parameter:
105 + da    - The DMDA
106 - bx,by,bz - One of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC
107 
108   Level: intermediate
109 
110 .keywords:  distributed array, periodicity
111 .seealso: DMDACreate(), DMDestroy(), DMDA, DMDABoundaryType
112 @*/
113 PetscErrorCode  DMDASetBoundaryType(DM da,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz)
114 {
115   DM_DA *dd = (DM_DA*)da->data;
116 
117   PetscFunctionBegin;
118   PetscValidHeaderSpecific(da,DM_CLASSID,1);
119   PetscValidLogicalCollectiveEnum(da,bx,2);
120   PetscValidLogicalCollectiveEnum(da,by,3);
121   PetscValidLogicalCollectiveEnum(da,bz,4);
122   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
123   dd->bx = bx;
124   dd->by = by;
125   dd->bz = bz;
126   PetscFunctionReturn(0);
127 }
128 
129 #undef __FUNCT__
130 #define __FUNCT__ "DMDASetDof"
131 /*@
132   DMDASetDof - Sets the number of degrees of freedom per vertex
133 
134   Not collective
135 
136   Input Parameter:
137 + da  - The DMDA
138 - dof - Number of degrees of freedom
139 
140   Level: intermediate
141 
142 .keywords:  distributed array, degrees of freedom
143 .seealso: DMDACreate(), DMDestroy(), DMDA
144 @*/
145 PetscErrorCode  DMDASetDof(DM da, PetscInt dof)
146 {
147   DM_DA *dd = (DM_DA*)da->data;
148 
149   PetscFunctionBegin;
150   PetscValidHeaderSpecific(da,DM_CLASSID,1);
151   PetscValidLogicalCollectiveInt(da,dof,2);
152   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
153   dd->w = dof;
154   da->bs = dof;
155   PetscFunctionReturn(0);
156 }
157 
158 #undef __FUNCT__
159 #define __FUNCT__ "DMDASetOverlap"
160 /*@
161   DMDASetOverlap - Sets the size of the per-processor overlap.
162 
163   Not collective
164 
165   Input Parameter:
166 + da  - The DMDA
167 - dof - Number of degrees of freedom
168 
169   Level: intermediate
170 
171 .keywords:  distributed array, degrees of freedom
172 .seealso: DMDACreate(), DMDestroy(), DMDA
173 @*/
174 PetscErrorCode  DMDASetOverlap(DM da, PetscInt overlap)
175 {
176   DM_DA *dd = (DM_DA*)da->data;
177 
178   PetscFunctionBegin;
179   PetscValidHeaderSpecific(da,DM_CLASSID,1);
180   PetscValidLogicalCollectiveInt(da,overlap,2);
181   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
182   dd->overlap = overlap;
183   PetscFunctionReturn(0);
184 }
185 
186 
187 #undef __FUNCT__
188 #define __FUNCT__ "DMDASetStencilType"
189 /*@
190   DMDASetStencilType - Sets the type of the communication stencil
191 
192   Logically Collective on DMDA
193 
194   Input Parameter:
195 + da    - The DMDA
196 - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
197 
198   Level: intermediate
199 
200 .keywords:  distributed array, stencil
201 .seealso: DMDACreate(), DMDestroy(), DMDA
202 @*/
203 PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
204 {
205   DM_DA *dd = (DM_DA*)da->data;
206 
207   PetscFunctionBegin;
208   PetscValidHeaderSpecific(da,DM_CLASSID,1);
209   PetscValidLogicalCollectiveEnum(da,stype,2);
210   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
211   dd->stencil_type = stype;
212   PetscFunctionReturn(0);
213 }
214 
215 #undef __FUNCT__
216 #define __FUNCT__ "DMDASetStencilWidth"
217 /*@
218   DMDASetStencilWidth - Sets the width of the communication stencil
219 
220   Logically Collective on DMDA
221 
222   Input Parameter:
223 + da    - The DMDA
224 - width - The stencil width
225 
226   Level: intermediate
227 
228 .keywords:  distributed array, stencil
229 .seealso: DMDACreate(), DMDestroy(), DMDA
230 @*/
231 PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
232 {
233   DM_DA *dd = (DM_DA*)da->data;
234 
235   PetscFunctionBegin;
236   PetscValidHeaderSpecific(da,DM_CLASSID,1);
237   PetscValidLogicalCollectiveInt(da,width,2);
238   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
239   dd->s = width;
240   PetscFunctionReturn(0);
241 }
242 
243 #undef __FUNCT__
244 #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
245 static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
246 {
247   PetscInt i,sum;
248 
249   PetscFunctionBegin;
250   if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
251   for (i=sum=0; i<m; i++) sum += lx[i];
252   if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
253   PetscFunctionReturn(0);
254 }
255 
256 #undef __FUNCT__
257 #define __FUNCT__ "DMDASetOwnershipRanges"
258 /*@
259   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
260 
261   Logically Collective on DMDA
262 
263   Input Parameter:
264 + da - The DMDA
265 . 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
266 . 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
267 - 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.
268 
269   Level: intermediate
270 
271 .keywords:  distributed array
272 .seealso: DMDACreate(), DMDestroy(), DMDA
273 @*/
274 PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
275 {
276   PetscErrorCode ierr;
277   DM_DA          *dd = (DM_DA*)da->data;
278 
279   PetscFunctionBegin;
280   PetscValidHeaderSpecific(da,DM_CLASSID,1);
281   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
282   if (lx) {
283     if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
284     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
285     if (!dd->lx) {
286       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
287     }
288     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
289   }
290   if (ly) {
291     if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
292     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
293     if (!dd->ly) {
294       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
295     }
296     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
297   }
298   if (lz) {
299     if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
300     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
301     if (!dd->lz) {
302       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
303     }
304     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));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 DMCreateInterpolation()
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 DMCreateInterpolation()
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 DMCreateInterpolation()
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(), DMCreateInterpolation()
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: DMCreateMatrix(), 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->creatematrix = 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 || ((startf+want)/ratio < nextc - stencil_width)
576           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
577     }
578     lf[i] = want;
579     startc += lc[i];
580     startf += lf[i];
581     remaining -= lf[i];
582   }
583   PetscFunctionReturn(0);
584 }
585 
586 #undef __FUNCT__
587 #define __FUNCT__ "DMDACoarsenOwnershipRanges"
588 /*
589   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
590   fine grid boundaries to fall within one stencil width of the coarse partition.
591 
592   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
593 */
594 static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
595 {
596   PetscInt       i,totalf,remaining,startc,startf;
597   PetscErrorCode ierr;
598 
599   PetscFunctionBegin;
600   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
601   if (ratio == 1) {
602     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
603     PetscFunctionReturn(0);
604   }
605   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
606   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
607   for (i=0,startc=0,startf=0; i<m; i++) {
608     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
609     if (i == m-1) lc[i] = want;
610     else {
611       const PetscInt nextf = startf+lf[i];
612       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
613        * node is within one stencil width. */
614       while (nextf/ratio < startc+want-stencil_width) want--;
615       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
616        * fine node is within one stencil width. */
617       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
618       if (want < 0 || want > remaining
619           || (nextf/ratio < startc+want-stencil_width) || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
620     }
621     lc[i] = want;
622     startc += lc[i];
623     startf += lf[i];
624     remaining -= lc[i];
625   }
626   PetscFunctionReturn(0);
627 }
628 
629 #undef __FUNCT__
630 #define __FUNCT__ "DMRefine_DA"
631 PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
632 {
633   PetscErrorCode ierr;
634   PetscInt       M,N,P,i;
635   DM             da2;
636   DM_DA          *dd = (DM_DA*)da->data,*dd2;
637 
638   PetscFunctionBegin;
639   PetscValidHeaderSpecific(da,DM_CLASSID,1);
640   PetscValidPointer(daref,3);
641 
642   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
643     M = dd->refine_x*dd->M;
644   } else {
645     M = 1 + dd->refine_x*(dd->M - 1);
646   }
647   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
648     if (dd->dim > 1) {
649       N = dd->refine_y*dd->N;
650     } else {
651       N = 1;
652     }
653   } else {
654     N = 1 + dd->refine_y*(dd->N - 1);
655   }
656   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
657     if (dd->dim > 2) {
658       P = dd->refine_z*dd->P;
659     } else {
660       P = 1;
661     }
662   } else {
663     P = 1 + dd->refine_z*(dd->P - 1);
664   }
665   ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr);
666   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
667   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
668   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
669   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
670   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
671   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
672   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
673   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
674   if (dd->dim == 3) {
675     PetscInt *lx,*ly,*lz;
676     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
677     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);
678     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);
679     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);
680     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
681     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
682   } else if (dd->dim == 2) {
683     PetscInt *lx,*ly;
684     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
685     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);
686     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);
687     ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr);
688     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
689   } else if (dd->dim == 1) {
690     PetscInt *lx;
691     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
692     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);
693     ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
694     ierr = PetscFree(lx);CHKERRQ(ierr);
695   }
696   dd2 = (DM_DA*)da2->data;
697 
698   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
699   da2->ops->creatematrix        = da->ops->creatematrix;
700   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
701   da2->ops->getcoloring      = da->ops->getcoloring;
702   dd2->interptype            = dd->interptype;
703 
704   /* copy fill information if given */
705   if (dd->dfill) {
706     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
707     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
708   }
709   if (dd->ofill) {
710     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
711     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
712   }
713   /* copy the refine information */
714   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
715   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
716   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
717 
718   /* copy vector type information */
719   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
720   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
721 
722   dd2->lf = dd->lf;
723   dd2->lj = dd->lj;
724 
725   da2->leveldown = da->leveldown;
726   da2->levelup   = da->levelup + 1;
727   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
728   ierr = DMSetUp(da2);CHKERRQ(ierr);
729   ierr = DMView_DA_Private(da2);CHKERRQ(ierr);
730 
731   /* interpolate coordinates if they are set on the coarse grid */
732   if (dd->coordinates) {
733     DM  cdaf,cdac;
734     Vec coordsc,coordsf;
735     Mat II;
736 
737     ierr = DMDAGetCoordinateDA(da,&cdac);CHKERRQ(ierr);
738     ierr = DMDAGetCoordinates(da,&coordsc);CHKERRQ(ierr);
739     ierr = DMDAGetCoordinateDA(da2,&cdaf);CHKERRQ(ierr);
740     /* force creation of the coordinate vector */
741     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
742     ierr = DMDAGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
743     ierr = DMCreateInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr);
744     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
745     ierr = MatDestroy(&II);CHKERRQ(ierr);
746   }
747 
748   for (i=0; i<da->bs; i++) {
749     const char *fieldname;
750     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
751     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
752   }
753 
754   *daref = da2;
755   PetscFunctionReturn(0);
756 }
757 
758 
759 #undef __FUNCT__
760 #define __FUNCT__ "DMCoarsen_DA"
761 PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
762 {
763   PetscErrorCode ierr;
764   PetscInt       M,N,P,i;
765   DM             da2;
766   DM_DA          *dd = (DM_DA*)da->data,*dd2;
767 
768   PetscFunctionBegin;
769   PetscValidHeaderSpecific(da,DM_CLASSID,1);
770   PetscValidPointer(daref,3);
771 
772   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
773     M = dd->M / dd->coarsen_x;
774   } else {
775     M = 1 + (dd->M - 1) / dd->coarsen_x;
776   }
777   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
778     if (dd->dim > 1) {
779       N = dd->N / dd->coarsen_y;
780     } else {
781       N = 1;
782     }
783   } else {
784     N = 1 + (dd->N - 1) / dd->coarsen_y;
785   }
786   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
787     if (dd->dim > 2) {
788       P = dd->P / dd->coarsen_z;
789     } else {
790       P = 1;
791     }
792   } else {
793     P = 1 + (dd->P - 1) / dd->coarsen_z;
794   }
795   ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr);
796   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
797   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
798   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
799   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
800   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
801   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
802   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
803   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
804   if (dd->dim == 3) {
805     PetscInt *lx,*ly,*lz;
806     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
807     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
808     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
809     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
810     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
811     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
812   } else if (dd->dim == 2) {
813     PetscInt *lx,*ly;
814     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
815     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
816     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
817     ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr);
818     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
819   } else if (dd->dim == 1) {
820     PetscInt *lx;
821     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
822     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
823     ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
824     ierr = PetscFree(lx);CHKERRQ(ierr);
825   }
826   dd2 = (DM_DA*)da2->data;
827 
828   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
829   /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */
830   da2->ops->creatematrix        = da->ops->creatematrix;
831   da2->ops->getcoloring      = da->ops->getcoloring;
832   dd2->interptype            = dd->interptype;
833 
834   /* copy fill information if given */
835   if (dd->dfill) {
836     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
837     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
838   }
839   if (dd->ofill) {
840     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
841     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
842   }
843   /* copy the refine information */
844   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
845   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
846   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
847 
848   /* copy vector type information */
849   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
850   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
851 
852   dd2->lf = dd->lf;
853   dd2->lj = dd->lj;
854 
855   da2->leveldown = da->leveldown + 1;
856   da2->levelup   = da->levelup;
857   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
858   ierr = DMSetUp(da2);CHKERRQ(ierr);
859   ierr = DMView_DA_Private(da2);CHKERRQ(ierr);
860 
861   /* inject coordinates if they are set on the fine grid */
862   if (dd->coordinates) {
863     DM         cdaf,cdac;
864     Vec        coordsc,coordsf;
865     VecScatter inject;
866 
867     ierr = DMDAGetCoordinateDA(da,&cdaf);CHKERRQ(ierr);
868     ierr = DMDAGetCoordinates(da,&coordsf);CHKERRQ(ierr);
869     ierr = DMDAGetCoordinateDA(da2,&cdac);CHKERRQ(ierr);
870     /* force creation of the coordinate vector */
871     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
872     ierr = DMDAGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
873 
874     ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
875     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
876     ierr = VecScatterEnd(inject  ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
877     ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
878   }
879 
880   for (i=0; i<da->bs; i++) {
881     const char *fieldname;
882     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
883     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
884   }
885 
886   *daref = da2;
887   PetscFunctionReturn(0);
888 }
889 
890 #undef __FUNCT__
891 #define __FUNCT__ "DMRefineHierarchy_DA"
892 PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
893 {
894   PetscErrorCode ierr;
895   PetscInt       i,n,*refx,*refy,*refz;
896 
897   PetscFunctionBegin;
898   PetscValidHeaderSpecific(da,DM_CLASSID,1);
899   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
900   if (nlevels == 0) PetscFunctionReturn(0);
901   PetscValidPointer(daf,3);
902 
903   /* Get refinement factors, defaults taken from the coarse DMDA */
904   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
905   for (i=0; i<nlevels; i++) {
906     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
907   }
908   n = nlevels;
909   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
910   n = nlevels;
911   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
912   n = nlevels;
913   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
914 
915   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
916   ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr);
917   for (i=1; i<nlevels; i++) {
918     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
919     ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr);
920   }
921   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
922   PetscFunctionReturn(0);
923 }
924 
925 #undef __FUNCT__
926 #define __FUNCT__ "DMCoarsenHierarchy_DA"
927 PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
928 {
929   PetscErrorCode ierr;
930   PetscInt       i;
931 
932   PetscFunctionBegin;
933   PetscValidHeaderSpecific(da,DM_CLASSID,1);
934   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
935   if (nlevels == 0) PetscFunctionReturn(0);
936   PetscValidPointer(dac,3);
937   ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr);
938   for (i=1; i<nlevels; i++) {
939     ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr);
940   }
941   PetscFunctionReturn(0);
942 }
943