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