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