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