Actual source code: plexcgns.c

petsc-3.4.5 2014-06-29
  1: #define PETSCDM_DLL
  2: #include <petsc-private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/

  4: #if defined(PETSC_HAVE_CGNS)
  5: #include <cgnslib.h>
  6: #include <cgns_io.h>
  7: #endif

 11: /*@
 12:   DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file.

 14:   Collective on comm

 16:   Input Parameters:
 17: + comm  - The MPI communicator
 18: . cgid - The CG id associated with a file and obtained using cg_open
 19: - interpolate - Create faces and edges in the mesh

 21:   Output Parameter:
 22: . dm  - The DM object representing the mesh

 24:   Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html

 26:   Level: beginner

 28: .keywords: mesh,CGNS
 29: .seealso: DMPlexCreate(), DMPlexCreateExodus()
 30: @*/
 31: PetscErrorCode DMPlexCreateCGNS(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm)
 32: {
 33: #if defined(PETSC_HAVE_CGNS)
 34:   PetscMPIInt    num_proc, rank;
 35:   PetscSection   coordSection;
 36:   Vec            coordinates;
 37:   PetscScalar    *coords;
 38:   PetscInt       coordSize, v;
 40:   /* Read from file */
 41:   char basename[CGIO_MAX_NAME_LENGTH+1];
 42:   char buffer[CGIO_MAX_NAME_LENGTH+1];
 43:   int  dim    = 0, physDim = 0, numVertices = 0, numCells = 0;
 44:   int  nzones = 0;
 45: #endif

 48: #if defined(PETSC_HAVE_CGNS)
 49:   MPI_Comm_rank(comm, &rank);
 50:   MPI_Comm_size(comm, &num_proc);
 51:   DMCreate(comm, dm);
 52:   DMSetType(*dm, DMPLEX);
 53:   /* Open CGNS II file and read basic informations on rank 0, then broadcast to all processors */
 54:   if (!rank) {
 55:     int nbases, z;

 57:     cg_nbases(cgid, &nbases);
 58:     if (nbases > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single base, not %d\n",nbases);
 59:     cg_base_read(cgid, 1, basename, &dim, &physDim);
 60:     cg_nzones(cgid, 1, &nzones);
 61:     for (z = 1; z <= nzones; ++z) {
 62:       cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */

 64:       cg_zone_read(cgid, 1, z, buffer, sizes);
 65:       numVertices += sizes[0];
 66:       numCells    += sizes[1];
 67:     }
 68:   }
 69:   MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH+1, MPI_CHAR, 0, comm);
 70:   MPI_Bcast(&dim, 1, MPI_INT, 0, comm);
 71:   MPI_Bcast(&nzones, 1, MPI_INT, 0, comm);
 72:   PetscObjectSetName((PetscObject) *dm, basename);
 73:   DMPlexSetDimension(*dm, dim);
 74:   DMPlexSetChart(*dm, 0, numCells+numVertices);

 76:   /* Read zone information */
 77:   if (!rank) {
 78:     int z, c, c_loc, v, v_loc;

 80:     /* Read the cell set connectivity table and build mesh topology
 81:        CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */
 82:     /* First set sizes */
 83:     for (z = 1, c = 0; z <= nzones; ++z) {
 84:       ZoneType_t    zonetype;
 85:       int           nsections;
 86:       ElementType_t cellType;
 87:       cgsize_t      start, end;
 88:       int           nbndry, parentFlag;
 89:       PetscInt      numCorners;

 91:       cg_zone_type(cgid, 1, z, &zonetype);
 92:       if (zonetype == Structured) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Can only handle Unstructured zones for CGNS");
 93:       cg_nsections(cgid, 1, z, &nsections);
 94:       if (nsections > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single section, not %d\n",nsections);
 95:       cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);
 96:       /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
 97:       if (cellType == MIXED) {
 98:         cgsize_t elementDataSize, *elements;
 99:         PetscInt off;

101:         cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);
102:         PetscMalloc(elementDataSize * sizeof(cgsize_t), &elements);
103:         cg_elements_read(cgid, 1, z, 1, elements, NULL);
104:         for (c_loc = start, off = 0; c_loc < end; ++c_loc, ++c) {
105:           cellType = elements[off];
106:           switch (cellType) {
107:           case TRI_3:   numCorners = 3;break;
108:           case QUAD_4:  numCorners = 4;break;
109:           case TETRA_4: numCorners = 4;break;
110:           case HEXA_8:  numCorners = 8;break;
111:           default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
112:           }
113:           DMPlexSetConeSize(*dm, c, numCorners);
114:           off += numCorners+1;
115:         }
116:         PetscFree(elements);
117:       } else {
118:         switch (cellType) {
119:         case TRI_3:   numCorners = 3;break;
120:         case QUAD_4:  numCorners = 4;break;
121:         case TETRA_4: numCorners = 4;break;
122:         case HEXA_8:  numCorners = 8;break;
123:         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
124:         }
125:         for (c_loc = start; c_loc < end; ++c_loc, ++c) {
126:           DMPlexSetConeSize(*dm, c, numCorners);
127:         }
128:       }
129:     }
130:     DMSetUp(*dm);
131:     for (z = 1, c = 0; z <= nzones; ++z) {
132:       ElementType_t cellType;
133:       cgsize_t     *elements, elementDataSize, start, end;
134:       int           nbndry, parentFlag;
135:       PetscInt     *cone, numc, numCorners, maxCorners = 27;

137:       cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);
138:       numc = end - start;
139:       /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
140:       cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);
141:       PetscMalloc2(elementDataSize,cgsize_t,&elements,maxCorners,PetscInt,&cone);
142:       cg_elements_read(cgid, 1, z, 1, elements, NULL);
143:       if (cellType == MIXED) {
144:         /* CGNS uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */
145:         for (c_loc = 0, v = 0; c_loc < numc; ++c_loc, ++c) {
146:           cellType = elements[v]; ++v;
147:           switch (cellType) {
148:           case TRI_3:   numCorners = 3;break;
149:           case QUAD_4:  numCorners = 4;break;
150:           case TETRA_4: numCorners = 4;break;
151:           case HEXA_8:  numCorners = 8;break;
152:           default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
153:           }
154:           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) {
155:             cone[v_loc] = elements[v]+numCells-1;
156:           }
157:           /* Tetrahedra are inverted */
158:           if (cellType == TETRA_4) {
159:             PetscInt tmp = cone[0];
160:             cone[0] = cone[1];
161:             cone[1] = tmp;
162:           }
163:           /* Hexahedra are inverted */
164:           if (cellType == HEXA_8) {
165:             PetscInt tmp = cone[1];
166:             cone[1] = cone[3];
167:             cone[3] = tmp;
168:           }
169:           DMPlexSetCone(*dm, c, cone);
170:           DMPlexSetLabelValue(*dm, "zone", c, z);
171:         }
172:       } else {
173:         switch (cellType) {
174:         case TRI_3:   numCorners = 3;break;
175:         case QUAD_4:  numCorners = 4;break;
176:         case TETRA_4: numCorners = 4;break;
177:         case HEXA_8:  numCorners = 8;break;
178:         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
179:         }

181:         /* CGNS uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */
182:         for (c_loc = 0, v = 0; c_loc < numc; ++c_loc, ++c) {
183:           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) {
184:             cone[v_loc] = elements[v]+numCells-1;
185:           }
186:           /* Tetrahedra are inverted */
187:           if (cellType == TETRA_4) {
188:             PetscInt tmp = cone[0];
189:             cone[0] = cone[1];
190:             cone[1] = tmp;
191:           }
192:           /* Hexahedra are inverted */
193:           if (cellType == HEXA_8) {
194:             PetscInt tmp = cone[1];
195:             cone[1] = cone[3];
196:             cone[3] = tmp;
197:           }
198:           DMPlexSetCone(*dm, c, cone);
199:           DMPlexSetLabelValue(*dm, "zone", c, z);
200:         }
201:       }
202:       PetscFree2(elements,cone);
203:     }
204:   }
205:   DMPlexSymmetrize(*dm);
206:   DMPlexStratify(*dm);
207:   if (interpolate) {
208:     DM idm;

210:     DMPlexInterpolate(*dm, &idm);
211:     /* Maintain zone label */
212:     {
213:       DMLabel label;

215:       DMPlexRemoveLabel(*dm, "zone", &label);
216:       if (label) {DMPlexAddLabel(idm, label);}
217:     }
218:     DMDestroy(dm);
219:     *dm  = idm;
220:   }

222:   /* Read coordinates */
223:   DMPlexGetCoordinateSection(*dm, &coordSection);
224:   PetscSectionSetNumFields(coordSection, 1);
225:   PetscSectionSetFieldComponents(coordSection, 0, dim);
226:   PetscSectionSetChart(coordSection, numCells, numCells + numVertices);
227:   for (v = numCells; v < numCells+numVertices; ++v) {
228:     PetscSectionSetDof(coordSection, v, dim);
229:     PetscSectionSetFieldDof(coordSection, v, 0, dim);
230:   }
231:   PetscSectionSetUp(coordSection);
232:   PetscSectionGetStorageSize(coordSection, &coordSize);
233:   VecCreate(comm, &coordinates);
234:   PetscObjectSetName((PetscObject) coordinates, "coordinates");
235:   VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
236:   VecSetFromOptions(coordinates);
237:   VecGetArray(coordinates, &coords);
238:   if (!rank) {
239:     PetscInt off = 0;
240:     float   *x[3];
241:     int      z, c, d;

243:     PetscMalloc3(numVertices,float,&x[0],numVertices,float,&x[1],numVertices,float,&x[2]);
244:     for (z = 1, c = 0; z <= nzones; ++z) {
245:       DataType_t datatype;
246:       cgsize_t   sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
247:       cgsize_t   range_min[3] = {1, 1, 1};
248:       cgsize_t   range_max[3] = {1, 1, 1};
249:       int        ngrids, ncoords;


252:       cg_zone_read(cgid, 1, z, buffer, sizes);
253:       range_max[0] = sizes[0];
254:       cg_ngrids(cgid, 1, z, &ngrids);
255:       if (ngrids > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single grid, not %d\n",ngrids);
256:       cg_ncoords(cgid, 1, z, &ncoords);
257:       if (ncoords != dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a coordinate array for each dimension, not %d\n",ncoords);
258:       for (d = 0; d < dim; ++d) {
259:         cg_coord_info(cgid, 1, z, 1, &datatype, buffer);
260:         cg_coord_read(cgid, 1, z, buffer, RealSingle, range_min, range_max, x[d]);
261:       }
262:       if (dim > 0) {
263:         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+0] = x[0][v];
264:       }
265:       if (dim > 1) {
266:         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+1] = x[1][v];
267:       }
268:       if (dim > 2) {
269:         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+2] = x[2][v];
270:       }
271:       off += sizes[0];
272:     }
273:     PetscFree3(x[0],x[1],x[2]);
274:   }
275:   VecRestoreArray(coordinates, &coords);
276:   DMSetCoordinatesLocal(*dm, coordinates);
277:   VecDestroy(&coordinates);
278: #else
279:   SETERRQ(comm, PETSC_ERR_SUP, "This method requires CGNS support. Reconfigure using --with-cgns-dir");
280: #endif
281:   return(0);
282: }