Actual source code: dacorn.c

  1: #define PETSCDM_DLL

  3: /*
  4:   Code for manipulating distributed regular arrays in parallel.
  5: */

 7:  #include src/dm/da/daimpl.h

 11: /*@
 12:    DASetCoordinates - Sets into the DA a vector that indicates the 
 13:       coordinates of the local nodes (NOT including ghost nodes).

 15:    Not Collective

 17:    Input Parameter:
 18: +  da - the distributed array
 19: -  c - coordinate vector

 21:    Note:
 22:     The coordinates should NOT include those for all ghost points

 24:      Does NOT increase the reference count of this vector, so caller should NOT
 25:   destroy the vector.

 27:   Level: intermediate

 29: .keywords: distributed array, get, corners, nodes, local indices, coordinates

 31: .seealso: DAGetGhostCorners(), DAGetCoordinates(), DASetUniformCoordinates(). DAGetGhostCoordinates(), DAGetCoordinateDA()
 32: @*/
 33: PetscErrorCode  DASetCoordinates(DA da,Vec c)
 34: {

 40:   if (da->coordinates) {
 41:     VecDestroy(da->coordinates);
 42:   }
 43:   da->coordinates = c;
 44:   VecSetBlockSize(c,da->dim);
 45:   return(0);
 46: }

 50: /*@
 51:    DAGetCoordinates - Gets the node coordinates associated with a DA.

 53:    Not Collective

 55:    Input Parameter:
 56: .  da - the distributed array

 58:    Output Parameter:
 59: .  c - coordinate vector

 61:    Note:
 62:     Each process has only the coordinates for its local nodes (does NOT have the
 63:   coordinates for the ghost nodes).

 65:     For two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
 66:     and (x_0,y_0,z_0,x_1,y_1,z_1...)

 68:     You should not destroy or keep around this vector after the DA is destroyed.

 70:   Level: intermediate

 72: .keywords: distributed array, get, corners, nodes, local indices, coordinates

 74: .seealso: DAGetGhostCorners(), DASetCoordinates(), DASetUniformCoordinates(), DAGetGhostedCoordinates(), DAGetCoordinateDA()
 75: @*/
 76: PetscErrorCode  DAGetCoordinates(DA da,Vec *c)
 77: {
 79: 
 82:   *c = da->coordinates;
 83:   return(0);
 84: }

 88: /*@
 89:    DAGetCoordinateDA - Gets the DA that scatters between global and local DA coordinates

 91:    Collective on DA

 93:    Input Parameter:
 94: .  da - the distributed array

 96:    Output Parameter:
 97: .  dac - coordinate DA

 99:    Note: You should not destroy or keep around this vector after the DA is destroyed.

101:   Level: intermediate

103: .keywords: distributed array, get, corners, nodes, local indices, coordinates

105: .seealso: DAGetGhostCorners(), DASetCoordinates(), DASetUniformCoordinates(), DAGetCoordinates(), DAGetGhostedCoordinates()
106: @*/
107: PetscErrorCode  DAGetCoordinateDA(DA da,DA *cda)
108: {
109:   DAStencilType  st;
111:   PetscMPIInt    size;

114:   if (da->da_coordinates) {
115:     *cda = da->da_coordinates;
116:     return(0);
117:   }
118:   MPI_Comm_size(da->comm,&size);
119:   if (da->dim == 1) {
120:     if (da->w == 1) {
121:       da->da_coordinates = da;
122:     } else {
123:       PetscInt            s,m,*lc,l;
124:       DAPeriodicType pt;
125:       DAGetInfo(da,0,&m,0,0,0,0,0,0,&s,&pt,0);
126:       DAGetCorners(da,0,0,0,&l,0,0);
127:       PetscMalloc(size*sizeof(PetscInt),&lc);
128:       MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,da->comm);
129:       DACreate1d(da->comm,pt,m,1,s,lc,&da->da_coordinates);
130:       PetscFree(lc);
131:     }
132:   } else if (da->dim == 2) {
133:     DAGetInfo(da,0,0,0,0,0,0,0,0,0,0,&st);
134:     if (da->w == 2 && st == DA_STENCIL_BOX) {
135:       da->da_coordinates = da;
136:     } else {
137:       PetscInt            i,s,m,*lc,*ld,l,k,n,M,N;
138:       DAPeriodicType pt;
139:       DAGetInfo(da,0,&m,&n,0,&M,&N,0,0,&s,&pt,0);
140:       DAGetCorners(da,0,0,0,&l,&k,0);
141:       PetscMalloc(size*sizeof(PetscInt),&lc);
142:       PetscMalloc(size*sizeof(PetscInt),&ld);
143:       /* only first M values in lc matter */
144:       MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,da->comm);
145:       /* every Mth value in ld matters */
146:       MPI_Allgather(&k,1,MPIU_INT,ld,1,MPIU_INT,da->comm);
147:       for ( i=0; i<N; i++) {
148:         ld[i] = ld[M*i];
149:       }
150:       DACreate2d(da->comm,pt,DA_STENCIL_BOX,m,n,M,N,2,s,lc,ld,&da->da_coordinates);
151:       PetscFree(lc);
152:       PetscFree(ld);
153:     }
154:   } else if (da->dim == 3) {
155:     DAGetInfo(da,0,0,0,0,0,0,0,0,0,0,&st);
156:     if (da->w == 3 && st == DA_STENCIL_BOX) {
157:       da->da_coordinates = da;
158:     } else {
159:       PetscInt            i,s,m,*lc,*ld,*le,l,k,q,n,M,N,P,p;
160:       DAPeriodicType pt;
161:       DAGetInfo(da,0,&m,&n,&p,&M,&N,&P,0,&s,&pt,0);
162:       DAGetCorners(da,0,0,0,&l,&k,&q);
163:       PetscMalloc(size*sizeof(PetscInt),&lc);
164:       PetscMalloc(size*sizeof(PetscInt),&ld);
165:       PetscMalloc(size*sizeof(PetscInt),&le);
166:       /* only first M values in lc matter */
167:       MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,da->comm);
168:       /* every Mth value in ld matters */
169:       MPI_Allgather(&k,1,MPIU_INT,ld,1,MPIU_INT,da->comm);
170:       for ( i=0; i<N; i++) {
171:         ld[i] = ld[M*i];
172:       }
173:       MPI_Allgather(&q,1,MPIU_INT,le,1,MPIU_INT,da->comm);
174:       for ( i=0; i<P; i++) {
175:         le[i] = le[M*N*i];
176:       }
177:       DACreate3d(da->comm,pt,DA_STENCIL_BOX,m,n,p,M,N,P,3,s,lc,ld,le,&da->da_coordinates);
178:       PetscFree(lc);
179:       PetscFree(ld);
180:       PetscFree(le);
181:     }
182:   }
183:   *cda = da->da_coordinates;
184:   return(0);
185: }


190: /*@
191:    DAGetGhostedCoordinates - Gets the node coordinates associated with a DA.

193:    Collective on DA the first time it is called

195:    Input Parameter:
196: .  da - the distributed array

198:    Output Parameter:
199: .  c - coordinate vector

201:    Note:
202:     Each process has only the coordinates for its local AND ghost nodes

204:     For two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
205:     and (x_0,y_0,z_0,x_1,y_1,z_1...)

207:     You should not destroy or keep around this vector after the DA is destroyed.

209:   Level: intermediate

211: .keywords: distributed array, get, corners, nodes, local indices, coordinates

213: .seealso: DAGetGhostCorners(), DASetCoordinates(), DASetUniformCoordinates(), DAGetCoordinates(), DAGetCoordinateDA()
214: @*/
215: PetscErrorCode  DAGetGhostedCoordinates(DA da,Vec *c)
216: {
218: 
221:   if (!da->coordinates) SETERRQ(PETSC_ERR_ORDER,"You must call DASetCoordinates() before this call");
222:   if (!da->ghosted_coordinates) {
223:     DA  dac;
225:     DAGetCoordinateDA(da,&dac);
226:     DACreateLocalVector(dac,&da->ghosted_coordinates);
227:     if (dac == da) {PetscObjectDereference((PetscObject)dac);}
228:     DAGlobalToLocalBegin(dac,da->coordinates,INSERT_VALUES,da->ghosted_coordinates);
229:     DAGlobalToLocalEnd(dac,da->coordinates,INSERT_VALUES,da->ghosted_coordinates);
230:   }
231:   *c = da->ghosted_coordinates;
232:   return(0);
233: }

237: /*@C
238:    DASetFieldName - Sets the names of individual field components in multicomponent
239:    vectors associated with a DA.

241:    Not Collective

243:    Input Parameters:
244: +  da - the distributed array
245: .  nf - field number for the DA (0, 1, ... dof-1), where dof indicates the 
246:         number of degrees of freedom per node within the DA
247: -  names - the name of the field (component)

249:   Level: intermediate

251: .keywords: distributed array, get, component name

253: .seealso: DAGetFieldName()
254: @*/
255: PetscErrorCode  DASetFieldName(DA da,PetscInt nf,const char name[])
256: {

260: 
262:   if (nf < 0 || nf >= da->w) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
263:   if (da->fieldname[nf]) {PetscFree(da->fieldname[nf]);}
264: 
265:   PetscStrallocpy(name,&da->fieldname[nf]);
266:   return(0);
267: }

271: /*@C
272:    DAGetFieldName - Gets the names of individual field components in multicomponent
273:    vectors associated with a DA.

275:    Not Collective

277:    Input Parameter:
278: +  da - the distributed array
279: -  nf - field number for the DA (0, 1, ... dof-1), where dof indicates the 
280:         number of degrees of freedom per node within the DA

282:    Output Parameter:
283: .  names - the name of the field (component)

285:   Level: intermediate

287: .keywords: distributed array, get, component name

289: .seealso: DASetFieldName()
290: @*/
291: PetscErrorCode  DAGetFieldName(DA da,PetscInt nf,char **name)
292: {
294: 
297:   if (nf < 0 || nf >= da->w) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
298:   *name = da->fieldname[nf];
299:   return(0);
300: }

304: /*@
305:    DAGetCorners - Returns the global (x,y,z) indices of the lower left
306:    corner of the local region, excluding ghost points.

308:    Not Collective

310:    Input Parameter:
311: .  da - the distributed array

313:    Output Parameters:
314: +  x,y,z - the corner indices (where y and z are optional; these are used
315:            for 2D and 3D problems)
316: -  m,n,p - widths in the corresponding directions (where n and p are optional;
317:            these are used for 2D and 3D problems)

319:    Note:
320:    The corner information is independent of the number of degrees of 
321:    freedom per node set with the DACreateXX() routine. Thus the x, y, z, and
322:    m, n, p can be thought of as coordinates on a logical grid, where each
323:    grid point has (potentially) several degrees of freedom.
324:    Any of y, z, n, and p can be passed in as PETSC_NULL if not needed.

326:   Level: beginner

328: .keywords: distributed array, get, corners, nodes, local indices

330: .seealso: DAGetGhostCorners()
331: @*/
332: PetscErrorCode  DAGetCorners(DA da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p)
333: {
334:   PetscInt w;

338:   /* since the xs, xe ... have all been multiplied by the number of degrees 
339:      of freedom per cell, w = da->w, we divide that out before returning.*/
340:   w = da->w;
341:   if (x) *x = da->xs/w; if(m) *m = (da->xe - da->xs)/w;
342:   /* the y and z have NOT been multiplied by w */
343:   if (y) *y = da->ys;   if (n) *n = (da->ye - da->ys);
344:   if (z) *z = da->zs;   if (p) *p = (da->ze - da->zs);
345:   return(0);
346: }