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