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