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