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