Actual source code: dacorn.c
petsc-3.6.1 2015-08-06
2: /*
3: Code for manipulating distributed regular arrays in parallel.
4: */
6: #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/
10: PetscErrorCode DMCreateCoordinateDM_DA(DM dm, DM *cdm)
11: {
14: DMDAGetReducedDMDA(dm,dm->dim,cdm);
15: return(0);
16: }
20: /*@C
21: DMDASetFieldName - Sets the names of individual field components in multicomponent
22: vectors associated with a DMDA.
24: Not Collective
26: Input Parameters:
27: + da - the distributed array
28: . nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the
29: number of degrees of freedom per node within the DMDA
30: - names - the name of the field (component)
32: Level: intermediate
34: .keywords: distributed array, get, component name
36: .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldNames()
37: @*/
38: PetscErrorCode DMDASetFieldName(DM da,PetscInt nf,const char name[])
39: {
41: DM_DA *dd = (DM_DA*)da->data;
45: if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
46: PetscFree(dd->fieldname[nf]);
47: PetscStrallocpy(name,&dd->fieldname[nf]);
48: return(0);
49: }
53: /*@C
54: DMDAGetFieldNames - Gets the name of each component in the vector associated with the DMDA
56: Collective on TS
58: Input Parameter:
59: . dm - the DMDA object
61: Output Parameter:
62: . names - the names of the components, final string is NULL, will have the same number of entries as the dof used in creating the DMDA
64: Level: intermediate
66: .keywords: distributed array, get, component name
68: .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldName(), DMDASetFieldNames()
69: @*/
70: PetscErrorCode DMDAGetFieldNames(DM da,const char * const **names)
71: {
72: DM_DA *dd = (DM_DA*)da->data;
75: *names = (const char * const *) dd->fieldname;
76: return(0);
77: }
81: /*@C
82: DMDASetFieldNames - Sets the name of each component in the vector associated with the DMDA
84: Collective on TS
86: Input Parameters:
87: + dm - the DMDA object
88: - names - the names of the components, final string must be NULL, must have the same number of entries as the dof used in creating the DMDA
90: Level: intermediate
92: .keywords: distributed array, get, component name
94: .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldName()
95: @*/
96: PetscErrorCode DMDASetFieldNames(DM da,const char * const *names)
97: {
98: PetscErrorCode ierr;
99: DM_DA *dd = (DM_DA*)da->data;
102: PetscStrArrayDestroy(&dd->fieldname);
103: PetscStrArrayallocpy(names,&dd->fieldname);
104: return(0);
105: }
109: /*@C
110: DMDAGetFieldName - Gets the names of individual field components in multicomponent
111: vectors associated with a DMDA.
113: Not Collective
115: Input Parameter:
116: + da - the distributed array
117: - nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the
118: number of degrees of freedom per node within the DMDA
120: Output Parameter:
121: . names - the name of the field (component)
123: Level: intermediate
125: .keywords: distributed array, get, component name
127: .seealso: DMDASetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName()
128: @*/
129: PetscErrorCode DMDAGetFieldName(DM da,PetscInt nf,const char **name)
130: {
131: DM_DA *dd = (DM_DA*)da->data;
136: if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
137: *name = dd->fieldname[nf];
138: return(0);
139: }
143: /*@C
144: DMDASetCoordinateName - Sets the name of the coordinate directions associated with a DMDA, for example "x" or "y"
146: Not Collective
148: Input Parameters:
149: + dm - the DM
150: . nf - coordinate number for the DMDA (0, 1, ... dim-1),
151: - name - the name of the coordinate
153: Level: intermediate
155: .keywords: distributed array, get, component name
157: .seealso: DMDAGetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName()
158: @*/
159: PetscErrorCode DMDASetCoordinateName(DM dm,PetscInt nf,const char name[])
160: {
162: DM_DA *dd = (DM_DA*)dm->data;
166: if (nf < 0 || nf >= dm->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf);
167: PetscFree(dd->coordinatename[nf]);
168: PetscStrallocpy(name,&dd->coordinatename[nf]);
169: return(0);
170: }
174: /*@C
175: DMDAGetCoordinateName - Gets the name of a coodinate direction associated with a DMDA.
177: Not Collective
179: Input Parameter:
180: + dm - the DM
181: - nf - number for the DMDA (0, 1, ... dim-1)
183: Output Parameter:
184: . names - the name of the coordinate direction
186: Level: intermediate
188: .keywords: distributed array, get, component name
190: .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName()
191: @*/
192: PetscErrorCode DMDAGetCoordinateName(DM dm,PetscInt nf,const char **name)
193: {
194: DM_DA *dd = (DM_DA*)dm->data;
199: if (nf < 0 || nf >= dm->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf);
200: *name = dd->coordinatename[nf];
201: return(0);
202: }
206: /*@
207: DMDAGetCorners - Returns the global (x,y,z) indices of the lower left
208: corner of the local region, excluding ghost points.
210: Not Collective
212: Input Parameter:
213: . da - the distributed array
215: Output Parameters:
216: + x,y,z - the corner indices (where y and z are optional; these are used
217: for 2D and 3D problems)
218: - m,n,p - widths in the corresponding directions (where n and p are optional;
219: these are used for 2D and 3D problems)
221: Note:
222: The corner information is independent of the number of degrees of
223: freedom per node set with the DMDACreateXX() routine. Thus the x, y, z, and
224: m, n, p can be thought of as coordinates on a logical grid, where each
225: grid point has (potentially) several degrees of freedom.
226: Any of y, z, n, and p can be passed in as NULL if not needed.
228: Level: beginner
230: .keywords: distributed array, get, corners, nodes, local indices
232: .seealso: DMDAGetGhostCorners(), DMDAGetOwnershipRanges()
233: @*/
234: PetscErrorCode DMDAGetCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p)
235: {
236: PetscInt w;
237: DM_DA *dd = (DM_DA*)da->data;
241: /* since the xs, xe ... have all been multiplied by the number of degrees
242: of freedom per cell, w = dd->w, we divide that out before returning.*/
243: w = dd->w;
244: if (x) *x = dd->xs/w + dd->xo; if (m) *m = (dd->xe - dd->xs)/w;
245: /* the y and z have NOT been multiplied by w */
246: if (y) *y = dd->ys + dd->yo; if (n) *n = (dd->ye - dd->ys);
247: if (z) *z = dd->zs + dd->zo; if (p) *p = (dd->ze - dd->zs);
248: return(0);
249: }
253: /*@
254: DMDAGetLocalBoundingBox - Returns the local bounding box for the DMDA.
256: Not Collective
258: Input Parameter:
259: . dm - the DM
261: Output Parameters:
262: + lmin - local minimum coordinates (length dim, optional)
263: - lmax - local maximim coordinates (length dim, optional)
265: Level: beginner
267: .keywords: distributed array, get, coordinates
269: .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetBoundingBox()
270: @*/
271: PetscErrorCode DMDAGetLocalBoundingBox(DM dm,PetscReal lmin[],PetscReal lmax[])
272: {
273: PetscErrorCode ierr;
274: Vec coords = NULL;
275: PetscInt dim,i,j;
276: const PetscScalar *local_coords;
277: PetscReal min[3]={PETSC_MAX_REAL,PETSC_MAX_REAL,PETSC_MAX_REAL},max[3]={PETSC_MIN_REAL,PETSC_MIN_REAL,PETSC_MIN_REAL};
278: PetscInt N,Ni;
282: dim = dm->dim;
283: DMGetCoordinates(dm,&coords);
284: if (coords) {
285: VecGetArrayRead(coords,&local_coords);
286: VecGetLocalSize(coords,&N);
287: Ni = N/dim;
288: for (i=0; i<Ni; i++) {
289: for (j=0; j<3; j++) {
290: min[j] = j < dim ? PetscMin(min[j],PetscRealPart(local_coords[i*dim+j])) : 0;
291: max[j] = j < dim ? PetscMax(max[j],PetscRealPart(local_coords[i*dim+j])) : 0;
292: }
293: }
294: VecRestoreArrayRead(coords,&local_coords);
295: } else { /* Just use grid indices */
296: DMDALocalInfo info;
297: DMDAGetLocalInfo(dm,&info);
298: min[0] = info.xs;
299: min[1] = info.ys;
300: min[2] = info.zs;
301: max[0] = info.xs + info.xm-1;
302: max[1] = info.ys + info.ym-1;
303: max[2] = info.zs + info.zm-1;
304: }
305: if (lmin) {PetscMemcpy(lmin,min,dim*sizeof(PetscReal));}
306: if (lmax) {PetscMemcpy(lmax,max,dim*sizeof(PetscReal));}
307: return(0);
308: }
312: /*@
313: DMDAGetBoundingBox - Returns the global bounding box for the DMDA.
315: Collective on DMDA
317: Input Parameter:
318: . dm - the DM
320: Output Parameters:
321: + gmin - global minimum coordinates (length dim, optional)
322: - gmax - global maximim coordinates (length dim, optional)
324: Level: beginner
326: .keywords: distributed array, get, coordinates
328: .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetLocalBoundingBox()
329: @*/
330: PetscErrorCode DMDAGetBoundingBox(DM dm,PetscReal gmin[],PetscReal gmax[])
331: {
333: PetscMPIInt count;
334: PetscReal lmin[3],lmax[3];
338: PetscMPIIntCast(dm->dim,&count);
339: DMDAGetLocalBoundingBox(dm,lmin,lmax);
340: if (gmin) {MPI_Allreduce(lmin,gmin,count,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)dm));}
341: if (gmax) {MPI_Allreduce(lmax,gmax,count,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)dm));}
342: return(0);
343: }
347: /*@
348: DMDAGetReducedDMDA - Gets the DMDA with the same layout but with fewer or more fields
350: Collective on DMDA
352: Input Parameters:
353: + da - the distributed array
354: - nfields - number of fields in new DMDA
356: Output Parameter:
357: . nda - the new DMDA
359: Level: intermediate
361: .keywords: distributed array, get, corners, nodes, local indices, coordinates
363: .seealso: DMDAGetGhostCorners(), DMSetCoordinates(), DMDASetUniformCoordinates(), DMGetCoordinates(), DMDAGetGhostedCoordinates()
364: @*/
365: PetscErrorCode DMDAGetReducedDMDA(DM da,PetscInt nfields,DM *nda)
366: {
367: PetscErrorCode ierr;
368: DM_DA *dd = (DM_DA*)da->data;
369: PetscInt s,m,n,p,M,N,P,dim,Mo,No,Po;
370: const PetscInt *lx,*ly,*lz;
371: DMBoundaryType bx,by,bz;
372: DMDAStencilType stencil_type;
373: PetscInt ox,oy,oz;
374: PetscInt cl,rl;
377: dim = da->dim;
378: M = dd->M;
379: N = dd->N;
380: P = dd->P;
381: m = dd->m;
382: n = dd->n;
383: p = dd->p;
384: s = dd->s;
385: bx = dd->bx;
386: by = dd->by;
387: bz = dd->bz;
389: stencil_type = dd->stencil_type;
391: DMDAGetOwnershipRanges(da,&lx,&ly,&lz);
392: if (dim == 1) {
393: DMDACreate1d(PetscObjectComm((PetscObject)da),bx,M,nfields,s,dd->lx,nda);
394: } else if (dim == 2) {
395: DMDACreate2d(PetscObjectComm((PetscObject)da),bx,by,stencil_type,M,N,m,n,nfields,s,lx,ly,nda);
396: } else if (dim == 3) {
397: DMDACreate3d(PetscObjectComm((PetscObject)da),bx,by,bz,stencil_type,M,N,P,m,n,p,nfields,s,lx,ly,lz,nda);
398: }
399: if (da->coordinates) {
400: PetscObjectReference((PetscObject)da->coordinates);
402: (*nda)->coordinates = da->coordinates;
403: }
405: /* allow for getting a reduced DA corresponding to a domain decomposition */
406: DMDAGetOffset(da,&ox,&oy,&oz,&Mo,&No,&Po);
407: DMDASetOffset(*nda,ox,oy,oz,Mo,No,Po);
409: /* allow for getting a reduced DA corresponding to a coarsened DA */
410: DMGetCoarsenLevel(da,&cl);
411: DMGetRefineLevel(da,&rl);
413: (*nda)->levelup = rl;
414: (*nda)->leveldown = cl;
415: return(0);
416: }
420: /*@C
421: DMDAGetCoordinateArray - Gets an array containing the coordinates of the DMDA
423: Not Collective
425: Input Parameter:
426: . dm - the DM
428: Output Parameter:
429: . xc - the coordinates
431: Level: intermediate
433: .keywords: distributed array, get, component name
435: .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMDARestoreCoordinateArray()
436: @*/
437: PetscErrorCode DMDAGetCoordinateArray(DM dm,void *xc)
438: {
440: DM cdm;
441: Vec x;
445: DMGetCoordinates(dm,&x);
446: DMGetCoordinateDM(dm,&cdm);
447: DMDAVecGetArray(cdm,x,xc);
448: return(0);
449: }
453: /*@C
454: DMDARestoreCoordinateArray - Sets an array containing the coordinates of the DMDA
456: Not Collective
458: Input Parameter:
459: + dm - the DM
460: - xc - the coordinates
462: Level: intermediate
464: .keywords: distributed array, get, component name
466: .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMDAGetCoordinateArray()
467: @*/
468: PetscErrorCode DMDARestoreCoordinateArray(DM dm,void *xc)
469: {
471: DM cdm;
472: Vec x;
476: DMGetCoordinates(dm,&x);
477: DMGetCoordinateDM(dm,&cdm);
478: DMDAVecRestoreArray(cdm,x,xc);
479: return(0);
480: }