1552f7358SJed Brown #define PETSCDM_DLL 234541f0dSBarry Smith #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3552f7358SJed Brown 4519f805aSKarl Rupp #if defined(PETSC_HAVE_EXODUSII) 5552f7358SJed Brown #include <netcdf.h> 6552f7358SJed Brown #include <exodusII.h> 7552f7358SJed Brown #endif 8552f7358SJed Brown 9552f7358SJed Brown #undef __FUNCT__ 10552f7358SJed Brown #define __FUNCT__ "DMPlexCreateExodus" 11552f7358SJed Brown /*@ 12552f7358SJed Brown DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file. 13552f7358SJed Brown 14552f7358SJed Brown Collective on comm 15552f7358SJed Brown 16552f7358SJed Brown Input Parameters: 17552f7358SJed Brown + comm - The MPI communicator 18552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 19552f7358SJed Brown - interpolate - Create faces and edges in the mesh 20552f7358SJed Brown 21552f7358SJed Brown Output Parameter: 22552f7358SJed Brown . dm - The DM object representing the mesh 23552f7358SJed Brown 24552f7358SJed Brown Level: beginner 25552f7358SJed Brown 26552f7358SJed Brown .keywords: mesh,ExodusII 27e84e7769SJed Brown .seealso: DMPLEX, DMCreate() 28552f7358SJed Brown @*/ 29552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 30552f7358SJed Brown { 31552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 32552f7358SJed Brown PetscMPIInt num_proc, rank; 33552f7358SJed Brown PetscSection coordSection; 34552f7358SJed Brown Vec coordinates; 35552f7358SJed Brown PetscScalar *coords; 36552f7358SJed Brown PetscInt coordSize, v; 37552f7358SJed Brown PetscErrorCode ierr; 38552f7358SJed Brown /* Read from ex_get_init() */ 39552f7358SJed Brown char title[PETSC_MAX_PATH_LEN+1]; 40552f7358SJed Brown int dim = 0, numVertices = 0, numCells = 0; 41552f7358SJed Brown int num_cs = 0, num_vs = 0, num_fs = 0; 42552f7358SJed Brown #endif 43552f7358SJed Brown 44552f7358SJed Brown PetscFunctionBegin; 45552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 46552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 47552f7358SJed Brown ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 48552f7358SJed Brown ierr = DMCreate(comm, dm);CHKERRQ(ierr); 49552f7358SJed Brown ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 50552f7358SJed Brown /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */ 51552f7358SJed Brown if (!rank) { 5239bba695SBarry Smith ierr = PetscMemzero(title,(PETSC_MAX_PATH_LEN+1)*sizeof(char));CHKERRQ(ierr); 53552f7358SJed Brown ierr = ex_get_init(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs); 54552f7358SJed Brown if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ExodusII ex_get_init() failed with error code %D\n",ierr); 55552f7358SJed Brown if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n"); 56552f7358SJed Brown } 57552f7358SJed Brown ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 58552f7358SJed Brown ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 59552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr); 60552f7358SJed Brown ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 61552f7358SJed Brown ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 62552f7358SJed Brown 63552f7358SJed Brown /* Read cell sets information */ 64552f7358SJed Brown if (!rank) { 65552f7358SJed Brown PetscInt *cone; 66552f7358SJed Brown int c, cs, c_loc, v, v_loc; 67552f7358SJed Brown /* Read from ex_get_elem_blk_ids() */ 68552f7358SJed Brown int *cs_id; 69552f7358SJed Brown /* Read from ex_get_elem_block() */ 70552f7358SJed Brown char buffer[PETSC_MAX_PATH_LEN+1]; 71552f7358SJed Brown int num_cell_in_set, num_vertex_per_cell, num_attr; 72552f7358SJed Brown /* Read from ex_get_elem_conn() */ 73552f7358SJed Brown int *cs_connect; 74552f7358SJed Brown 75552f7358SJed Brown /* Get cell sets IDs */ 76552f7358SJed Brown ierr = PetscMalloc(num_cs * sizeof(int), &cs_id);CHKERRQ(ierr); 77552f7358SJed Brown ierr = ex_get_elem_blk_ids(exoid, cs_id);CHKERRQ(ierr); 78552f7358SJed Brown /* Read the cell set connectivity table and build mesh topology 79552f7358SJed Brown EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 80552f7358SJed Brown /* First set sizes */ 81552f7358SJed Brown for (cs = 0, c = 0; cs < num_cs; cs++) { 82552f7358SJed Brown ierr = ex_get_elem_block(exoid, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, &num_attr);CHKERRQ(ierr); 83552f7358SJed Brown for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 84552f7358SJed Brown ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr); 85552f7358SJed Brown } 86552f7358SJed Brown } 87552f7358SJed Brown ierr = DMSetUp(*dm);CHKERRQ(ierr); 88552f7358SJed Brown for (cs = 0, c = 0; cs < num_cs; cs++) { 89552f7358SJed Brown ierr = ex_get_elem_block(exoid, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, &num_attr);CHKERRQ(ierr); 90*dcca6d9dSJed Brown ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr); 91552f7358SJed Brown ierr = ex_get_elem_conn(exoid, cs_id[cs], cs_connect);CHKERRQ(ierr); 92552f7358SJed Brown /* EXO uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */ 93552f7358SJed Brown for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 94552f7358SJed Brown for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) { 95552f7358SJed Brown cone[v_loc] = cs_connect[v]+numCells-1; 96552f7358SJed Brown } 97731dcdf4SMatthew G. Knepley if (dim == 3) { 982e1b13c2SMatthew G. Knepley /* Tetrahedra are inverted */ 992e1b13c2SMatthew G. Knepley if (num_vertex_per_cell == 4) { 1002e1b13c2SMatthew G. Knepley PetscInt tmp = cone[0]; 1012e1b13c2SMatthew G. Knepley cone[0] = cone[1]; 1022e1b13c2SMatthew G. Knepley cone[1] = tmp; 1032e1b13c2SMatthew G. Knepley } 1042e1b13c2SMatthew G. Knepley /* Hexahedra are inverted */ 1052e1b13c2SMatthew G. Knepley if (num_vertex_per_cell == 8) { 1062e1b13c2SMatthew G. Knepley PetscInt tmp = cone[1]; 1072e1b13c2SMatthew G. Knepley cone[1] = cone[3]; 1082e1b13c2SMatthew G. Knepley cone[3] = tmp; 1092e1b13c2SMatthew G. Knepley } 110731dcdf4SMatthew G. Knepley } 111552f7358SJed Brown ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 112552f7358SJed Brown ierr = DMPlexSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr); 113552f7358SJed Brown } 114552f7358SJed Brown ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr); 115552f7358SJed Brown } 116552f7358SJed Brown ierr = PetscFree(cs_id);CHKERRQ(ierr); 117552f7358SJed Brown } 118552f7358SJed Brown ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 119552f7358SJed Brown ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 120552f7358SJed Brown if (interpolate) { 121552f7358SJed Brown DM idm; 122552f7358SJed Brown 1239f074e33SMatthew G Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 124552f7358SJed Brown /* Maintain Cell Sets label */ 125552f7358SJed Brown { 126552f7358SJed Brown DMLabel label; 127552f7358SJed Brown 128552f7358SJed Brown ierr = DMPlexRemoveLabel(*dm, "Cell Sets", &label);CHKERRQ(ierr); 129552f7358SJed Brown if (label) {ierr = DMPlexAddLabel(idm, label);CHKERRQ(ierr);} 130552f7358SJed Brown } 131552f7358SJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 132552f7358SJed Brown *dm = idm; 133552f7358SJed Brown } 134552f7358SJed Brown 135552f7358SJed Brown /* Create vertex set label */ 136552f7358SJed Brown if (!rank && (num_vs > 0)) { 137552f7358SJed Brown int vs, v; 138552f7358SJed Brown /* Read from ex_get_node_set_ids() */ 139552f7358SJed Brown int *vs_id; 140552f7358SJed Brown /* Read from ex_get_node_set_param() */ 141552f7358SJed Brown int num_vertex_in_set, num_attr; 142552f7358SJed Brown /* Read from ex_get_node_set() */ 143552f7358SJed Brown int *vs_vertex_list; 144552f7358SJed Brown 145552f7358SJed Brown /* Get vertex set ids */ 146552f7358SJed Brown ierr = PetscMalloc(num_vs * sizeof(int), &vs_id);CHKERRQ(ierr); 147552f7358SJed Brown ierr = ex_get_node_set_ids(exoid, vs_id);CHKERRQ(ierr); 148552f7358SJed Brown for (vs = 0; vs < num_vs; ++vs) { 149552f7358SJed Brown ierr = ex_get_node_set_param(exoid, vs_id[vs], &num_vertex_in_set, &num_attr);CHKERRQ(ierr); 150552f7358SJed Brown ierr = PetscMalloc(num_vertex_in_set * sizeof(int), &vs_vertex_list);CHKERRQ(ierr); 151552f7358SJed Brown ierr = ex_get_node_set(exoid, vs_id[vs], vs_vertex_list);CHKERRQ(ierr); 152552f7358SJed Brown for (v = 0; v < num_vertex_in_set; ++v) { 153552f7358SJed Brown ierr = DMPlexSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr); 154552f7358SJed Brown } 155552f7358SJed Brown ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr); 156552f7358SJed Brown } 157552f7358SJed Brown ierr = PetscFree(vs_id);CHKERRQ(ierr); 158552f7358SJed Brown } 159552f7358SJed Brown /* Read coordinates */ 160552f7358SJed Brown ierr = DMPlexGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 161552f7358SJed Brown ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 162552f7358SJed Brown ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 163552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 164552f7358SJed Brown for (v = numCells; v < numCells+numVertices; ++v) { 165552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 166552f7358SJed Brown ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 167552f7358SJed Brown } 168552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 169552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 170552f7358SJed Brown ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 171552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 172552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1732eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 174552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1750aba08f6SKarl Rupp if (!rank) { 176552f7358SJed Brown float *x, *y, *z; 177552f7358SJed Brown 178*dcca6d9dSJed Brown ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr); 179552f7358SJed Brown ierr = ex_get_coord(exoid, x, y, z);CHKERRQ(ierr); 1800d644c17SKarl Rupp if (dim > 0) { 1810d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v]; 1820d644c17SKarl Rupp } 1830d644c17SKarl Rupp if (dim > 1) { 1840d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v]; 1850d644c17SKarl Rupp } 1860d644c17SKarl Rupp if (dim > 2) { 1870d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v]; 1880d644c17SKarl Rupp } 189552f7358SJed Brown ierr = PetscFree3(x,y,z);CHKERRQ(ierr); 190552f7358SJed Brown } 191552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 192552f7358SJed Brown ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 193552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 194552f7358SJed Brown 195552f7358SJed Brown /* Create side set label */ 196552f7358SJed Brown if (!rank && interpolate && (num_fs > 0)) { 197552f7358SJed Brown int fs, f, voff; 198552f7358SJed Brown /* Read from ex_get_side_set_ids() */ 199552f7358SJed Brown int *fs_id; 200552f7358SJed Brown /* Read from ex_get_side_set_param() */ 201552f7358SJed Brown int num_side_in_set, num_dist_fact_in_set; 202552f7358SJed Brown /* Read from ex_get_side_set_node_list() */ 203552f7358SJed Brown int *fs_vertex_count_list, *fs_vertex_list; 204552f7358SJed Brown 205552f7358SJed Brown /* Get side set ids */ 206552f7358SJed Brown ierr = PetscMalloc(num_fs * sizeof(int), &fs_id);CHKERRQ(ierr); 207552f7358SJed Brown ierr = ex_get_side_set_ids(exoid, fs_id);CHKERRQ(ierr); 208552f7358SJed Brown for (fs = 0; fs < num_fs; ++fs) { 209552f7358SJed Brown ierr = ex_get_side_set_param(exoid, fs_id[fs], &num_side_in_set, &num_dist_fact_in_set);CHKERRQ(ierr); 210*dcca6d9dSJed Brown ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr); 211552f7358SJed Brown ierr = ex_get_side_set_node_list(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);CHKERRQ(ierr); 212552f7358SJed Brown for (f = 0, voff = 0; f < num_side_in_set; ++f) { 2130298fd71SBarry Smith const PetscInt *faces = NULL; 214552f7358SJed Brown PetscInt faceSize = fs_vertex_count_list[f], numFaces; 215552f7358SJed Brown PetscInt faceVertices[4], v; 216552f7358SJed Brown 217552f7358SJed Brown if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize); 218552f7358SJed Brown for (v = 0; v < faceSize; ++v, ++voff) { 219552f7358SJed Brown faceVertices[v] = fs_vertex_list[voff]+numCells-1; 220552f7358SJed Brown } 221552f7358SJed Brown ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 222552f7358SJed Brown if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces); 223552f7358SJed Brown ierr = DMPlexSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr); 224552f7358SJed Brown ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 225552f7358SJed Brown } 226552f7358SJed Brown ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr); 227552f7358SJed Brown } 228552f7358SJed Brown ierr = PetscFree(fs_id);CHKERRQ(ierr); 229552f7358SJed Brown } 230552f7358SJed Brown #else 231552f7358SJed Brown SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 232552f7358SJed Brown #endif 233552f7358SJed Brown PetscFunctionReturn(0); 234552f7358SJed Brown } 235