xref: /petsc/src/dm/impls/da/dasub.c (revision 9a42bb27a39f0cdf3306a1e22d33cd9809484eaa)
1 #define PETSCDM_DLL
2 
3 /*
4   Code for manipulating distributed regular arrays in parallel.
5 */
6 
7 #include "private/daimpl.h"    /*I   "petscda.h"   I*/
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "DAGetProcessorSubset"
11 /*@C
12    DAGetProcessorSubset - Returns a communicator consisting only of the
13    processors in a DA that own a particular global x, y, or z grid point
14    (corresponding to a logical plane in a 3D grid or a line in a 2D grid).
15 
16    Collective on DA
17 
18    Input Parameters:
19 +  da - the distributed array
20 .  dir - Cartesian direction, either DA_X, DA_Y, or DA_Z
21 -  gp - global grid point number in this direction
22 
23    Output Parameters:
24 .  comm - new communicator
25 
26    Level: advanced
27 
28    Notes:
29    All processors that share the DA must call this with the same gp value
30 
31    This routine is particularly useful to compute boundary conditions
32    or other application-specific calculations that require manipulating
33    sets of data throughout a logical plane of grid points.
34 
35 .keywords: distributed array, get, processor subset
36 @*/
37 PetscErrorCode PETSCDM_DLLEXPORT DAGetProcessorSubset(DM da,DADirection dir,PetscInt gp,MPI_Comm *comm)
38 {
39   MPI_Group      group,subgroup;
40   PetscErrorCode ierr;
41   PetscInt       i,ict,flag,*owners,xs,xm,ys,ym,zs,zm;
42   PetscMPIInt    size,*ranks = PETSC_NULL;
43   DM_DA          *dd = (DM_DA*)da->data;
44 
45   PetscFunctionBegin;
46   PetscValidHeaderSpecific(da,DM_CLASSID,1);
47   flag = 0;
48   ierr = DAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
49   ierr = MPI_Comm_size(((PetscObject)da)->comm,&size);CHKERRQ(ierr);
50   if (dir == DA_Z) {
51     if (dd->dim < 3) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"DA_Z invalid for DA dim < 3");
52     if (gp < 0 || gp > dd->P) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
53     if (gp >= zs && gp < zs+zm) flag = 1;
54   } else if (dir == DA_Y) {
55     if (dd->dim == 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"DA_Y invalid for DA dim = 1");
56     if (gp < 0 || gp > dd->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
57     if (gp >= ys && gp < ys+ym) flag = 1;
58   } else if (dir == DA_X) {
59     if (gp < 0 || gp > dd->M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
60     if (gp >= xs && gp < xs+xm) flag = 1;
61   } else SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");
62 
63   ierr = PetscMalloc2(size,PetscInt,&owners,size,PetscMPIInt,&ranks);CHKERRQ(ierr);
64   ierr = MPI_Allgather(&flag,1,MPIU_INT,owners,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
65   ict  = 0;
66   ierr = PetscInfo2(da,"DAGetProcessorSubset: dim=%D, direction=%d, procs: ",dd->dim,(int)dir);CHKERRQ(ierr);
67   for (i=0; i<size; i++) {
68     if (owners[i]) {
69       ranks[ict] = i; ict++;
70       ierr = PetscInfo1(da,"%D ",i);CHKERRQ(ierr);
71     }
72   }
73   ierr = PetscInfo(da,"\n");CHKERRQ(ierr);
74   ierr = MPI_Comm_group(((PetscObject)da)->comm,&group);CHKERRQ(ierr);
75   ierr = MPI_Group_incl(group,ict,ranks,&subgroup);CHKERRQ(ierr);
76   ierr = MPI_Comm_create(((PetscObject)da)->comm,subgroup,comm);CHKERRQ(ierr);
77   ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr);
78   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
79   ierr = PetscFree2(owners,ranks);CHKERRQ(ierr);
80   PetscFunctionReturn(0);
81 }
82 
83 #undef __FUNCT__
84 #define __FUNCT__ "DAGetProcessorSubsets"
85 /*@C
86    DAGetProcessorSubsets - Returns communicators consisting only of the
87    processors in a DA adjacent in a particular dimension,
88    corresponding to a logical plane in a 3D grid or a line in a 2D grid.
89 
90    Collective on DA
91 
92    Input Parameters:
93 +  da - the distributed array
94 -  dir - Cartesian direction, either DA_X, DA_Y, or DA_Z
95 
96    Output Parameters:
97 .  subcomm - new communicator
98 
99    Level: advanced
100 
101    Notes:
102    This routine is useful for distributing one-dimensional data in a tensor product grid.
103 
104 .keywords: distributed array, get, processor subset
105 @*/
106 PetscErrorCode PETSCDM_DLLEXPORT DAGetProcessorSubsets(DM da, DADirection dir, MPI_Comm *subcomm)
107 {
108   MPI_Comm       comm;
109   MPI_Group      group, subgroup;
110   PetscInt       subgroupSize = 0;
111   PetscInt      *firstPoints;
112   PetscMPIInt    size, *subgroupRanks = PETSC_NULL;
113   PetscInt       xs, xm, ys, ym, zs, zm, firstPoint, p;
114   PetscErrorCode ierr;
115   DM_DA          *dd = (DM_DA*)da->data;
116 
117   PetscFunctionBegin;
118   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
119   comm = ((PetscObject) da)->comm;
120   ierr = DAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm);CHKERRQ(ierr);
121   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
122   if (dir == DA_Z) {
123     if (dd->dim < 3) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DA_Z invalid for DA dim < 3");
124     firstPoint = zs;
125   } else if (dir == DA_Y) {
126     if (dd->dim == 1) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DA_Y invalid for DA dim = 1");
127     firstPoint = ys;
128   } else if (dir == DA_X) {
129     firstPoint = xs;
130   } else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");
131 
132   ierr = PetscMalloc2(size, PetscInt, &firstPoints, size, PetscMPIInt, &subgroupRanks);CHKERRQ(ierr);
133   ierr = MPI_Allgather(&firstPoint, 1, MPIU_INT, firstPoints, 1, MPIU_INT, comm);CHKERRQ(ierr);
134   ierr = PetscInfo2(da,"DAGetProcessorSubset: dim=%D, direction=%d, procs: ",dd->dim,(int)dir);CHKERRQ(ierr);
135   for(p = 0; p < size; ++p) {
136     if (firstPoints[p] == firstPoint) {
137       subgroupRanks[subgroupSize++] = p;
138       ierr = PetscInfo1(da, "%D ", p);CHKERRQ(ierr);
139     }
140   }
141   ierr = PetscInfo(da, "\n");CHKERRQ(ierr);
142   ierr = MPI_Comm_group(comm, &group);CHKERRQ(ierr);
143   ierr = MPI_Group_incl(group, subgroupSize, subgroupRanks, &subgroup);CHKERRQ(ierr);
144   ierr = MPI_Comm_create(comm, subgroup, subcomm);CHKERRQ(ierr);
145   ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr);
146   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
147   ierr = PetscFree2(firstPoints, subgroupRanks);CHKERRQ(ierr);
148   PetscFunctionReturn(0);
149 }
150