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