xref: /petsc/src/dm/impls/plex/plexcgns.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1 #define PETSCDM_DLL
2 #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3 
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
15 
16 /*@C
17   DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file.
18 
19   Collective
20 
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
25 
26   Output Parameter:
27 . dm  - The DM object representing the mesh
28 
29   Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html
30 
31   Level: beginner
32 
33 .seealso: DMPlexCreate(), DMPlexCreateCGNS(), DMPlexCreateExodus()
34 @*/
35 PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
36 {
37   PetscMPIInt    rank;
38   PetscErrorCode ierr;
39 #if defined(PETSC_HAVE_CGNS)
40   int cgid = -1;
41 #endif
42 
43   PetscFunctionBegin;
44   PetscValidCharPointer(filename, 2);
45   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
46 #if defined(PETSC_HAVE_CGNS)
47   if (!rank) {
48     ierr = cg_open(filename, CG_MODE_READ, &cgid);CHKERRQ(ierr);
49     if (cgid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename);
50   }
51   ierr = DMPlexCreateCGNS(comm, cgid, interpolate, dm);CHKERRQ(ierr);
52   if (!rank) {ierr = cg_close(cgid);CHKERRQ(ierr);}
53 #else
54   SETERRQ(comm, PETSC_ERR_SUP, "Loading meshes requires CGNS support. Reconfigure using --with-cgns-dir");
55 #endif
56   PetscFunctionReturn(0);
57 }
58 
59 /*@
60   DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file ID.
61 
62   Collective
63 
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
68 
69   Output Parameter:
70 . dm  - The DM object representing the mesh
71 
72   Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html
73 
74   Level: beginner
75 
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;
87   PetscErrorCode ierr;
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
94 
95   PetscFunctionBegin;
96 #if defined(PETSC_HAVE_CGNS)
97   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
98   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
99   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
100   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
101   /* Open CGNS II file and read basic informations on rank 0, then broadcast to all processors */
102   if (!rank) {
103     int nbases, z;
104 
105     ierr = cg_nbases(cgid, &nbases);CHKERRQ(ierr);
106     if (nbases > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single base, not %d\n",nbases);
107     ierr = cg_base_read(cgid, 1, basename, &dim, &physDim);CHKERRQ(ierr);
108     ierr = cg_nzones(cgid, 1, &nzones);CHKERRQ(ierr);
109     ierr = PetscCalloc2(nzones+1, &cellStart, nzones+1, &vertStart);CHKERRQ(ierr);
110     for (z = 1; z <= nzones; ++z) {
111       cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
112 
113       ierr = cg_zone_read(cgid, 1, z, buffer, sizes);CHKERRQ(ierr);
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   ierr = MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
124   ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
125   ierr = MPI_Bcast(&nzones, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
126   ierr = PetscObjectSetName((PetscObject) *dm, basename);CHKERRQ(ierr);
127   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
128   ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
129 
130   /* Read zone information */
131   if (!rank) {
132     int z, c, c_loc, v, v_loc;
133 
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;
144 
145       ierr = cg_zone_type(cgid, 1, z, &zonetype);CHKERRQ(ierr);
146       if (zonetype == CGNS_ENUMV( Structured )) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Can only handle Unstructured zones for CGNS");
147       ierr = cg_nsections(cgid, 1, z, &nsections);CHKERRQ(ierr);
148       if (nsections > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single section, not %d\n",nsections);
149       ierr = cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);CHKERRQ(ierr);
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;
154 
155         ierr = cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);CHKERRQ(ierr);
156         ierr = PetscMalloc1(elementDataSize, &elements);CHKERRQ(ierr);
157         ierr = cg_elements_read(cgid, 1, z, 1, elements, NULL);CHKERRQ(ierr);
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           ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr);
167           off += numCorners+1;
168         }
169         ierr = PetscFree(elements);CHKERRQ(ierr);
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           ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr);
180         }
181       }
182     }
183     ierr = DMSetUp(*dm);CHKERRQ(ierr);
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;
189 
190       ierr = cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);CHKERRQ(ierr);
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       ierr = cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);CHKERRQ(ierr);
194       ierr = PetscMalloc2(elementDataSize,&elements,maxCorners,&cone);CHKERRQ(ierr);
195       ierr = cg_elements_read(cgid, 1, z, 1, elements, NULL);CHKERRQ(ierr);
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           ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
223           ierr = DMSetLabelValue(*dm, "zone", c, z);CHKERRQ(ierr);
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         }
233 
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           ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
252           ierr = DMSetLabelValue(*dm, "zone", c, z);CHKERRQ(ierr);
253         }
254       }
255       ierr = PetscFree2(elements,cone);CHKERRQ(ierr);
256     }
257   }
258   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
259   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
260   if (interpolate) {
261     DM idm;
262 
263     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
264     ierr = DMDestroy(dm);CHKERRQ(ierr);
265     *dm  = idm;
266   }
267 
268   /* Read coordinates */
269   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
270   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
271   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
272   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
273   for (v = numCells; v < numCells+numVertices; ++v) {
274     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
275     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
276   }
277   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
278   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
279   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
280   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
281   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
282   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
283   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
284   if (!rank) {
285     PetscInt off = 0;
286     float   *x[3];
287     int      z, d;
288 
289     ierr = PetscMalloc3(numVertices,&x[0],numVertices,&x[1],numVertices,&x[2]);CHKERRQ(ierr);
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;
296 
297 
298       ierr = cg_zone_read(cgid, 1, z, buffer, sizes);CHKERRQ(ierr);
299       range_max[0] = sizes[0];
300       ierr = cg_ngrids(cgid, 1, z, &ngrids);CHKERRQ(ierr);
301       if (ngrids > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single grid, not %d\n",ngrids);
302       ierr = cg_ncoords(cgid, 1, z, &ncoords);CHKERRQ(ierr);
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         ierr = cg_coord_info(cgid, 1, z, 1+d, &datatype, buffer);CHKERRQ(ierr);
306         ierr = cg_coord_read(cgid, 1, z, buffer, CGNS_ENUMV( RealSingle ), range_min, range_max, x[d]);CHKERRQ(ierr);
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     ierr = PetscFree3(x[0],x[1],x[2]);CHKERRQ(ierr);
320   }
321   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
322   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
323   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
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;
336 
337     for (z = 1; z <= nzones; ++z) {
338       ierr = cg_nbocos(cgid, 1, z, &nbc);CHKERRQ(ierr);
339       for (bc = 1; bc <= nbc; ++bc) {
340         ierr = cg_boco_info(cgid, 1, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets);CHKERRQ(ierr);
341         ierr = DMCreateLabel(*dm, bcname);CHKERRQ(ierr);
342         ierr = DMGetLabel(*dm, bcname, &label);CHKERRQ(ierr);
343         ierr = PetscMalloc2(npoints, &points, nnormals, &normals);CHKERRQ(ierr);
344         ierr = cg_boco_read(cgid, 1, z, bc, points, (void *) normals);CHKERRQ(ierr);
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             ierr = DMLabelSetValue(label, c - cellStart[z-1], 1);CHKERRQ(ierr);
349           }
350         } else if (pointtype == CGNS_ENUMV( ElementList )) {
351           /* List of cells */
352           for (c = 0; c < npoints; ++c) {
353             ierr = DMLabelSetValue(label, points[c] - cellStart[z-1], 1);CHKERRQ(ierr);
354           }
355         } else if (pointtype == CGNS_ENUMV( PointRange )) {
356           CGNS_ENUMT( GridLocation_t ) gridloc;
357 
358           /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
359           ierr = cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end");CHKERRQ(ierr);
360           ierr = cg_gridlocation_read(&gridloc);CHKERRQ(ierr);
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 )) {ierr = DMLabelSetValue(label, c - vertStart[z-1], 1);CHKERRQ(ierr);}
364             else                                 {ierr = DMLabelSetValue(label, c - cellStart[z-1], 1);CHKERRQ(ierr);}
365           }
366         } else if (pointtype == CGNS_ENUMV( PointList )) {
367           CGNS_ENUMT( GridLocation_t ) gridloc;
368 
369           /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
370           ierr = cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end");
371           ierr = cg_gridlocation_read(&gridloc);
372           for (c = 0; c < npoints; ++c) {
373             if (gridloc == CGNS_ENUMV( Vertex )) {ierr = DMLabelSetValue(label, points[c] - vertStart[z-1], 1);CHKERRQ(ierr);}
374             else                   {ierr = DMLabelSetValue(label, points[c] - cellStart[z-1], 1);CHKERRQ(ierr);}
375           }
376         } else SETERRQ1(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int) pointtype);
377         ierr = PetscFree2(points, normals);CHKERRQ(ierr);
378       }
379     }
380     ierr = PetscFree2(cellStart, vertStart);CHKERRQ(ierr);
381   }
382 #else
383   SETERRQ(comm, PETSC_ERR_SUP, "This method requires CGNS support. Reconfigure using --with-cgns-dir");
384 #endif
385   PetscFunctionReturn(0);
386 }
387