xref: /petsc/src/dm/impls/plex/plexcgns.c (revision 0c8c14c502e79aaeb7ec26e1cd1e476322d9010e)
1 #define PETSCDM_DLL
2 #include <petsc-private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3 
4 #if defined(PETSC_HAVE_CGNS)
5 #include <cgnslib.h>
6 #include <cgns_io.h>
7 #endif
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "DMPlexCreateCGNSFromFile"
11 /*@C
12   DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file.
13 
14   Collective on comm
15 
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
20 
21   Output Parameter:
22 . dm  - The DM object representing the mesh
23 
24   Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html
25 
26   Level: beginner
27 
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;
34   PetscErrorCode ierr;
35 #if defined(PETSC_HAVE_CGNS)
36   int cgid = -1;
37 #endif
38 
39   PetscFunctionBegin;
40   PetscValidCharPointer(filename, 2);
41   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
42 #if defined(PETSC_HAVE_CGNS)
43   if (!rank) {
44     ierr = cg_open(filename, CG_MODE_READ, &cgid);CHKERRQ(ierr);
45     if (cgid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename);
46   }
47   ierr = DMPlexCreateCGNS(comm, cgid, interpolate, dm);CHKERRQ(ierr);
48   if (!rank) {ierr = cg_close(cgid);CHKERRQ(ierr);}
49 #else
50   SETERRQ(comm, PETSC_ERR_SUP, "Loading meshes requires CGNS support. Reconfigure using --with-cgns-dir");
51 #endif
52   PetscFunctionReturn(0);
53 }
54 
55 #undef __FUNCT__
56 #define __FUNCT__ "DMPlexCreateCGNS"
57 /*@
58   DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file ID.
59 
60   Collective on comm
61 
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
66 
67   Output Parameter:
68 . dm  - The DM object representing the mesh
69 
70   Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html
71 
72   Level: beginner
73 
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       coordSize, v;
85   PetscErrorCode ierr;
86   /* Read from file */
87   char basename[CGIO_MAX_NAME_LENGTH+1];
88   char buffer[CGIO_MAX_NAME_LENGTH+1];
89   int  dim    = 0, physDim = 0, numVertices = 0, numCells = 0;
90   int  nzones = 0;
91 #endif
92 
93   PetscFunctionBegin;
94 #if defined(PETSC_HAVE_CGNS)
95   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
96   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
97   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
98   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
99   /* Open CGNS II file and read basic informations on rank 0, then broadcast to all processors */
100   if (!rank) {
101     int nbases, z;
102 
103     ierr = cg_nbases(cgid, &nbases);CHKERRQ(ierr);
104     if (nbases > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single base, not %d\n",nbases);
105     ierr = cg_base_read(cgid, 1, basename, &dim, &physDim);CHKERRQ(ierr);
106     ierr = cg_nzones(cgid, 1, &nzones);CHKERRQ(ierr);
107     for (z = 1; z <= nzones; ++z) {
108       cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
109 
110       ierr = cg_zone_read(cgid, 1, z, buffer, sizes);CHKERRQ(ierr);
111       numVertices += sizes[0];
112       numCells    += sizes[1];
113     }
114   }
115   ierr = MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
116   ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
117   ierr = MPI_Bcast(&nzones, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
118   ierr = PetscObjectSetName((PetscObject) *dm, basename);CHKERRQ(ierr);
119   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
120   ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
121 
122   /* Read zone information */
123   if (!rank) {
124     int z, c, c_loc, v, v_loc;
125 
126     /* Read the cell set connectivity table and build mesh topology
127        CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */
128     /* First set sizes */
129     for (z = 1, c = 0; z <= nzones; ++z) {
130       ZoneType_t    zonetype;
131       int           nsections;
132       ElementType_t cellType;
133       cgsize_t      start, end;
134       int           nbndry, parentFlag;
135       PetscInt      numCorners;
136 
137       ierr = cg_zone_type(cgid, 1, z, &zonetype);CHKERRQ(ierr);
138       if (zonetype == Structured) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Can only handle Unstructured zones for CGNS");
139       ierr = cg_nsections(cgid, 1, z, &nsections);CHKERRQ(ierr);
140       if (nsections > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single section, not %d\n",nsections);
141       ierr = cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);CHKERRQ(ierr);
142       /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
143       if (cellType == MIXED) {
144         cgsize_t elementDataSize, *elements;
145         PetscInt off;
146 
147         ierr = cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);CHKERRQ(ierr);
148         ierr = PetscMalloc1(elementDataSize, &elements);CHKERRQ(ierr);
149         ierr = cg_elements_read(cgid, 1, z, 1, elements, NULL);CHKERRQ(ierr);
150         for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) {
151           switch (elements[off]) {
152           case TRI_3:   numCorners = 3;break;
153           case QUAD_4:  numCorners = 4;break;
154           case TETRA_4: numCorners = 4;break;
155           case HEXA_8:  numCorners = 8;break;
156           default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[off]);
157           }
158           ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr);
159           off += numCorners+1;
160         }
161         ierr = PetscFree(elements);CHKERRQ(ierr);
162       } else {
163         switch (cellType) {
164         case TRI_3:   numCorners = 3;break;
165         case QUAD_4:  numCorners = 4;break;
166         case TETRA_4: numCorners = 4;break;
167         case HEXA_8:  numCorners = 8;break;
168         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
169         }
170         for (c_loc = start; c_loc <= end; ++c_loc, ++c) {
171           ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr);
172         }
173       }
174     }
175     ierr = DMSetUp(*dm);CHKERRQ(ierr);
176     for (z = 1, c = 0; z <= nzones; ++z) {
177       ElementType_t cellType;
178       cgsize_t     *elements, elementDataSize, start, end;
179       int           nbndry, parentFlag;
180       PetscInt     *cone, numc, numCorners, maxCorners = 27;
181 
182       ierr = cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);CHKERRQ(ierr);
183       numc = end - start;
184       /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
185       ierr = cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);CHKERRQ(ierr);
186       ierr = PetscMalloc2(elementDataSize,&elements,maxCorners,&cone);CHKERRQ(ierr);
187       ierr = cg_elements_read(cgid, 1, z, 1, elements, NULL);CHKERRQ(ierr);
188       if (cellType == MIXED) {
189         /* CGNS uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */
190         for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
191           switch (elements[v]) {
192           case TRI_3:   numCorners = 3;break;
193           case QUAD_4:  numCorners = 4;break;
194           case TETRA_4: numCorners = 4;break;
195           case HEXA_8:  numCorners = 8;break;
196           default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[v]);
197           }
198           ++v;
199           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) {
200             cone[v_loc] = elements[v]+numCells-1;
201           }
202           /* Tetrahedra are inverted */
203           if (elements[v] == TETRA_4) {
204             PetscInt tmp = cone[0];
205             cone[0] = cone[1];
206             cone[1] = tmp;
207           }
208           /* Hexahedra are inverted */
209           if (elements[v] == HEXA_8) {
210             PetscInt tmp = cone[1];
211             cone[1] = cone[3];
212             cone[3] = tmp;
213           }
214           ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
215           ierr = DMPlexSetLabelValue(*dm, "zone", c, z);CHKERRQ(ierr);
216         }
217       } else {
218         switch (cellType) {
219         case TRI_3:   numCorners = 3;break;
220         case QUAD_4:  numCorners = 4;break;
221         case TETRA_4: numCorners = 4;break;
222         case HEXA_8:  numCorners = 8;break;
223         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
224         }
225 
226         /* CGNS uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */
227         for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
228           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) {
229             cone[v_loc] = elements[v]+numCells-1;
230           }
231           /* Tetrahedra are inverted */
232           if (cellType == TETRA_4) {
233             PetscInt tmp = cone[0];
234             cone[0] = cone[1];
235             cone[1] = tmp;
236           }
237           /* Hexahedra are inverted */
238           if (cellType == HEXA_8) {
239             PetscInt tmp = cone[1];
240             cone[1] = cone[3];
241             cone[3] = tmp;
242           }
243           ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
244           ierr = DMPlexSetLabelValue(*dm, "zone", c, z);CHKERRQ(ierr);
245         }
246       }
247       ierr = PetscFree2(elements,cone);CHKERRQ(ierr);
248     }
249   }
250   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
251   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
252   if (interpolate) {
253     DM idm = NULL;
254 
255     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
256     /* Maintain zone label */
257     {
258       DMLabel label;
259 
260       ierr = DMPlexRemoveLabel(*dm, "zone", &label);CHKERRQ(ierr);
261       if (label) {ierr = DMPlexAddLabel(idm, label);CHKERRQ(ierr);}
262     }
263     ierr = DMDestroy(dm);CHKERRQ(ierr);
264     *dm  = idm;
265   }
266 
267   /* Read coordinates */
268   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
269   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
270   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
271   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
272   for (v = numCells; v < numCells+numVertices; ++v) {
273     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
274     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
275   }
276   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
277   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
278   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
279   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
280   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
281   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
282   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
283   if (!rank) {
284     PetscInt off = 0;
285     float   *x[3];
286     int      z, d;
287 
288     ierr = PetscMalloc3(numVertices,&x[0],numVertices,&x[1],numVertices,&x[2]);CHKERRQ(ierr);
289     for (z = 1; z <= nzones; ++z) {
290       DataType_t datatype;
291       cgsize_t   sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
292       cgsize_t   range_min[3] = {1, 1, 1};
293       cgsize_t   range_max[3] = {1, 1, 1};
294       int        ngrids, ncoords;
295 
296 
297       ierr = cg_zone_read(cgid, 1, z, buffer, sizes);CHKERRQ(ierr);
298       range_max[0] = sizes[0];
299       ierr = cg_ngrids(cgid, 1, z, &ngrids);CHKERRQ(ierr);
300       if (ngrids > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single grid, not %d\n",ngrids);
301       ierr = cg_ncoords(cgid, 1, z, &ncoords);CHKERRQ(ierr);
302       if (ncoords != dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a coordinate array for each dimension, not %d\n",ncoords);
303       for (d = 0; d < dim; ++d) {
304         ierr = cg_coord_info(cgid, 1, z, 1+d, &datatype, buffer);CHKERRQ(ierr);
305         ierr = cg_coord_read(cgid, 1, z, buffer, RealSingle, range_min, range_max, x[d]);CHKERRQ(ierr);
306       }
307       if (dim > 0) {
308         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+0] = x[0][v];
309       }
310       if (dim > 1) {
311         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+1] = x[1][v];
312       }
313       if (dim > 2) {
314         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+2] = x[2][v];
315       }
316       off += sizes[0];
317     }
318     ierr = PetscFree3(x[0],x[1],x[2]);CHKERRQ(ierr);
319   }
320   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
321   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
322   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
323 #else
324   SETERRQ(comm, PETSC_ERR_SUP, "This method requires CGNS support. Reconfigure using --with-cgns-dir");
325 #endif
326   PetscFunctionReturn(0);
327 }
328