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