1552f7358SJed Brown #define PETSCDM_DLL 2552f7358SJed Brown #include <petsc-private/pleximpl.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 PETSC_EXTERN PetscErrorCode DMPlexInterpolate_2D(DM, DM *); 10552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexInterpolate_3D(DM, DM *); 11552f7358SJed Brown 12552f7358SJed Brown #undef __FUNCT__ 13552f7358SJed Brown #define __FUNCT__ "DMPlexCreateExodus" 14552f7358SJed Brown /*@ 15552f7358SJed Brown DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file. 16552f7358SJed Brown 17552f7358SJed Brown Collective on comm 18552f7358SJed Brown 19552f7358SJed Brown Input Parameters: 20552f7358SJed Brown + comm - The MPI communicator 21552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 22552f7358SJed Brown - interpolate - Create faces and edges in the mesh 23552f7358SJed Brown 24552f7358SJed Brown Output Parameter: 25552f7358SJed Brown . dm - The DM object representing the mesh 26552f7358SJed Brown 27552f7358SJed Brown Level: beginner 28552f7358SJed Brown 29552f7358SJed Brown .keywords: mesh,ExodusII 30552f7358SJed Brown .seealso: MeshCreate(), MeshCreateExodus() 31552f7358SJed Brown @*/ 32552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 33552f7358SJed Brown { 34552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 35552f7358SJed Brown PetscMPIInt num_proc, rank; 36552f7358SJed Brown PetscSection coordSection; 37552f7358SJed Brown Vec coordinates; 38552f7358SJed Brown PetscScalar *coords; 39552f7358SJed Brown PetscInt coordSize, v; 40552f7358SJed Brown PetscErrorCode ierr; 41552f7358SJed Brown /* Read from ex_get_init() */ 42552f7358SJed Brown char title[PETSC_MAX_PATH_LEN+1]; 43552f7358SJed Brown int dim = 0, numVertices = 0, numCells = 0; 44552f7358SJed Brown int num_cs = 0, num_vs = 0, num_fs = 0; 45552f7358SJed Brown #endif 46552f7358SJed Brown 47552f7358SJed Brown PetscFunctionBegin; 48552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 49552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 50552f7358SJed Brown ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 51552f7358SJed Brown ierr = DMCreate(comm, dm);CHKERRQ(ierr); 52552f7358SJed Brown ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 53552f7358SJed Brown /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */ 54552f7358SJed Brown if (!rank) { 55552f7358SJed Brown ierr = ex_get_init(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs); 56552f7358SJed Brown if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ExodusII ex_get_init() failed with error code %D\n",ierr); 57552f7358SJed Brown if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n"); 58552f7358SJed Brown } 59552f7358SJed Brown ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 60552f7358SJed Brown ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 61552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr); 62552f7358SJed Brown ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 63552f7358SJed Brown ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 64552f7358SJed Brown 65552f7358SJed Brown /* Read cell sets information */ 66552f7358SJed Brown if (!rank) { 67552f7358SJed Brown PetscInt *cone; 68552f7358SJed Brown int c, cs, c_loc, v, v_loc; 69552f7358SJed Brown /* Read from ex_get_elem_blk_ids() */ 70552f7358SJed Brown int *cs_id; 71552f7358SJed Brown /* Read from ex_get_elem_block() */ 72552f7358SJed Brown char buffer[PETSC_MAX_PATH_LEN+1]; 73552f7358SJed Brown int num_cell_in_set, num_vertex_per_cell, num_attr; 74552f7358SJed Brown /* Read from ex_get_elem_conn() */ 75552f7358SJed Brown int *cs_connect; 76552f7358SJed Brown 77552f7358SJed Brown /* Get cell sets IDs */ 78552f7358SJed Brown ierr = PetscMalloc(num_cs * sizeof(int), &cs_id);CHKERRQ(ierr); 79552f7358SJed Brown ierr = ex_get_elem_blk_ids(exoid, cs_id);CHKERRQ(ierr); 80552f7358SJed Brown /* Read the cell set connectivity table and build mesh topology 81552f7358SJed Brown EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 82552f7358SJed Brown /* First set sizes */ 83552f7358SJed Brown for (cs = 0, c = 0; cs < num_cs; cs++) { 84552f7358SJed Brown ierr = ex_get_elem_block(exoid, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, &num_attr);CHKERRQ(ierr); 85552f7358SJed Brown for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 86552f7358SJed Brown ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr); 87552f7358SJed Brown } 88552f7358SJed Brown } 89552f7358SJed Brown ierr = DMSetUp(*dm);CHKERRQ(ierr); 90552f7358SJed Brown for (cs = 0, c = 0; cs < num_cs; cs++) { 91552f7358SJed Brown ierr = ex_get_elem_block(exoid, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, &num_attr);CHKERRQ(ierr); 92552f7358SJed Brown ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,int,&cs_connect,num_vertex_per_cell,PetscInt,&cone);CHKERRQ(ierr); 93552f7358SJed Brown ierr = ex_get_elem_conn(exoid, cs_id[cs], cs_connect);CHKERRQ(ierr); 94552f7358SJed Brown /* EXO uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */ 95552f7358SJed Brown for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 96552f7358SJed Brown for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) { 97552f7358SJed Brown cone[v_loc] = cs_connect[v]+numCells-1; 98552f7358SJed Brown } 99552f7358SJed Brown ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 100552f7358SJed Brown ierr = DMPlexSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr); 101552f7358SJed Brown } 102552f7358SJed Brown ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr); 103552f7358SJed Brown } 104552f7358SJed Brown ierr = PetscFree(cs_id);CHKERRQ(ierr); 105552f7358SJed Brown } 106552f7358SJed Brown ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 107552f7358SJed Brown ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 108552f7358SJed Brown if (interpolate) { 109552f7358SJed Brown DM idm; 110552f7358SJed Brown 111552f7358SJed Brown switch (dim) { 112552f7358SJed Brown case 2: 113552f7358SJed Brown ierr = DMPlexInterpolate_2D(*dm, &idm);CHKERRQ(ierr);break; 114552f7358SJed Brown case 3: 115552f7358SJed Brown ierr = DMPlexInterpolate_3D(*dm, &idm);CHKERRQ(ierr);break; 116552f7358SJed Brown default: 117552f7358SJed Brown SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No mesh interpolation support for dimension %D", dim); 118552f7358SJed Brown } 119552f7358SJed Brown /* Maintain Cell Sets label */ 120552f7358SJed Brown { 121552f7358SJed Brown DMLabel label; 122552f7358SJed Brown 123552f7358SJed Brown ierr = DMPlexRemoveLabel(*dm, "Cell Sets", &label);CHKERRQ(ierr); 124552f7358SJed Brown if (label) {ierr = DMPlexAddLabel(idm, label);CHKERRQ(ierr);} 125552f7358SJed Brown } 126552f7358SJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 127552f7358SJed Brown *dm = idm; 128552f7358SJed Brown } 129552f7358SJed Brown 130552f7358SJed Brown /* Create vertex set label */ 131552f7358SJed Brown if (!rank && (num_vs > 0)) { 132552f7358SJed Brown int vs, v; 133552f7358SJed Brown /* Read from ex_get_node_set_ids() */ 134552f7358SJed Brown int *vs_id; 135552f7358SJed Brown /* Read from ex_get_node_set_param() */ 136552f7358SJed Brown int num_vertex_in_set, num_attr; 137552f7358SJed Brown /* Read from ex_get_node_set() */ 138552f7358SJed Brown int *vs_vertex_list; 139552f7358SJed Brown 140552f7358SJed Brown /* Get vertex set ids */ 141552f7358SJed Brown ierr = PetscMalloc(num_vs * sizeof(int), &vs_id);CHKERRQ(ierr); 142552f7358SJed Brown ierr = ex_get_node_set_ids(exoid, vs_id);CHKERRQ(ierr); 143552f7358SJed Brown for (vs = 0; vs < num_vs; ++vs) { 144552f7358SJed Brown ierr = ex_get_node_set_param(exoid, vs_id[vs], &num_vertex_in_set, &num_attr);CHKERRQ(ierr); 145552f7358SJed Brown ierr = PetscMalloc(num_vertex_in_set * sizeof(int), &vs_vertex_list);CHKERRQ(ierr); 146552f7358SJed Brown ierr = ex_get_node_set(exoid, vs_id[vs], vs_vertex_list);CHKERRQ(ierr); 147552f7358SJed Brown for (v = 0; v < num_vertex_in_set; ++v) { 148552f7358SJed Brown ierr = DMPlexSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr); 149552f7358SJed Brown } 150552f7358SJed Brown ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr); 151552f7358SJed Brown } 152552f7358SJed Brown ierr = PetscFree(vs_id);CHKERRQ(ierr); 153552f7358SJed Brown } 154552f7358SJed Brown /* Read coordinates */ 155552f7358SJed Brown ierr = DMPlexGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 156552f7358SJed Brown ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 157552f7358SJed Brown ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 158552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 159552f7358SJed Brown for (v = numCells; v < numCells+numVertices; ++v) { 160552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 161552f7358SJed Brown ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 162552f7358SJed Brown } 163552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 164552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 165552f7358SJed Brown ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 166552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 167552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 168552f7358SJed Brown ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 169552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 170*0aba08f6SKarl Rupp if (!rank) { 171552f7358SJed Brown float *x, *y, *z; 172552f7358SJed Brown 173552f7358SJed Brown ierr = PetscMalloc3(numVertices,float,&x,numVertices,float,&y,numVertices,float,&z);CHKERRQ(ierr); 174552f7358SJed Brown ierr = ex_get_coord(exoid, x, y, z);CHKERRQ(ierr); 175552f7358SJed Brown if (dim > 0) {for (v = 0; v < numVertices; ++v) {coords[v*dim+0] = x[v];}} 176552f7358SJed Brown if (dim > 1) {for (v = 0; v < numVertices; ++v) {coords[v*dim+1] = y[v];}} 177552f7358SJed Brown if (dim > 2) {for (v = 0; v < numVertices; ++v) {coords[v*dim+2] = z[v];}} 178552f7358SJed Brown ierr = PetscFree3(x,y,z);CHKERRQ(ierr); 179552f7358SJed Brown } 180552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 181552f7358SJed Brown ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 182552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 183552f7358SJed Brown 184552f7358SJed Brown /* Create side set label */ 185552f7358SJed Brown if (!rank && interpolate && (num_fs > 0)) { 186552f7358SJed Brown int fs, f, voff; 187552f7358SJed Brown /* Read from ex_get_side_set_ids() */ 188552f7358SJed Brown int *fs_id; 189552f7358SJed Brown /* Read from ex_get_side_set_param() */ 190552f7358SJed Brown int num_side_in_set, num_dist_fact_in_set; 191552f7358SJed Brown /* Read from ex_get_side_set_node_list() */ 192552f7358SJed Brown int *fs_vertex_count_list, *fs_vertex_list; 193552f7358SJed Brown 194552f7358SJed Brown /* Get side set ids */ 195552f7358SJed Brown ierr = PetscMalloc(num_fs * sizeof(int), &fs_id);CHKERRQ(ierr); 196552f7358SJed Brown ierr = ex_get_side_set_ids(exoid, fs_id);CHKERRQ(ierr); 197552f7358SJed Brown for (fs = 0; fs < num_fs; ++fs) { 198552f7358SJed Brown ierr = ex_get_side_set_param(exoid, fs_id[fs], &num_side_in_set, &num_dist_fact_in_set);CHKERRQ(ierr); 199552f7358SJed Brown ierr = PetscMalloc2(num_side_in_set,int,&fs_vertex_count_list,num_side_in_set*4,int,&fs_vertex_list);CHKERRQ(ierr); 200552f7358SJed Brown ierr = ex_get_side_set_node_list(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);CHKERRQ(ierr); 201552f7358SJed Brown for (f = 0, voff = 0; f < num_side_in_set; ++f) { 202552f7358SJed Brown const PetscInt *faces = PETSC_NULL; 203552f7358SJed Brown PetscInt faceSize = fs_vertex_count_list[f], numFaces; 204552f7358SJed Brown PetscInt faceVertices[4], v; 205552f7358SJed Brown 206552f7358SJed Brown if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize); 207552f7358SJed Brown for (v = 0; v < faceSize; ++v, ++voff) { 208552f7358SJed Brown faceVertices[v] = fs_vertex_list[voff]+numCells-1; 209552f7358SJed Brown } 210552f7358SJed Brown ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 211552f7358SJed Brown if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces); 212552f7358SJed Brown ierr = DMPlexSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr); 213552f7358SJed Brown ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 214552f7358SJed Brown } 215552f7358SJed Brown ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr); 216552f7358SJed Brown } 217552f7358SJed Brown ierr = PetscFree(fs_id);CHKERRQ(ierr); 218552f7358SJed Brown } 219552f7358SJed Brown #else 220552f7358SJed Brown SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 221552f7358SJed Brown #endif 222552f7358SJed Brown PetscFunctionReturn(0); 223552f7358SJed Brown } 224