Actual source code: dasub.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  2: /*
  3:   Code for manipulating distributed regular arrays in parallel.
  4: */

  6: #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/

 10: /*@C
 11:    DMDAGetLogicalCoordinate - Returns a the i,j,k logical coordinate for the closest mesh point to a x,y,z point in the coordinates of the DMDA

 13:    Collective on DMDA

 15:    Input Parameters:
 16: +  da - the distributed array
 17: -  x,y,z - the physical coordinates

 19:    Output Parameters:
 20: +   II, JJ, KK - the logical coordinate (-1 on processes that do not contain that point)
 21: -   X, Y, Z, - (optional) the coordinates of the located grid point

 23:    Level: advanced

 25:    Notes:
 26:    All processors that share the DMDA must call this with the same coordinate value

 28: .keywords: distributed array, get, processor subset
 29: @*/
 30: PetscErrorCode  DMDAGetLogicalCoordinate(DM da,PetscScalar x,PetscScalar y,PetscScalar z,PetscInt *II,PetscInt *JJ,PetscInt *KK,PetscScalar *X,PetscScalar *Y,PetscScalar *Z)
 31: {
 33:   Vec            coors;
 34:   DM             dacoors;
 35:   DMDACoor2d     **c;
 36:   PetscInt       i,j,xs,xm,ys,ym;
 37:   PetscReal      d,D = PETSC_MAX_REAL,Dv;
 38:   PetscMPIInt    rank,root;

 41:   if (da->dim == 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 1d DMDA");
 42:   if (da->dim == 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 3d DMDA");

 44:   *II = -1;
 45:   *JJ = -1;

 47:   DMGetCoordinateDM(da,&dacoors);
 48:   DMDAGetCorners(dacoors,&xs,&ys,NULL,&xm,&ym,NULL);
 49:   DMGetCoordinates(da,&coors);
 50:   DMDAVecGetArrayRead(dacoors,coors,&c);
 51:   for (j=ys; j<ys+ym; j++) {
 52:     for (i=xs; i<xs+xm; i++) {
 53:       d = PetscSqrtReal(PetscRealPart( (c[j][i].x - x)*(c[j][i].x - x) + (c[j][i].y - y)*(c[j][i].y - y) ));
 54:       if (d < D) {
 55:         D   = d;
 56:         *II = i;
 57:         *JJ = j;
 58:       }
 59:     }
 60:   }
 61:   MPIU_Allreduce(&D,&Dv,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da));
 62:   if (D != Dv) {
 63:     *II  = -1;
 64:     *JJ  = -1;
 65:     rank = 0;
 66:   } else {
 67:     *X = c[*JJ][*II].x;
 68:     *Y = c[*JJ][*II].y;
 69:     MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);
 70:     rank++;
 71:   }
 72:   MPIU_Allreduce(&rank,&root,1,MPI_INT,MPI_SUM,PetscObjectComm((PetscObject)da));
 73:   root--;
 74:   MPI_Bcast(X,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da));
 75:   MPI_Bcast(Y,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da));
 76:   DMDAVecRestoreArrayRead(dacoors,coors,&c);
 77:   return(0);
 78: }

 82: /*@C
 83:    DMDAGetRay - Returns a vector on process zero that contains a row or column of the values in a DMDA vector

 85:    Collective on DMDA

 87:    Input Parameters:
 88: +  da - the distributed array
 89: .  vec - the vector
 90: .  dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z
 91: -  gp - global grid point number in this direction

 93:    Output Parameters:
 94: +  newvec - the new vector that can hold the values (size zero on all processes except process 0)
 95: -  scatter - the VecScatter that will map from the original vector to the slice

 97:    Level: advanced

 99:    Notes:
100:    All processors that share the DMDA must call this with the same gp value

102: .keywords: distributed array, get, processor subset
103: @*/
104: PetscErrorCode  DMDAGetRay(DM da,DMDADirection dir,PetscInt gp,Vec *newvec,VecScatter *scatter)
105: {
106:   PetscMPIInt    rank;
107:   DM_DA          *dd = (DM_DA*)da->data;
109:   IS             is;
110:   AO             ao;
111:   Vec            vec;
112:   PetscInt       *indices,i,j;

115:   if (da->dim == 3) SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get slice from 3d DMDA");
116:   MPI_Comm_rank(PetscObjectComm((PetscObject) da), &rank);
117:   DMDAGetAO(da, &ao);
118:   if (!rank) {
119:     if (da->dim == 1) {
120:       if (dir == DMDA_X) {
121:         PetscMalloc1(dd->w, &indices);
122:         indices[0] = dd->w*gp;
123:         for (i = 1; i < dd->w; ++i) indices[i] = indices[i-1] + 1;
124:         AOApplicationToPetsc(ao, dd->w, indices);
125:         VecCreate(PETSC_COMM_SELF, newvec);
126:         VecSetBlockSize(*newvec, dd->w);
127:         VecSetSizes(*newvec, dd->w, PETSC_DETERMINE);
128:         VecSetType(*newvec, VECSEQ);
129:         ISCreateGeneral(PETSC_COMM_SELF, dd->w, indices, PETSC_OWN_POINTER, &is);
130:       } else if (dir == DMDA_Y) SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get Y slice from 1d DMDA");
131:       else SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_ARG_OUTOFRANGE, "Unknown DMDADirection");
132:     } else {
133:       if (dir == DMDA_Y) {
134:         PetscMalloc1(dd->w*dd->M,&indices);
135:         indices[0] = gp*dd->M*dd->w;
136:         for (i=1; i<dd->M*dd->w; i++) indices[i] = indices[i-1] + 1;

138:         AOApplicationToPetsc(ao,dd->M*dd->w,indices);
139:         VecCreate(PETSC_COMM_SELF,newvec);
140:         VecSetBlockSize(*newvec,dd->w);
141:         VecSetSizes(*newvec,dd->M*dd->w,PETSC_DETERMINE);
142:         VecSetType(*newvec,VECSEQ);
143:         ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->M,indices,PETSC_OWN_POINTER,&is);
144:       } else if (dir == DMDA_X) {
145:         PetscMalloc1(dd->w*dd->N,&indices);
146:         indices[0] = dd->w*gp;
147:         for (j=1; j<dd->w; j++) indices[j] = indices[j-1] + 1;
148:         for (i=1; i<dd->N; i++) {
149:           indices[i*dd->w] = indices[i*dd->w-1] + dd->w*dd->M - dd->w + 1;
150:           for (j=1; j<dd->w; j++) indices[i*dd->w + j] = indices[i*dd->w + j - 1] + 1;
151:         }
152:         AOApplicationToPetsc(ao,dd->w*dd->N,indices);
153:         VecCreate(PETSC_COMM_SELF,newvec);
154:         VecSetBlockSize(*newvec,dd->w);
155:         VecSetSizes(*newvec,dd->N*dd->w,PETSC_DETERMINE);
156:         VecSetType(*newvec,VECSEQ);
157:         ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->N,indices,PETSC_OWN_POINTER,&is);
158:       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown DMDADirection");
159:     }
160:   } else {
161:     VecCreateSeq(PETSC_COMM_SELF, 0, newvec);
162:     ISCreateGeneral(PETSC_COMM_SELF, 0, 0, PETSC_COPY_VALUES, &is);
163:   }
164:   DMGetGlobalVector(da, &vec);
165:   VecScatterCreate(vec, is, *newvec, NULL, scatter);
166:   DMRestoreGlobalVector(da, &vec);
167:   ISDestroy(&is);
168:   return(0);
169: }

173: /*@C
174:    DMDAGetProcessorSubset - Returns a communicator consisting only of the
175:    processors in a DMDA that own a particular global x, y, or z grid point
176:    (corresponding to a logical plane in a 3D grid or a line in a 2D grid).

178:    Collective on DMDA

180:    Input Parameters:
181: +  da - the distributed array
182: .  dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z
183: -  gp - global grid point number in this direction

185:    Output Parameters:
186: .  comm - new communicator

188:    Level: advanced

190:    Notes:
191:    All processors that share the DMDA must call this with the same gp value

193:    This routine is particularly useful to compute boundary conditions
194:    or other application-specific calculations that require manipulating
195:    sets of data throughout a logical plane of grid points.

197: .keywords: distributed array, get, processor subset
198: @*/
199: PetscErrorCode  DMDAGetProcessorSubset(DM da,DMDADirection dir,PetscInt gp,MPI_Comm *comm)
200: {
201:   MPI_Group      group,subgroup;
203:   PetscInt       i,ict,flag,*owners,xs,xm,ys,ym,zs,zm;
204:   PetscMPIInt    size,*ranks = NULL;
205:   DM_DA          *dd = (DM_DA*)da->data;

209:   flag = 0;
210:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
211:   MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);
212:   if (dir == DMDA_Z) {
213:     if (da->dim < 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Z invalid for DMDA dim < 3");
214:     if (gp < 0 || gp > dd->P) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
215:     if (gp >= zs && gp < zs+zm) flag = 1;
216:   } else if (dir == DMDA_Y) {
217:     if (da->dim == 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Y invalid for DMDA dim = 1");
218:     if (gp < 0 || gp > dd->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
219:     if (gp >= ys && gp < ys+ym) flag = 1;
220:   } else if (dir == DMDA_X) {
221:     if (gp < 0 || gp > dd->M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
222:     if (gp >= xs && gp < xs+xm) flag = 1;
223:   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");

225:   PetscMalloc2(size,&owners,size,&ranks);
226:   MPI_Allgather(&flag,1,MPIU_INT,owners,1,MPIU_INT,PetscObjectComm((PetscObject)da));
227:   ict  = 0;
228:   PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",da->dim,(int)dir);
229:   for (i=0; i<size; i++) {
230:     if (owners[i]) {
231:       ranks[ict] = i; ict++;
232:       PetscInfo1(da,"%D ",i);
233:     }
234:   }
235:   PetscInfo(da,"\n");
236:   MPI_Comm_group(PetscObjectComm((PetscObject)da),&group);
237:   MPI_Group_incl(group,ict,ranks,&subgroup);
238:   MPI_Comm_create(PetscObjectComm((PetscObject)da),subgroup,comm);
239:   MPI_Group_free(&subgroup);
240:   MPI_Group_free(&group);
241:   PetscFree2(owners,ranks);
242:   return(0);
243: }

247: /*@C
248:    DMDAGetProcessorSubsets - Returns communicators consisting only of the
249:    processors in a DMDA adjacent in a particular dimension,
250:    corresponding to a logical plane in a 3D grid or a line in a 2D grid.

252:    Collective on DMDA

254:    Input Parameters:
255: +  da - the distributed array
256: -  dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z

258:    Output Parameters:
259: .  subcomm - new communicator

261:    Level: advanced

263:    Notes:
264:    This routine is useful for distributing one-dimensional data in a tensor product grid.

266: .keywords: distributed array, get, processor subset
267: @*/
268: PetscErrorCode  DMDAGetProcessorSubsets(DM da, DMDADirection dir, MPI_Comm *subcomm)
269: {
270:   MPI_Comm       comm;
271:   MPI_Group      group, subgroup;
272:   PetscInt       subgroupSize = 0;
273:   PetscInt       *firstPoints;
274:   PetscMPIInt    size, *subgroupRanks = NULL;
275:   PetscInt       xs, xm, ys, ym, zs, zm, firstPoint, p;

280:   PetscObjectGetComm((PetscObject)da,&comm);
281:   DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm);
282:   MPI_Comm_size(comm, &size);
283:   if (dir == DMDA_Z) {
284:     if (da->dim < 3) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Z invalid for DMDA dim < 3");
285:     firstPoint = zs;
286:   } else if (dir == DMDA_Y) {
287:     if (da->dim == 1) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Y invalid for DMDA dim = 1");
288:     firstPoint = ys;
289:   } else if (dir == DMDA_X) {
290:     firstPoint = xs;
291:   } else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");

293:   PetscMalloc2(size, &firstPoints, size, &subgroupRanks);
294:   MPI_Allgather(&firstPoint, 1, MPIU_INT, firstPoints, 1, MPIU_INT, comm);
295:   PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",da->dim,(int)dir);
296:   for (p = 0; p < size; ++p) {
297:     if (firstPoints[p] == firstPoint) {
298:       subgroupRanks[subgroupSize++] = p;
299:       PetscInfo1(da, "%D ", p);
300:     }
301:   }
302:   PetscInfo(da, "\n");
303:   MPI_Comm_group(comm, &group);
304:   MPI_Group_incl(group, subgroupSize, subgroupRanks, &subgroup);
305:   MPI_Comm_create(comm, subgroup, subcomm);
306:   MPI_Group_free(&subgroup);
307:   MPI_Group_free(&group);
308:   PetscFree2(firstPoints, subgroupRanks);
309:   return(0);
310: }