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