xref: /petsc/src/dm/impls/plex/plexcgns.c (revision 9566063d113dddea24716c546802770db7481bc0)
1e1b39ce3SMatthew G. Knepley #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3e1b39ce3SMatthew G. Knepley 
41b7ec9c1SMatthew G. Knepley #undef I /* Very old CGNS stupidly uses I as a variable, which fails when using complex. Curse you idiot package managers */
5e1b39ce3SMatthew G. Knepley #if defined(PETSC_HAVE_CGNS)
6e1b39ce3SMatthew G. Knepley #include <cgnslib.h>
7e1b39ce3SMatthew G. Knepley #include <cgns_io.h>
8e1b39ce3SMatthew G. Knepley #endif
9c4cbe8f1SToby Isaac #if !defined(CGNS_ENUMT)
10c4cbe8f1SToby Isaac #define CGNS_ENUMT(a) a
11c4cbe8f1SToby Isaac #endif
12c4cbe8f1SToby Isaac #if !defined(CGNS_ENUMV)
13c4cbe8f1SToby Isaac #define CGNS_ENUMV(a) a
14c4cbe8f1SToby Isaac #endif
15e1b39ce3SMatthew G. Knepley 
16*9566063dSJacob Faibussowitsch #define PetscCallCGNS(ierr) \
17c261946bSLisandro Dalcin do { \
18c261946bSLisandro Dalcin   int _cgns_ier = (ierr); \
1928b400f6SJacob Faibussowitsch   PetscCheck(!_cgns_ier,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS error %d %s",_cgns_ier,cg_get_error()); \
20c261946bSLisandro Dalcin } while (0)
21c261946bSLisandro Dalcin 
2244cd5272SMichael Lange /*@C
2344cd5272SMichael Lange   DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file.
2444cd5272SMichael Lange 
25d083f849SBarry Smith   Collective
2644cd5272SMichael Lange 
2744cd5272SMichael Lange   Input Parameters:
2844cd5272SMichael Lange + comm  - The MPI communicator
2944cd5272SMichael Lange . filename - The name of the CGNS file
3044cd5272SMichael Lange - interpolate - Create faces and edges in the mesh
3144cd5272SMichael Lange 
3244cd5272SMichael Lange   Output Parameter:
3344cd5272SMichael Lange . dm  - The DM object representing the mesh
3444cd5272SMichael Lange 
3544cd5272SMichael Lange   Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html
3644cd5272SMichael Lange 
3744cd5272SMichael Lange   Level: beginner
3844cd5272SMichael Lange 
3944cd5272SMichael Lange .seealso: DMPlexCreate(), DMPlexCreateCGNS(), DMPlexCreateExodus()
4044cd5272SMichael Lange @*/
4144cd5272SMichael Lange PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
4244cd5272SMichael Lange {
4344cd5272SMichael Lange   PetscMPIInt    rank;
4444cd5272SMichael Lange #if defined(PETSC_HAVE_CGNS)
4544cd5272SMichael Lange   int cgid = -1;
4644cd5272SMichael Lange #endif
4744cd5272SMichael Lange 
4844cd5272SMichael Lange   PetscFunctionBegin;
4944cd5272SMichael Lange   PetscValidCharPointer(filename, 2);
50*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5144cd5272SMichael Lange #if defined(PETSC_HAVE_CGNS)
52dd400576SPatrick Sanan   if (rank == 0) {
53*9566063dSJacob Faibussowitsch     PetscCallCGNS(cg_open(filename, CG_MODE_READ, &cgid));
542c71b3e2SJacob Faibussowitsch     PetscCheckFalse(cgid <= 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename);
5544cd5272SMichael Lange   }
56*9566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm));
57*9566063dSJacob Faibussowitsch   if (rank == 0) PetscCallCGNS(cg_close(cgid));
58b458e8f1SJose E. Roman   PetscFunctionReturn(0);
5944cd5272SMichael Lange #else
6044cd5272SMichael Lange   SETERRQ(comm, PETSC_ERR_SUP, "Loading meshes requires CGNS support. Reconfigure using --with-cgns-dir");
6144cd5272SMichael Lange #endif
6244cd5272SMichael Lange }
6344cd5272SMichael Lange 
64e1b39ce3SMatthew G. Knepley /*@
6544cd5272SMichael Lange   DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file ID.
66e1b39ce3SMatthew G. Knepley 
67d083f849SBarry Smith   Collective
68e1b39ce3SMatthew G. Knepley 
69e1b39ce3SMatthew G. Knepley   Input Parameters:
70e1b39ce3SMatthew G. Knepley + comm  - The MPI communicator
71e1b39ce3SMatthew G. Knepley . cgid - The CG id associated with a file and obtained using cg_open
72e1b39ce3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
73e1b39ce3SMatthew G. Knepley 
74e1b39ce3SMatthew G. Knepley   Output Parameter:
75e1b39ce3SMatthew G. Knepley . dm  - The DM object representing the mesh
76e1b39ce3SMatthew G. Knepley 
77e1b39ce3SMatthew G. Knepley   Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html
78e1b39ce3SMatthew G. Knepley 
79e1b39ce3SMatthew G. Knepley   Level: beginner
80e1b39ce3SMatthew G. Knepley 
81e1b39ce3SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexCreateExodus()
82e1b39ce3SMatthew G. Knepley @*/
83e1b39ce3SMatthew G. Knepley PetscErrorCode DMPlexCreateCGNS(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm)
84e1b39ce3SMatthew G. Knepley {
85e1b39ce3SMatthew G. Knepley #if defined(PETSC_HAVE_CGNS)
86e1b39ce3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
87c261946bSLisandro Dalcin   DM             cdm;
884e6cc61cSLisandro Dalcin   DMLabel        label;
89e1b39ce3SMatthew G. Knepley   PetscSection   coordSection;
90e1b39ce3SMatthew G. Knepley   Vec            coordinates;
91e1b39ce3SMatthew G. Knepley   PetscScalar   *coords;
92c261946bSLisandro Dalcin   PetscInt      *cellStart, *vertStart, v;
934e6cc61cSLisandro Dalcin   PetscInt       labelIdRange[2], labelId;
94e1b39ce3SMatthew G. Knepley   /* Read from file */
95e1b39ce3SMatthew G. Knepley   char basename[CGIO_MAX_NAME_LENGTH+1];
96e1b39ce3SMatthew G. Knepley   char buffer[CGIO_MAX_NAME_LENGTH+1];
97c261946bSLisandro Dalcin   int  dim    = 0, physDim = 0, coordDim =0, numVertices = 0, numCells = 0;
98e1b39ce3SMatthew G. Knepley   int  nzones = 0;
99e1b39ce3SMatthew G. Knepley #endif
100e1b39ce3SMatthew G. Knepley 
101e1b39ce3SMatthew G. Knepley   PetscFunctionBegin;
102e1b39ce3SMatthew G. Knepley #if defined(PETSC_HAVE_CGNS)
103*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
104*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &num_proc));
105*9566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
106*9566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
107c261946bSLisandro Dalcin 
108a5b23f4aSJose E. Roman   /* Open CGNS II file and read basic information on rank 0, then broadcast to all processors */
109dd400576SPatrick Sanan   if (rank == 0) {
110e1b39ce3SMatthew G. Knepley     int nbases, z;
111e1b39ce3SMatthew G. Knepley 
112*9566063dSJacob Faibussowitsch     PetscCallCGNS(cg_nbases(cgid, &nbases));
1132c71b3e2SJacob Faibussowitsch     PetscCheckFalse(nbases > 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single base, not %d",nbases);
114*9566063dSJacob Faibussowitsch     PetscCallCGNS(cg_base_read(cgid, 1, basename, &dim, &physDim));
115*9566063dSJacob Faibussowitsch     PetscCallCGNS(cg_nzones(cgid, 1, &nzones));
116*9566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(nzones+1, &cellStart, nzones+1, &vertStart));
117e1b39ce3SMatthew G. Knepley     for (z = 1; z <= nzones; ++z) {
118e1b39ce3SMatthew G. Knepley       cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
119e1b39ce3SMatthew G. Knepley 
120*9566063dSJacob Faibussowitsch       PetscCallCGNS(cg_zone_read(cgid, 1, z, buffer, sizes));
121e1b39ce3SMatthew G. Knepley       numVertices += sizes[0];
122e1b39ce3SMatthew G. Knepley       numCells    += sizes[1];
123f8716870SMatthew G. Knepley       cellStart[z] += sizes[1] + cellStart[z-1];
124f8716870SMatthew G. Knepley       vertStart[z] += sizes[0] + vertStart[z-1];
125f8716870SMatthew G. Knepley     }
126f8716870SMatthew G. Knepley     for (z = 1; z <= nzones; ++z) {
127f8716870SMatthew G. Knepley       vertStart[z] += numCells;
128e1b39ce3SMatthew G. Knepley     }
129c261946bSLisandro Dalcin     coordDim = dim;
130e1b39ce3SMatthew G. Knepley   }
131*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH+1, MPI_CHAR, 0, comm));
132*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm));
133*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm));
134*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&nzones, 1, MPI_INT, 0, comm));
135c261946bSLisandro Dalcin 
136*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *dm, basename));
137*9566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*dm, dim));
138*9566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dm, "celltype"));
139*9566063dSJacob Faibussowitsch   PetscCall(DMPlexSetChart(*dm, 0, numCells+numVertices));
140e1b39ce3SMatthew G. Knepley 
141e1b39ce3SMatthew G. Knepley   /* Read zone information */
142dd400576SPatrick Sanan   if (rank == 0) {
143c261946bSLisandro Dalcin     int z, c, c_loc;
144e1b39ce3SMatthew G. Knepley 
145e1b39ce3SMatthew G. Knepley     /* Read the cell set connectivity table and build mesh topology
146e1b39ce3SMatthew G. Knepley        CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */
147e1b39ce3SMatthew G. Knepley     /* First set sizes */
148e1b39ce3SMatthew G. Knepley     for (z = 1, c = 0; z <= nzones; ++z) {
149c4cbe8f1SToby Isaac       CGNS_ENUMT(ZoneType_t)    zonetype;
150e1b39ce3SMatthew G. Knepley       int                       nsections;
151c4cbe8f1SToby Isaac       CGNS_ENUMT(ElementType_t) cellType;
152e1b39ce3SMatthew G. Knepley       cgsize_t                  start, end;
153e1b39ce3SMatthew G. Knepley       int                       nbndry, parentFlag;
154e1b39ce3SMatthew G. Knepley       PetscInt                  numCorners;
155c261946bSLisandro Dalcin       DMPolytopeType            ctype;
156e1b39ce3SMatthew G. Knepley 
157*9566063dSJacob Faibussowitsch       PetscCallCGNS(cg_zone_type(cgid, 1, z, &zonetype));
1582c71b3e2SJacob Faibussowitsch       PetscCheckFalse(zonetype == CGNS_ENUMV(Structured),PETSC_COMM_SELF,PETSC_ERR_LIB,"Can only handle Unstructured zones for CGNS");
159*9566063dSJacob Faibussowitsch       PetscCallCGNS(cg_nsections(cgid, 1, z, &nsections));
1602c71b3e2SJacob Faibussowitsch       PetscCheckFalse(nsections > 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single section, not %d",nsections);
161*9566063dSJacob Faibussowitsch       PetscCallCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag));
162e1b39ce3SMatthew G. Knepley       /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
163c4cbe8f1SToby Isaac       if (cellType == CGNS_ENUMV(MIXED)) {
164e1b39ce3SMatthew G. Knepley         cgsize_t elementDataSize, *elements;
165e1b39ce3SMatthew G. Knepley         PetscInt off;
166e1b39ce3SMatthew G. Knepley 
167*9566063dSJacob Faibussowitsch         PetscCallCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize));
168*9566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(elementDataSize, &elements));
169*9566063dSJacob Faibussowitsch         PetscCallCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL));
17059dfc837SStefan Lemvig Glimberg         for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) {
1714c9a562dSMatthew G. Knepley           switch (elements[off]) {
172c261946bSLisandro Dalcin           case CGNS_ENUMV(BAR_2):   numCorners = 2; ctype = DM_POLYTOPE_SEGMENT;       break;
173c261946bSLisandro Dalcin           case CGNS_ENUMV(TRI_3):   numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE;      break;
174c261946bSLisandro Dalcin           case CGNS_ENUMV(QUAD_4):  numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break;
175c261946bSLisandro Dalcin           case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON;   break;
176c261946bSLisandro Dalcin           case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM;     break;
177c261946bSLisandro Dalcin           case CGNS_ENUMV(HEXA_8):  numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON;    break;
17898921bdaSJacob Faibussowitsch           default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[off]);
179e1b39ce3SMatthew G. Knepley           }
180*9566063dSJacob Faibussowitsch           PetscCall(DMPlexSetConeSize(*dm, c, numCorners));
181*9566063dSJacob Faibussowitsch           PetscCall(DMPlexSetCellType(*dm, c, ctype));
182e1b39ce3SMatthew G. Knepley           off += numCorners+1;
183e1b39ce3SMatthew G. Knepley         }
184*9566063dSJacob Faibussowitsch         PetscCall(PetscFree(elements));
185e1b39ce3SMatthew G. Knepley       } else {
186e1b39ce3SMatthew G. Knepley         switch (cellType) {
187c261946bSLisandro Dalcin         case CGNS_ENUMV(BAR_2):   numCorners = 2; ctype = DM_POLYTOPE_SEGMENT;       break;
188c261946bSLisandro Dalcin         case CGNS_ENUMV(TRI_3):   numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE;      break;
189c261946bSLisandro Dalcin         case CGNS_ENUMV(QUAD_4):  numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break;
190c261946bSLisandro Dalcin         case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON;   break;
191c261946bSLisandro Dalcin         case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM;     break;
192c261946bSLisandro Dalcin         case CGNS_ENUMV(HEXA_8):  numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON;    break;
19398921bdaSJacob Faibussowitsch         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
194e1b39ce3SMatthew G. Knepley         }
19559dfc837SStefan Lemvig Glimberg         for (c_loc = start; c_loc <= end; ++c_loc, ++c) {
196*9566063dSJacob Faibussowitsch           PetscCall(DMPlexSetConeSize(*dm, c, numCorners));
197*9566063dSJacob Faibussowitsch           PetscCall(DMPlexSetCellType(*dm, c, ctype));
198e1b39ce3SMatthew G. Knepley         }
199e1b39ce3SMatthew G. Knepley       }
200e1b39ce3SMatthew G. Knepley     }
201c261946bSLisandro Dalcin     for (v = numCells; v < numCells+numVertices; ++v) {
202*9566063dSJacob Faibussowitsch       PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
203c261946bSLisandro Dalcin     }
204c261946bSLisandro Dalcin   }
205c261946bSLisandro Dalcin 
206*9566063dSJacob Faibussowitsch   PetscCall(DMSetUp(*dm));
207c261946bSLisandro Dalcin 
208*9566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dm, "zone"));
209dd400576SPatrick Sanan   if (rank == 0) {
210c261946bSLisandro Dalcin     int z, c, c_loc, v_loc;
211c261946bSLisandro Dalcin 
212*9566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "zone", &label));
213e1b39ce3SMatthew G. Knepley     for (z = 1, c = 0; z <= nzones; ++z) {
214c4cbe8f1SToby Isaac       CGNS_ENUMT(ElementType_t)   cellType;
215c261946bSLisandro Dalcin       cgsize_t                    elementDataSize, *elements, start, end;
216e1b39ce3SMatthew G. Knepley       int                          nbndry, parentFlag;
217e1b39ce3SMatthew G. Knepley       PetscInt                    *cone, numc, numCorners, maxCorners = 27;
218e1b39ce3SMatthew G. Knepley 
219*9566063dSJacob Faibussowitsch       PetscCallCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag));
220e1b39ce3SMatthew G. Knepley       numc = end - start;
221e1b39ce3SMatthew G. Knepley       /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
222*9566063dSJacob Faibussowitsch       PetscCallCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize));
223*9566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(elementDataSize,&elements,maxCorners,&cone));
224*9566063dSJacob Faibussowitsch       PetscCallCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL));
225c4cbe8f1SToby Isaac       if (cellType == CGNS_ENUMV(MIXED)) {
226eaf898f9SPatrick Sanan         /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
22759dfc837SStefan Lemvig Glimberg         for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
2284c9a562dSMatthew G. Knepley           switch (elements[v]) {
229c261946bSLisandro Dalcin           case CGNS_ENUMV(BAR_2):   numCorners = 2; break;
230c4cbe8f1SToby Isaac           case CGNS_ENUMV(TRI_3):   numCorners = 3; break;
231c4cbe8f1SToby Isaac           case CGNS_ENUMV(QUAD_4):  numCorners = 4; break;
232c4cbe8f1SToby Isaac           case CGNS_ENUMV(TETRA_4): numCorners = 4; break;
233c261946bSLisandro Dalcin           case CGNS_ENUMV(PENTA_6): numCorners = 6; break;
234c4cbe8f1SToby Isaac           case CGNS_ENUMV(HEXA_8):  numCorners = 8; break;
23598921bdaSJacob Faibussowitsch           default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[v]);
236e1b39ce3SMatthew G. Knepley           }
2374c9a562dSMatthew G. Knepley           ++v;
238e1b39ce3SMatthew G. Knepley           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) {
239e1b39ce3SMatthew G. Knepley             cone[v_loc] = elements[v]+numCells-1;
240e1b39ce3SMatthew G. Knepley           }
241*9566063dSJacob Faibussowitsch           PetscCall(DMPlexReorderCell(*dm, c, cone));
242*9566063dSJacob Faibussowitsch           PetscCall(DMPlexSetCone(*dm, c, cone));
243*9566063dSJacob Faibussowitsch           PetscCall(DMLabelSetValue(label, c, z));
244e1b39ce3SMatthew G. Knepley         }
245e1b39ce3SMatthew G. Knepley       } else {
246e1b39ce3SMatthew G. Knepley         switch (cellType) {
247c261946bSLisandro Dalcin         case CGNS_ENUMV(BAR_2):   numCorners = 2; break;
248c4cbe8f1SToby Isaac         case CGNS_ENUMV(TRI_3):   numCorners = 3; break;
249c4cbe8f1SToby Isaac         case CGNS_ENUMV(QUAD_4):  numCorners = 4; break;
250c4cbe8f1SToby Isaac         case CGNS_ENUMV(TETRA_4): numCorners = 4; break;
251c261946bSLisandro Dalcin         case CGNS_ENUMV(PENTA_6): numCorners = 6; break;
252c4cbe8f1SToby Isaac         case CGNS_ENUMV(HEXA_8):  numCorners = 8; break;
25398921bdaSJacob Faibussowitsch         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
254e1b39ce3SMatthew G. Knepley         }
255eaf898f9SPatrick Sanan         /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
25659dfc837SStefan Lemvig Glimberg         for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
257e1b39ce3SMatthew G. Knepley           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) {
258e1b39ce3SMatthew G. Knepley             cone[v_loc] = elements[v]+numCells-1;
259e1b39ce3SMatthew G. Knepley           }
260*9566063dSJacob Faibussowitsch           PetscCall(DMPlexReorderCell(*dm, c, cone));
261*9566063dSJacob Faibussowitsch           PetscCall(DMPlexSetCone(*dm, c, cone));
262*9566063dSJacob Faibussowitsch           PetscCall(DMLabelSetValue(label, c, z));
263e1b39ce3SMatthew G. Knepley         }
264e1b39ce3SMatthew G. Knepley       }
265*9566063dSJacob Faibussowitsch       PetscCall(PetscFree2(elements,cone));
266e1b39ce3SMatthew G. Knepley     }
267e1b39ce3SMatthew G. Knepley   }
268c261946bSLisandro Dalcin 
269*9566063dSJacob Faibussowitsch   PetscCall(DMPlexSymmetrize(*dm));
270*9566063dSJacob Faibussowitsch   PetscCall(DMPlexStratify(*dm));
271e1b39ce3SMatthew G. Knepley   if (interpolate) {
2725fd9971aSMatthew G. Knepley     DM idm;
273e1b39ce3SMatthew G. Knepley 
274*9566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
275*9566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
276e1b39ce3SMatthew G. Knepley     *dm  = idm;
277e1b39ce3SMatthew G. Knepley   }
278e1b39ce3SMatthew G. Knepley 
279e1b39ce3SMatthew G. Knepley   /* Read coordinates */
280*9566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateDim(*dm, coordDim));
281*9566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(*dm, &cdm));
282*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(cdm, &coordSection));
283*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(coordSection, 1));
284*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, coordDim));
285*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
286e1b39ce3SMatthew G. Knepley   for (v = numCells; v < numCells+numVertices; ++v) {
287*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(coordSection, v, dim));
288*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, coordDim));
289e1b39ce3SMatthew G. Knepley   }
290*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(coordSection));
291c261946bSLisandro Dalcin 
292*9566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(cdm, &coordinates));
293*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &coords));
294dd400576SPatrick Sanan   if (rank == 0) {
295e1b39ce3SMatthew G. Knepley     PetscInt off = 0;
296e1b39ce3SMatthew G. Knepley     float   *x[3];
29759dfc837SStefan Lemvig Glimberg     int      z, d;
298e1b39ce3SMatthew G. Knepley 
299*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(numVertices,&x[0],numVertices,&x[1],numVertices,&x[2]));
30059dfc837SStefan Lemvig Glimberg     for (z = 1; z <= nzones; ++z) {
301c4cbe8f1SToby Isaac       CGNS_ENUMT(DataType_t) datatype;
302e1b39ce3SMatthew G. Knepley       cgsize_t               sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
303e1b39ce3SMatthew G. Knepley       cgsize_t               range_min[3] = {1, 1, 1};
304e1b39ce3SMatthew G. Knepley       cgsize_t               range_max[3] = {1, 1, 1};
305e1b39ce3SMatthew G. Knepley       int                    ngrids, ncoords;
306e1b39ce3SMatthew G. Knepley 
307*9566063dSJacob Faibussowitsch       PetscCallCGNS(cg_zone_read(cgid, 1, z, buffer, sizes));
308e1b39ce3SMatthew G. Knepley       range_max[0] = sizes[0];
309*9566063dSJacob Faibussowitsch       PetscCallCGNS(cg_ngrids(cgid, 1, z, &ngrids));
3102c71b3e2SJacob Faibussowitsch       PetscCheckFalse(ngrids > 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single grid, not %d",ngrids);
311*9566063dSJacob Faibussowitsch       PetscCallCGNS(cg_ncoords(cgid, 1, z, &ncoords));
3122c71b3e2SJacob Faibussowitsch       PetscCheckFalse(ncoords != coordDim,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a coordinate array for each dimension, not %d",ncoords);
313c261946bSLisandro Dalcin       for (d = 0; d < coordDim; ++d) {
314*9566063dSJacob Faibussowitsch         PetscCallCGNS(cg_coord_info(cgid, 1, z, 1+d, &datatype, buffer));
315*9566063dSJacob Faibussowitsch         PetscCallCGNS(cg_coord_read(cgid, 1, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d]));
316e1b39ce3SMatthew G. Knepley       }
317c261946bSLisandro Dalcin       if (coordDim >= 1) {
318c261946bSLisandro Dalcin         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+0] = x[0][v];
319e1b39ce3SMatthew G. Knepley       }
320c261946bSLisandro Dalcin       if (coordDim >= 2) {
321c261946bSLisandro Dalcin         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+1] = x[1][v];
322e1b39ce3SMatthew G. Knepley       }
323c261946bSLisandro Dalcin       if (coordDim >= 3) {
324c261946bSLisandro Dalcin         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+2] = x[2][v];
325e1b39ce3SMatthew G. Knepley       }
326e1b39ce3SMatthew G. Knepley       off += sizes[0];
327e1b39ce3SMatthew G. Knepley     }
328*9566063dSJacob Faibussowitsch     PetscCall(PetscFree3(x[0],x[1],x[2]));
329e1b39ce3SMatthew G. Knepley   }
330*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &coords));
331c261946bSLisandro Dalcin 
332*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates"));
333*9566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(coordinates, coordDim));
334*9566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
335*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coordinates));
336c261946bSLisandro Dalcin 
337f8716870SMatthew G. Knepley   /* Read boundary conditions */
338*9566063dSJacob Faibussowitsch   PetscCall(DMGetNumLabels(*dm, &labelIdRange[0]));
339dd400576SPatrick Sanan   if (rank == 0) {
340c4cbe8f1SToby Isaac     CGNS_ENUMT(BCType_t)        bctype;
341c4cbe8f1SToby Isaac     CGNS_ENUMT(DataType_t)      datatype;
342c4cbe8f1SToby Isaac     CGNS_ENUMT(PointSetType_t)  pointtype;
343f8716870SMatthew G. Knepley     cgsize_t                    *points;
344f8716870SMatthew G. Knepley     PetscReal                   *normals;
345f8716870SMatthew G. Knepley     int                         normal[3];
346f8716870SMatthew G. Knepley     char                        *bcname = buffer;
347f8716870SMatthew G. Knepley     cgsize_t                    npoints, nnormals;
348f8716870SMatthew G. Knepley     int                         z, nbc, bc, c, ndatasets;
349f8716870SMatthew G. Knepley 
350f8716870SMatthew G. Knepley     for (z = 1; z <= nzones; ++z) {
351*9566063dSJacob Faibussowitsch       PetscCallCGNS(cg_nbocos(cgid, 1, z, &nbc));
352f8716870SMatthew G. Knepley       for (bc = 1; bc <= nbc; ++bc) {
353*9566063dSJacob Faibussowitsch         PetscCallCGNS(cg_boco_info(cgid, 1, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets));
354*9566063dSJacob Faibussowitsch         PetscCall(DMCreateLabel(*dm, bcname));
355*9566063dSJacob Faibussowitsch         PetscCall(DMGetLabel(*dm, bcname, &label));
356*9566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(npoints, &points, nnormals, &normals));
357*9566063dSJacob Faibussowitsch         PetscCallCGNS(cg_boco_read(cgid, 1, z, bc, points, (void *) normals));
358c4cbe8f1SToby Isaac         if (pointtype == CGNS_ENUMV(ElementRange)) {
359f8716870SMatthew G. Knepley           /* Range of cells: assuming half-open interval since the documentation sucks */
360f8716870SMatthew G. Knepley           for (c = points[0]; c < points[1]; ++c) {
361*9566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(label, c - cellStart[z-1], 1));
362f8716870SMatthew G. Knepley           }
363c4cbe8f1SToby Isaac         } else if (pointtype == CGNS_ENUMV(ElementList)) {
364f8716870SMatthew G. Knepley           /* List of cells */
365f8716870SMatthew G. Knepley           for (c = 0; c < npoints; ++c) {
366*9566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(label, points[c] - cellStart[z-1], 1));
367f8716870SMatthew G. Knepley           }
368c4cbe8f1SToby Isaac         } else if (pointtype == CGNS_ENUMV(PointRange)) {
369c4cbe8f1SToby Isaac           CGNS_ENUMT(GridLocation_t) gridloc;
370f8716870SMatthew G. Knepley 
371f8716870SMatthew G. Knepley           /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
372*9566063dSJacob Faibussowitsch           PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end"));
373*9566063dSJacob Faibussowitsch           PetscCallCGNS(cg_gridlocation_read(&gridloc));
374f8716870SMatthew G. Knepley           /* Range of points: assuming half-open interval since the documentation sucks */
375f8716870SMatthew G. Knepley           for (c = points[0]; c < points[1]; ++c) {
376*9566063dSJacob Faibussowitsch             if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, c - vertStart[z-1], 1));
377*9566063dSJacob Faibussowitsch             else                               PetscCall(DMLabelSetValue(label, c - cellStart[z-1], 1));
378f8716870SMatthew G. Knepley           }
379c4cbe8f1SToby Isaac         } else if (pointtype == CGNS_ENUMV(PointList)) {
380c4cbe8f1SToby Isaac           CGNS_ENUMT(GridLocation_t) gridloc;
381f8716870SMatthew G. Knepley 
382f8716870SMatthew G. Knepley           /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
383*9566063dSJacob Faibussowitsch           PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end"));
384*9566063dSJacob Faibussowitsch           PetscCallCGNS(cg_gridlocation_read(&gridloc));
385f8716870SMatthew G. Knepley           for (c = 0; c < npoints; ++c) {
386*9566063dSJacob Faibussowitsch             if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, points[c] - vertStart[z-1], 1));
387*9566063dSJacob Faibussowitsch             else                               PetscCall(DMLabelSetValue(label, points[c] - cellStart[z-1], 1));
388f8716870SMatthew G. Knepley           }
38998921bdaSJacob Faibussowitsch         } else SETERRQ(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int) pointtype);
390*9566063dSJacob Faibussowitsch         PetscCall(PetscFree2(points, normals));
391f8716870SMatthew G. Knepley       }
392f8716870SMatthew G. Knepley     }
393*9566063dSJacob Faibussowitsch     PetscCall(PetscFree2(cellStart, vertStart));
394f8716870SMatthew G. Knepley   }
395*9566063dSJacob Faibussowitsch   PetscCall(DMGetNumLabels(*dm, &labelIdRange[1]));
396*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(labelIdRange, 2, MPIU_INT, 0, comm));
3974e6cc61cSLisandro Dalcin 
3984e6cc61cSLisandro Dalcin   /* Create BC labels at all processes */
3994e6cc61cSLisandro Dalcin   for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) {
4004e6cc61cSLisandro Dalcin     char *labelName = buffer;
4014e6cc61cSLisandro Dalcin     size_t len = sizeof(buffer);
4024e6cc61cSLisandro Dalcin     const char *locName;
4034e6cc61cSLisandro Dalcin 
404dd400576SPatrick Sanan     if (rank == 0) {
405*9566063dSJacob Faibussowitsch       PetscCall(DMGetLabelByNum(*dm, labelId, &label));
406*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)label, &locName));
407*9566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(labelName, locName, len));
4084e6cc61cSLisandro Dalcin     }
409*9566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm));
410*9566063dSJacob Faibussowitsch     PetscCallMPI(DMCreateLabel(*dm, labelName));
4114e6cc61cSLisandro Dalcin   }
412b458e8f1SJose E. Roman   PetscFunctionReturn(0);
413e1b39ce3SMatthew G. Knepley #else
414e1b39ce3SMatthew G. Knepley   SETERRQ(comm, PETSC_ERR_SUP, "This method requires CGNS support. Reconfigure using --with-cgns-dir");
415e1b39ce3SMatthew G. Knepley #endif
416e1b39ce3SMatthew G. Knepley }
417