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