1 #define PETSCDM_DLL 2 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3 4 #if defined(PETSC_HAVE_EXODUSII) 5 #include <netcdf.h> 6 #include <exodusII.h> 7 #endif 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "DMPlexCreateExodus" 11 /*@ 12 DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file. 13 14 Collective on comm 15 16 Input Parameters: 17 + comm - The MPI communicator 18 . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 19 - interpolate - Create faces and edges in the mesh 20 21 Output Parameter: 22 . dm - The DM object representing the mesh 23 24 Level: beginner 25 26 .keywords: mesh,ExodusII 27 .seealso: MeshCreate(), MeshCreateExodus() 28 @*/ 29 PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 30 { 31 #if defined(PETSC_HAVE_EXODUSII) 32 PetscMPIInt num_proc, rank; 33 PetscSection coordSection; 34 Vec coordinates; 35 PetscScalar *coords; 36 PetscInt coordSize, v; 37 PetscErrorCode ierr; 38 /* Read from ex_get_init() */ 39 char title[PETSC_MAX_PATH_LEN+1]; 40 int dim = 0, numVertices = 0, numCells = 0; 41 int num_cs = 0, num_vs = 0, num_fs = 0; 42 #endif 43 44 PetscFunctionBegin; 45 #if defined(PETSC_HAVE_EXODUSII) 46 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 47 ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 48 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 49 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 50 /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */ 51 if (!rank) { 52 ierr = ex_get_init(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs); 53 if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ExodusII ex_get_init() failed with error code %D\n",ierr); 54 if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n"); 55 } 56 ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 57 ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 58 ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr); 59 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 60 ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 61 62 /* Read cell sets information */ 63 if (!rank) { 64 PetscInt *cone; 65 int c, cs, c_loc, v, v_loc; 66 /* Read from ex_get_elem_blk_ids() */ 67 int *cs_id; 68 /* Read from ex_get_elem_block() */ 69 char buffer[PETSC_MAX_PATH_LEN+1]; 70 int num_cell_in_set, num_vertex_per_cell, num_attr; 71 /* Read from ex_get_elem_conn() */ 72 int *cs_connect; 73 74 /* Get cell sets IDs */ 75 ierr = PetscMalloc(num_cs * sizeof(int), &cs_id);CHKERRQ(ierr); 76 ierr = ex_get_elem_blk_ids(exoid, cs_id);CHKERRQ(ierr); 77 /* Read the cell set connectivity table and build mesh topology 78 EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 79 /* First set sizes */ 80 for (cs = 0, c = 0; cs < num_cs; cs++) { 81 ierr = ex_get_elem_block(exoid, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, &num_attr);CHKERRQ(ierr); 82 for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 83 ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr); 84 } 85 } 86 ierr = DMSetUp(*dm);CHKERRQ(ierr); 87 for (cs = 0, c = 0; cs < num_cs; cs++) { 88 ierr = ex_get_elem_block(exoid, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, &num_attr);CHKERRQ(ierr); 89 ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,int,&cs_connect,num_vertex_per_cell,PetscInt,&cone);CHKERRQ(ierr); 90 ierr = ex_get_elem_conn(exoid, cs_id[cs], cs_connect);CHKERRQ(ierr); 91 /* EXO uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */ 92 for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 93 for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) { 94 cone[v_loc] = cs_connect[v]+numCells-1; 95 } 96 /* Tetrahedra are inverted */ 97 if (num_vertex_per_cell == 4) { 98 PetscInt tmp = cone[0]; 99 cone[0] = cone[1]; 100 cone[1] = tmp; 101 } 102 /* Hexahedra are inverted */ 103 if (num_vertex_per_cell == 8) { 104 PetscInt tmp = cone[1]; 105 cone[1] = cone[3]; 106 cone[3] = tmp; 107 } 108 ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 109 ierr = DMPlexSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr); 110 } 111 ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr); 112 } 113 ierr = PetscFree(cs_id);CHKERRQ(ierr); 114 } 115 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 116 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 117 if (interpolate) { 118 DM idm; 119 120 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 121 /* Maintain Cell Sets label */ 122 { 123 DMLabel label; 124 125 ierr = DMPlexRemoveLabel(*dm, "Cell Sets", &label);CHKERRQ(ierr); 126 if (label) {ierr = DMPlexAddLabel(idm, label);CHKERRQ(ierr);} 127 } 128 ierr = DMDestroy(dm);CHKERRQ(ierr); 129 *dm = idm; 130 } 131 132 /* Create vertex set label */ 133 if (!rank && (num_vs > 0)) { 134 int vs, v; 135 /* Read from ex_get_node_set_ids() */ 136 int *vs_id; 137 /* Read from ex_get_node_set_param() */ 138 int num_vertex_in_set, num_attr; 139 /* Read from ex_get_node_set() */ 140 int *vs_vertex_list; 141 142 /* Get vertex set ids */ 143 ierr = PetscMalloc(num_vs * sizeof(int), &vs_id);CHKERRQ(ierr); 144 ierr = ex_get_node_set_ids(exoid, vs_id);CHKERRQ(ierr); 145 for (vs = 0; vs < num_vs; ++vs) { 146 ierr = ex_get_node_set_param(exoid, vs_id[vs], &num_vertex_in_set, &num_attr);CHKERRQ(ierr); 147 ierr = PetscMalloc(num_vertex_in_set * sizeof(int), &vs_vertex_list);CHKERRQ(ierr); 148 ierr = ex_get_node_set(exoid, vs_id[vs], vs_vertex_list);CHKERRQ(ierr); 149 for (v = 0; v < num_vertex_in_set; ++v) { 150 ierr = DMPlexSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr); 151 } 152 ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr); 153 } 154 ierr = PetscFree(vs_id);CHKERRQ(ierr); 155 } 156 /* Read coordinates */ 157 ierr = DMPlexGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 158 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 159 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 160 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 161 for (v = numCells; v < numCells+numVertices; ++v) { 162 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 163 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 164 } 165 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 166 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 167 ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 168 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 169 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 170 ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 171 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 172 if (!rank) { 173 float *x, *y, *z; 174 175 ierr = PetscMalloc3(numVertices,float,&x,numVertices,float,&y,numVertices,float,&z);CHKERRQ(ierr); 176 ierr = ex_get_coord(exoid, x, y, z);CHKERRQ(ierr); 177 if (dim > 0) { 178 for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v]; 179 } 180 if (dim > 1) { 181 for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v]; 182 } 183 if (dim > 2) { 184 for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v]; 185 } 186 ierr = PetscFree3(x,y,z);CHKERRQ(ierr); 187 } 188 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 189 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 190 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 191 192 /* Create side set label */ 193 if (!rank && interpolate && (num_fs > 0)) { 194 int fs, f, voff; 195 /* Read from ex_get_side_set_ids() */ 196 int *fs_id; 197 /* Read from ex_get_side_set_param() */ 198 int num_side_in_set, num_dist_fact_in_set; 199 /* Read from ex_get_side_set_node_list() */ 200 int *fs_vertex_count_list, *fs_vertex_list; 201 202 /* Get side set ids */ 203 ierr = PetscMalloc(num_fs * sizeof(int), &fs_id);CHKERRQ(ierr); 204 ierr = ex_get_side_set_ids(exoid, fs_id);CHKERRQ(ierr); 205 for (fs = 0; fs < num_fs; ++fs) { 206 ierr = ex_get_side_set_param(exoid, fs_id[fs], &num_side_in_set, &num_dist_fact_in_set);CHKERRQ(ierr); 207 ierr = PetscMalloc2(num_side_in_set,int,&fs_vertex_count_list,num_side_in_set*4,int,&fs_vertex_list);CHKERRQ(ierr); 208 ierr = ex_get_side_set_node_list(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);CHKERRQ(ierr); 209 for (f = 0, voff = 0; f < num_side_in_set; ++f) { 210 const PetscInt *faces = NULL; 211 PetscInt faceSize = fs_vertex_count_list[f], numFaces; 212 PetscInt faceVertices[4], v; 213 214 if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize); 215 for (v = 0; v < faceSize; ++v, ++voff) { 216 faceVertices[v] = fs_vertex_list[voff]+numCells-1; 217 } 218 ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 219 if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces); 220 ierr = DMPlexSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr); 221 ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 222 } 223 ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr); 224 } 225 ierr = PetscFree(fs_id);CHKERRQ(ierr); 226 } 227 #else 228 SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 229 #endif 230 PetscFunctionReturn(0); 231 } 232