xref: /petsc/src/dm/impls/da/dasub.c (revision db05f41b88540d61f4bd88b8a8012f1ff01b70f8)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6b45d2f2cSJed Brown #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
747c6ae99SBarry Smith 
847c6ae99SBarry Smith #undef __FUNCT__
9*db05f41bSBarry Smith #define __FUNCT__ "DMDAGetRay"
10*db05f41bSBarry Smith /*@C
11*db05f41bSBarry Smith    DMDAGetRay - Returns a vector on process zero that contains a row or column of the values in a DMDA vector
12*db05f41bSBarry Smith 
13*db05f41bSBarry Smith    Collective on DMDA
14*db05f41bSBarry Smith 
15*db05f41bSBarry Smith    Input Parameters:
16*db05f41bSBarry Smith +  da - the distributed array
17*db05f41bSBarry Smith .  vec - the vector
18*db05f41bSBarry Smith .  dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z
19*db05f41bSBarry Smith -  gp - global grid point number in this direction
20*db05f41bSBarry Smith 
21*db05f41bSBarry Smith    Output Parameters:
22*db05f41bSBarry Smith +  newvec - the new vector that can hold the values (size zero on all processes except process 0)
23*db05f41bSBarry Smith -  scatter - the VecScatter that will map from the original vector to the slice
24*db05f41bSBarry Smith 
25*db05f41bSBarry Smith    Level: advanced
26*db05f41bSBarry Smith 
27*db05f41bSBarry Smith    Notes:
28*db05f41bSBarry Smith    All processors that share the DMDA must call this with the same gp value
29*db05f41bSBarry Smith 
30*db05f41bSBarry Smith .keywords: distributed array, get, processor subset
31*db05f41bSBarry Smith @*/
32*db05f41bSBarry Smith PetscErrorCode  DMDAGetRay(DM da,DMDADirection dir,PetscInt gp,Vec *newvec,VecScatter *scatter)
33*db05f41bSBarry Smith {
34*db05f41bSBarry Smith   PetscMPIInt    rank;
35*db05f41bSBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
36*db05f41bSBarry Smith   PetscErrorCode ierr;
37*db05f41bSBarry Smith   IS             is;
38*db05f41bSBarry Smith   AO             ao;
39*db05f41bSBarry Smith   Vec            vec;
40*db05f41bSBarry Smith 
41*db05f41bSBarry Smith   PetscFunctionBegin;
42*db05f41bSBarry Smith   if (dd->dim == 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Cannot get slice from 1d DMDA");
43*db05f41bSBarry Smith   if (dd->dim == 3) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Cannot get slice from 3d DMDA");
44*db05f41bSBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
45*db05f41bSBarry Smith   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
46*db05f41bSBarry Smith   if (!rank) {
47*db05f41bSBarry Smith     if (dir == DMDA_Y) {
48*db05f41bSBarry Smith       PetscInt *indices,i;
49*db05f41bSBarry Smith       ierr = PetscMalloc(dd->M*sizeof(PetscInt),&indices);CHKERRQ(ierr);
50*db05f41bSBarry Smith       indices[0] = gp*dd->M;
51*db05f41bSBarry Smith       for (i=1; i<dd->M; i++) {indices[i] = indices[i-1] + 1;}
52*db05f41bSBarry Smith       ierr = AOApplicationToPetsc(ao,dd->M,indices);CHKERRQ(ierr);
53*db05f41bSBarry Smith       ierr = VecCreate(PETSC_COMM_SELF,newvec);CHKERRQ(ierr);
54*db05f41bSBarry Smith       ierr = VecSetBlockSize(*newvec,dd->w);CHKERRQ(ierr);
55*db05f41bSBarry Smith       ierr = VecSetSizes(*newvec,dd->M*dd->w,PETSC_DETERMINE);CHKERRQ(ierr);
56*db05f41bSBarry Smith       ierr = VecSetType(*newvec,VECSEQ);CHKERRQ(ierr);
57*db05f41bSBarry Smith       ierr = ISCreateBlock(PETSC_COMM_SELF,dd->w,dd->M,indices,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
58*db05f41bSBarry Smith     } else if (dir == DMDA_X) {
59*db05f41bSBarry Smith       PetscInt *indices,i;
60*db05f41bSBarry Smith       ierr = PetscMalloc(dd->N*sizeof(PetscInt),&indices);CHKERRQ(ierr);
61*db05f41bSBarry Smith       indices[0] = gp;
62*db05f41bSBarry Smith       for (i=1; i<dd->N; i++) {indices[i] = indices[i-1] + dd->M;}
63*db05f41bSBarry Smith       ierr = AOApplicationToPetsc(ao,dd->M,indices);CHKERRQ(ierr);
64*db05f41bSBarry Smith       ierr = VecCreate(PETSC_COMM_SELF,newvec);CHKERRQ(ierr);
65*db05f41bSBarry Smith       ierr = VecSetBlockSize(*newvec,dd->w);CHKERRQ(ierr);
66*db05f41bSBarry Smith       ierr = VecSetSizes(*newvec,dd->N*dd->w,PETSC_DETERMINE);CHKERRQ(ierr);
67*db05f41bSBarry Smith       ierr = VecSetType(*newvec,VECSEQ);CHKERRQ(ierr);
68*db05f41bSBarry Smith       ierr = ISCreateBlock(PETSC_COMM_SELF,dd->w,dd->N,indices,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
69*db05f41bSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown DMDADirection");
70*db05f41bSBarry Smith   } else {
71*db05f41bSBarry Smith     ierr = VecCreateSeq(PETSC_COMM_SELF,0,newvec);CHKERRQ(ierr);
72*db05f41bSBarry Smith     ierr = ISCreateBlock(PETSC_COMM_SELF,dd->w,0,0,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
73*db05f41bSBarry Smith   }
74*db05f41bSBarry Smith   ierr = DMGetGlobalVector(da,&vec);CHKERRQ(ierr);
75*db05f41bSBarry Smith   ierr = VecScatterCreate(vec,is,*newvec,PETSC_NULL,scatter);CHKERRQ(ierr);
76*db05f41bSBarry Smith   ierr = DMRestoreGlobalVector(da,&vec);CHKERRQ(ierr);
77*db05f41bSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
78*db05f41bSBarry Smith   PetscFunctionReturn(0);
79*db05f41bSBarry Smith }
80*db05f41bSBarry Smith 
81*db05f41bSBarry Smith #undef __FUNCT__
82aa219208SBarry Smith #define __FUNCT__ "DMDAGetProcessorSubset"
8347c6ae99SBarry Smith /*@C
84aa219208SBarry Smith    DMDAGetProcessorSubset - Returns a communicator consisting only of the
85aa219208SBarry Smith    processors in a DMDA that own a particular global x, y, or z grid point
8647c6ae99SBarry Smith    (corresponding to a logical plane in a 3D grid or a line in a 2D grid).
8747c6ae99SBarry Smith 
88aa219208SBarry Smith    Collective on DMDA
8947c6ae99SBarry Smith 
9047c6ae99SBarry Smith    Input Parameters:
9147c6ae99SBarry Smith +  da - the distributed array
92aa219208SBarry Smith .  dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z
9347c6ae99SBarry Smith -  gp - global grid point number in this direction
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith    Output Parameters:
9647c6ae99SBarry Smith .  comm - new communicator
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith    Level: advanced
9947c6ae99SBarry Smith 
10047c6ae99SBarry Smith    Notes:
101aa219208SBarry Smith    All processors that share the DMDA must call this with the same gp value
10247c6ae99SBarry Smith 
10347c6ae99SBarry Smith    This routine is particularly useful to compute boundary conditions
10447c6ae99SBarry Smith    or other application-specific calculations that require manipulating
10547c6ae99SBarry Smith    sets of data throughout a logical plane of grid points.
10647c6ae99SBarry Smith 
10747c6ae99SBarry Smith .keywords: distributed array, get, processor subset
10847c6ae99SBarry Smith @*/
1097087cfbeSBarry Smith PetscErrorCode  DMDAGetProcessorSubset(DM da,DMDADirection dir,PetscInt gp,MPI_Comm *comm)
11047c6ae99SBarry Smith {
11147c6ae99SBarry Smith   MPI_Group      group,subgroup;
11247c6ae99SBarry Smith   PetscErrorCode ierr;
11347c6ae99SBarry Smith   PetscInt       i,ict,flag,*owners,xs,xm,ys,ym,zs,zm;
11447c6ae99SBarry Smith   PetscMPIInt    size,*ranks = PETSC_NULL;
11547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith   PetscFunctionBegin;
11847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
11947c6ae99SBarry Smith   flag = 0;
120aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
12147c6ae99SBarry Smith   ierr = MPI_Comm_size(((PetscObject)da)->comm,&size);CHKERRQ(ierr);
122aa219208SBarry Smith   if (dir == DMDA_Z) {
123aa219208SBarry Smith     if (dd->dim < 3) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Z invalid for DMDA dim < 3");
12447c6ae99SBarry Smith     if (gp < 0 || gp > dd->P) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
12547c6ae99SBarry Smith     if (gp >= zs && gp < zs+zm) flag = 1;
126aa219208SBarry Smith   } else if (dir == DMDA_Y) {
127aa219208SBarry Smith     if (dd->dim == 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Y invalid for DMDA dim = 1");
12847c6ae99SBarry Smith     if (gp < 0 || gp > dd->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
12947c6ae99SBarry Smith     if (gp >= ys && gp < ys+ym) flag = 1;
130aa219208SBarry Smith   } else if (dir == DMDA_X) {
13147c6ae99SBarry Smith     if (gp < 0 || gp > dd->M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
13247c6ae99SBarry Smith     if (gp >= xs && gp < xs+xm) flag = 1;
13347c6ae99SBarry Smith   } else SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");
13447c6ae99SBarry Smith 
13547c6ae99SBarry Smith   ierr = PetscMalloc2(size,PetscInt,&owners,size,PetscMPIInt,&ranks);CHKERRQ(ierr);
13647c6ae99SBarry Smith   ierr = MPI_Allgather(&flag,1,MPIU_INT,owners,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
13747c6ae99SBarry Smith   ict  = 0;
138aa219208SBarry Smith   ierr = PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",dd->dim,(int)dir);CHKERRQ(ierr);
13947c6ae99SBarry Smith   for (i=0; i<size; i++) {
14047c6ae99SBarry Smith     if (owners[i]) {
14147c6ae99SBarry Smith       ranks[ict] = i; ict++;
14247c6ae99SBarry Smith       ierr = PetscInfo1(da,"%D ",i);CHKERRQ(ierr);
14347c6ae99SBarry Smith     }
14447c6ae99SBarry Smith   }
14547c6ae99SBarry Smith   ierr = PetscInfo(da,"\n");CHKERRQ(ierr);
14647c6ae99SBarry Smith   ierr = MPI_Comm_group(((PetscObject)da)->comm,&group);CHKERRQ(ierr);
14747c6ae99SBarry Smith   ierr = MPI_Group_incl(group,ict,ranks,&subgroup);CHKERRQ(ierr);
14847c6ae99SBarry Smith   ierr = MPI_Comm_create(((PetscObject)da)->comm,subgroup,comm);CHKERRQ(ierr);
14947c6ae99SBarry Smith   ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr);
15047c6ae99SBarry Smith   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
15147c6ae99SBarry Smith   ierr = PetscFree2(owners,ranks);CHKERRQ(ierr);
15247c6ae99SBarry Smith   PetscFunctionReturn(0);
15347c6ae99SBarry Smith }
15447c6ae99SBarry Smith 
15547c6ae99SBarry Smith #undef __FUNCT__
156aa219208SBarry Smith #define __FUNCT__ "DMDAGetProcessorSubsets"
15747c6ae99SBarry Smith /*@C
158aa219208SBarry Smith    DMDAGetProcessorSubsets - Returns communicators consisting only of the
159aa219208SBarry Smith    processors in a DMDA adjacent in a particular dimension,
16047c6ae99SBarry Smith    corresponding to a logical plane in a 3D grid or a line in a 2D grid.
16147c6ae99SBarry Smith 
162aa219208SBarry Smith    Collective on DMDA
16347c6ae99SBarry Smith 
16447c6ae99SBarry Smith    Input Parameters:
16547c6ae99SBarry Smith +  da - the distributed array
166aa219208SBarry Smith -  dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z
16747c6ae99SBarry Smith 
16847c6ae99SBarry Smith    Output Parameters:
16947c6ae99SBarry Smith .  subcomm - new communicator
17047c6ae99SBarry Smith 
17147c6ae99SBarry Smith    Level: advanced
17247c6ae99SBarry Smith 
17347c6ae99SBarry Smith    Notes:
17447c6ae99SBarry Smith    This routine is useful for distributing one-dimensional data in a tensor product grid.
17547c6ae99SBarry Smith 
17647c6ae99SBarry Smith .keywords: distributed array, get, processor subset
17747c6ae99SBarry Smith @*/
1787087cfbeSBarry Smith PetscErrorCode  DMDAGetProcessorSubsets(DM da, DMDADirection dir, MPI_Comm *subcomm)
17947c6ae99SBarry Smith {
18047c6ae99SBarry Smith   MPI_Comm       comm;
18147c6ae99SBarry Smith   MPI_Group      group, subgroup;
18247c6ae99SBarry Smith   PetscInt       subgroupSize = 0;
18347c6ae99SBarry Smith   PetscInt      *firstPoints;
18447c6ae99SBarry Smith   PetscMPIInt    size, *subgroupRanks = PETSC_NULL;
18547c6ae99SBarry Smith   PetscInt       xs, xm, ys, ym, zs, zm, firstPoint, p;
18647c6ae99SBarry Smith   PetscErrorCode ierr;
18747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
18847c6ae99SBarry Smith 
18947c6ae99SBarry Smith   PetscFunctionBegin;
19047c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
19147c6ae99SBarry Smith   comm = ((PetscObject) da)->comm;
192aa219208SBarry Smith   ierr = DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm);CHKERRQ(ierr);
19347c6ae99SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
194aa219208SBarry Smith   if (dir == DMDA_Z) {
195aa219208SBarry Smith     if (dd->dim < 3) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Z invalid for DMDA dim < 3");
19647c6ae99SBarry Smith     firstPoint = zs;
197aa219208SBarry Smith   } else if (dir == DMDA_Y) {
198aa219208SBarry Smith     if (dd->dim == 1) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Y invalid for DMDA dim = 1");
19947c6ae99SBarry Smith     firstPoint = ys;
200aa219208SBarry Smith   } else if (dir == DMDA_X) {
20147c6ae99SBarry Smith     firstPoint = xs;
20247c6ae99SBarry Smith   } else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");
20347c6ae99SBarry Smith 
20447c6ae99SBarry Smith   ierr = PetscMalloc2(size, PetscInt, &firstPoints, size, PetscMPIInt, &subgroupRanks);CHKERRQ(ierr);
20547c6ae99SBarry Smith   ierr = MPI_Allgather(&firstPoint, 1, MPIU_INT, firstPoints, 1, MPIU_INT, comm);CHKERRQ(ierr);
206aa219208SBarry Smith   ierr = PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",dd->dim,(int)dir);CHKERRQ(ierr);
20747c6ae99SBarry Smith   for (p = 0; p < size; ++p) {
20847c6ae99SBarry Smith     if (firstPoints[p] == firstPoint) {
20947c6ae99SBarry Smith       subgroupRanks[subgroupSize++] = p;
21047c6ae99SBarry Smith       ierr = PetscInfo1(da, "%D ", p);CHKERRQ(ierr);
21147c6ae99SBarry Smith     }
21247c6ae99SBarry Smith   }
21347c6ae99SBarry Smith   ierr = PetscInfo(da, "\n");CHKERRQ(ierr);
21447c6ae99SBarry Smith   ierr = MPI_Comm_group(comm, &group);CHKERRQ(ierr);
21547c6ae99SBarry Smith   ierr = MPI_Group_incl(group, subgroupSize, subgroupRanks, &subgroup);CHKERRQ(ierr);
21647c6ae99SBarry Smith   ierr = MPI_Comm_create(comm, subgroup, subcomm);CHKERRQ(ierr);
21747c6ae99SBarry Smith   ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr);
21847c6ae99SBarry Smith   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
21947c6ae99SBarry Smith   ierr = PetscFree2(firstPoints, subgroupRanks);CHKERRQ(ierr);
22047c6ae99SBarry Smith   PetscFunctionReturn(0);
22147c6ae99SBarry Smith }
222