Actual source code: plexcgns.c
petsc-3.7.7 2017-09-25
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: /*@C
12: DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file.
14: Collective on comm
16: Input Parameters:
17: + comm - The MPI communicator
18: . filename - The name of the CGNS file
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(), DMPlexCreateCGNS(), DMPlexCreateExodus()
30: @*/
31: PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
32: {
33: PetscMPIInt rank;
35: #if defined(PETSC_HAVE_CGNS)
36: int cgid = -1;
37: #endif
41: MPI_Comm_rank(comm, &rank);
42: #if defined(PETSC_HAVE_CGNS)
43: if (!rank) {
44: cg_open(filename, CG_MODE_READ, &cgid);
45: if (cgid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename);
46: }
47: DMPlexCreateCGNS(comm, cgid, interpolate, dm);
48: if (!rank) {cg_close(cgid);}
49: #else
50: SETERRQ(comm, PETSC_ERR_SUP, "Loading meshes requires CGNS support. Reconfigure using --with-cgns-dir");
51: #endif
52: return(0);
53: }
57: /*@
58: DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file ID.
60: Collective on comm
62: Input Parameters:
63: + comm - The MPI communicator
64: . cgid - The CG id associated with a file and obtained using cg_open
65: - interpolate - Create faces and edges in the mesh
67: Output Parameter:
68: . dm - The DM object representing the mesh
70: Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html
72: Level: beginner
74: .keywords: mesh,CGNS
75: .seealso: DMPlexCreate(), DMPlexCreateExodus()
76: @*/
77: PetscErrorCode DMPlexCreateCGNS(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm)
78: {
79: #if defined(PETSC_HAVE_CGNS)
80: PetscMPIInt num_proc, rank;
81: PetscSection coordSection;
82: Vec coordinates;
83: PetscScalar *coords;
84: PetscInt *cellStart, *vertStart;
85: PetscInt coordSize, v;
87: /* Read from file */
88: char basename[CGIO_MAX_NAME_LENGTH+1];
89: char buffer[CGIO_MAX_NAME_LENGTH+1];
90: int dim = 0, physDim = 0, numVertices = 0, numCells = 0;
91: int nzones = 0;
92: #endif
95: #if defined(PETSC_HAVE_CGNS)
96: MPI_Comm_rank(comm, &rank);
97: MPI_Comm_size(comm, &num_proc);
98: DMCreate(comm, dm);
99: DMSetType(*dm, DMPLEX);
100: /* Open CGNS II file and read basic informations on rank 0, then broadcast to all processors */
101: if (!rank) {
102: int nbases, z;
104: cg_nbases(cgid, &nbases);
105: if (nbases > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single base, not %d\n",nbases);
106: cg_base_read(cgid, 1, basename, &dim, &physDim);
107: cg_nzones(cgid, 1, &nzones);
108: PetscCalloc2(nzones+1, &cellStart, nzones+1, &vertStart);
109: for (z = 1; z <= nzones; ++z) {
110: cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
112: cg_zone_read(cgid, 1, z, buffer, sizes);
113: numVertices += sizes[0];
114: numCells += sizes[1];
115: cellStart[z] += sizes[1] + cellStart[z-1];
116: vertStart[z] += sizes[0] + vertStart[z-1];
117: }
118: for (z = 1; z <= nzones; ++z) {
119: vertStart[z] += numCells;
120: }
121: }
122: MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH+1, MPI_CHAR, 0, comm);
123: MPI_Bcast(&dim, 1, MPI_INT, 0, comm);
124: MPI_Bcast(&nzones, 1, MPI_INT, 0, comm);
125: PetscObjectSetName((PetscObject) *dm, basename);
126: DMSetDimension(*dm, dim);
127: DMPlexSetChart(*dm, 0, numCells+numVertices);
129: /* Read zone information */
130: if (!rank) {
131: int z, c, c_loc, v, v_loc;
133: /* Read the cell set connectivity table and build mesh topology
134: CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */
135: /* First set sizes */
136: for (z = 1, c = 0; z <= nzones; ++z) {
137: ZoneType_t zonetype;
138: int nsections;
139: ElementType_t cellType;
140: cgsize_t start, end;
141: int nbndry, parentFlag;
142: PetscInt numCorners;
144: cg_zone_type(cgid, 1, z, &zonetype);
145: if (zonetype == Structured) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Can only handle Unstructured zones for CGNS");
146: cg_nsections(cgid, 1, z, &nsections);
147: if (nsections > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single section, not %d\n",nsections);
148: cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);
149: /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
150: if (cellType == MIXED) {
151: cgsize_t elementDataSize, *elements;
152: PetscInt off;
154: cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);
155: PetscMalloc1(elementDataSize, &elements);
156: cg_elements_read(cgid, 1, z, 1, elements, NULL);
157: for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) {
158: switch (elements[off]) {
159: case TRI_3: numCorners = 3;break;
160: case QUAD_4: numCorners = 4;break;
161: case TETRA_4: numCorners = 4;break;
162: case HEXA_8: numCorners = 8;break;
163: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[off]);
164: }
165: DMPlexSetConeSize(*dm, c, numCorners);
166: off += numCorners+1;
167: }
168: PetscFree(elements);
169: } else {
170: switch (cellType) {
171: case TRI_3: numCorners = 3;break;
172: case QUAD_4: numCorners = 4;break;
173: case TETRA_4: numCorners = 4;break;
174: case HEXA_8: numCorners = 8;break;
175: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
176: }
177: for (c_loc = start; c_loc <= end; ++c_loc, ++c) {
178: DMPlexSetConeSize(*dm, c, numCorners);
179: }
180: }
181: }
182: DMSetUp(*dm);
183: for (z = 1, c = 0; z <= nzones; ++z) {
184: ElementType_t cellType;
185: cgsize_t *elements, elementDataSize, start, end;
186: int nbndry, parentFlag;
187: PetscInt *cone, numc, numCorners, maxCorners = 27;
189: cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);
190: numc = end - start;
191: /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
192: cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);
193: PetscMalloc2(elementDataSize,&elements,maxCorners,&cone);
194: cg_elements_read(cgid, 1, z, 1, elements, NULL);
195: if (cellType == MIXED) {
196: /* CGNS uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */
197: for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
198: switch (elements[v]) {
199: case TRI_3: numCorners = 3;break;
200: case QUAD_4: numCorners = 4;break;
201: case TETRA_4: numCorners = 4;break;
202: case HEXA_8: numCorners = 8;break;
203: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[v]);
204: }
205: ++v;
206: for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) {
207: cone[v_loc] = elements[v]+numCells-1;
208: }
209: /* Tetrahedra are inverted */
210: if (elements[v] == TETRA_4) {
211: PetscInt tmp = cone[0];
212: cone[0] = cone[1];
213: cone[1] = tmp;
214: }
215: /* Hexahedra are inverted */
216: if (elements[v] == HEXA_8) {
217: PetscInt tmp = cone[5];
218: cone[5] = cone[7];
219: cone[7] = tmp;
220: }
221: DMPlexSetCone(*dm, c, cone);
222: DMSetLabelValue(*dm, "zone", c, z);
223: }
224: } else {
225: switch (cellType) {
226: case TRI_3: numCorners = 3;break;
227: case QUAD_4: numCorners = 4;break;
228: case TETRA_4: numCorners = 4;break;
229: case HEXA_8: numCorners = 8;break;
230: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
231: }
233: /* CGNS uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */
234: for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
235: for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) {
236: cone[v_loc] = elements[v]+numCells-1;
237: }
238: /* Tetrahedra are inverted */
239: if (cellType == TETRA_4) {
240: PetscInt tmp = cone[0];
241: cone[0] = cone[1];
242: cone[1] = tmp;
243: }
244: /* Hexahedra are inverted, and they give the top first */
245: if (cellType == HEXA_8) {
246: PetscInt tmp = cone[5];
247: cone[5] = cone[7];
248: cone[7] = tmp;
249: }
250: DMPlexSetCone(*dm, c, cone);
251: DMSetLabelValue(*dm, "zone", c, z);
252: }
253: }
254: PetscFree2(elements,cone);
255: }
256: }
257: DMPlexSymmetrize(*dm);
258: DMPlexStratify(*dm);
259: if (interpolate) {
260: DM idm = NULL;
262: DMPlexInterpolate(*dm, &idm);
263: /* Maintain zone label */
264: {
265: DMLabel label;
267: DMRemoveLabel(*dm, "zone", &label);
268: if (label) {DMAddLabel(idm, label);}
269: }
270: DMDestroy(dm);
271: *dm = idm;
272: }
274: /* Read coordinates */
275: DMGetCoordinateSection(*dm, &coordSection);
276: PetscSectionSetNumFields(coordSection, 1);
277: PetscSectionSetFieldComponents(coordSection, 0, dim);
278: PetscSectionSetChart(coordSection, numCells, numCells + numVertices);
279: for (v = numCells; v < numCells+numVertices; ++v) {
280: PetscSectionSetDof(coordSection, v, dim);
281: PetscSectionSetFieldDof(coordSection, v, 0, dim);
282: }
283: PetscSectionSetUp(coordSection);
284: PetscSectionGetStorageSize(coordSection, &coordSize);
285: VecCreate(PETSC_COMM_SELF, &coordinates);
286: PetscObjectSetName((PetscObject) coordinates, "coordinates");
287: VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
288: VecSetType(coordinates,VECSTANDARD);
289: VecGetArray(coordinates, &coords);
290: if (!rank) {
291: PetscInt off = 0;
292: float *x[3];
293: int z, d;
295: PetscMalloc3(numVertices,&x[0],numVertices,&x[1],numVertices,&x[2]);
296: for (z = 1; z <= nzones; ++z) {
297: DataType_t datatype;
298: cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
299: cgsize_t range_min[3] = {1, 1, 1};
300: cgsize_t range_max[3] = {1, 1, 1};
301: int ngrids, ncoords;
304: cg_zone_read(cgid, 1, z, buffer, sizes);
305: range_max[0] = sizes[0];
306: cg_ngrids(cgid, 1, z, &ngrids);
307: if (ngrids > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single grid, not %d\n",ngrids);
308: cg_ncoords(cgid, 1, z, &ncoords);
309: if (ncoords != dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a coordinate array for each dimension, not %d\n",ncoords);
310: for (d = 0; d < dim; ++d) {
311: cg_coord_info(cgid, 1, z, 1+d, &datatype, buffer);
312: cg_coord_read(cgid, 1, z, buffer, RealSingle, range_min, range_max, x[d]);
313: }
314: if (dim > 0) {
315: for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+0] = x[0][v];
316: }
317: if (dim > 1) {
318: for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+1] = x[1][v];
319: }
320: if (dim > 2) {
321: for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+2] = x[2][v];
322: }
323: off += sizes[0];
324: }
325: PetscFree3(x[0],x[1],x[2]);
326: }
327: VecRestoreArray(coordinates, &coords);
328: DMSetCoordinatesLocal(*dm, coordinates);
329: VecDestroy(&coordinates);
330: /* Read boundary conditions */
331: if (!rank) {
332: DMLabel label;
333: BCType_t bctype;
334: DataType_t datatype;
335: PointSetType_t pointtype;
336: cgsize_t *points;
337: PetscReal *normals;
338: int normal[3];
339: char *bcname = buffer;
340: cgsize_t npoints, nnormals;
341: int z, nbc, bc, c, ndatasets;
343: for (z = 1; z <= nzones; ++z) {
344: cg_nbocos(cgid, 1, z, &nbc);
345: for (bc = 1; bc <= nbc; ++bc) {
346: cg_boco_info(cgid, 1, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets);
347: DMCreateLabel(*dm, bcname);
348: DMGetLabel(*dm, bcname, &label);
349: PetscMalloc2(npoints, &points, nnormals, &normals);
350: cg_boco_read(cgid, 1, z, bc, points, (void *) normals);
351: if (pointtype == ElementRange) {
352: /* Range of cells: assuming half-open interval since the documentation sucks */
353: for (c = points[0]; c < points[1]; ++c) {
354: DMLabelSetValue(label, c - cellStart[z-1], 1);
355: }
356: } else if (pointtype == ElementList) {
357: /* List of cells */
358: for (c = 0; c < npoints; ++c) {
359: DMLabelSetValue(label, points[c] - cellStart[z-1], 1);
360: }
361: } else if (pointtype == PointRange) {
362: GridLocation_t gridloc;
364: /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
365: cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end");
366: cg_gridlocation_read(&gridloc);
367: /* Range of points: assuming half-open interval since the documentation sucks */
368: for (c = points[0]; c < points[1]; ++c) {
369: if (gridloc == Vertex) {DMLabelSetValue(label, c - vertStart[z-1], 1);}
370: else {DMLabelSetValue(label, c - cellStart[z-1], 1);}
371: }
372: } else if (pointtype == PointList) {
373: GridLocation_t gridloc;
375: /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
376: cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end");
377: cg_gridlocation_read(&gridloc);
378: for (c = 0; c < npoints; ++c) {
379: if (gridloc == Vertex) {DMLabelSetValue(label, points[c] - vertStart[z-1], 1);}
380: else {DMLabelSetValue(label, points[c] - cellStart[z-1], 1);}
381: }
382: } else SETERRQ1(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int) pointtype);
383: PetscFree2(points, normals);
384: }
385: }
386: PetscFree2(cellStart, vertStart);
387: }
388: #else
389: SETERRQ(comm, PETSC_ERR_SUP, "This method requires CGNS support. Reconfigure using --with-cgns-dir");
390: #endif
391: return(0);
392: }