1552f7358SJed Brown #define PETSCDM_DLL 2af0996ceSBarry 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 933751fbdSMichael Lange /*@C 10eaf898f9SPatrick Sanan DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file. 1133751fbdSMichael Lange 1233751fbdSMichael Lange Collective on comm 1333751fbdSMichael Lange 1433751fbdSMichael Lange Input Parameters: 1533751fbdSMichael Lange + comm - The MPI communicator 1633751fbdSMichael Lange . filename - The name of the ExodusII file 1733751fbdSMichael Lange - interpolate - Create faces and edges in the mesh 1833751fbdSMichael Lange 1933751fbdSMichael Lange Output Parameter: 2033751fbdSMichael Lange . dm - The DM object representing the mesh 2133751fbdSMichael Lange 2233751fbdSMichael Lange Level: beginner 2333751fbdSMichael Lange 2433751fbdSMichael Lange .keywords: mesh,ExodusII 2533751fbdSMichael Lange .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus() 2633751fbdSMichael Lange @*/ 2733751fbdSMichael Lange PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 2833751fbdSMichael Lange { 2933751fbdSMichael Lange PetscMPIInt rank; 3033751fbdSMichael Lange PetscErrorCode ierr; 3133751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 32*e84f7031SBlaise Bourdin int CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1; 3333751fbdSMichael Lange float version; 3433751fbdSMichael Lange #endif 3533751fbdSMichael Lange 3633751fbdSMichael Lange PetscFunctionBegin; 3733751fbdSMichael Lange PetscValidCharPointer(filename, 2); 3833751fbdSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 3933751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 4033751fbdSMichael Lange if (!rank) { 4133751fbdSMichael Lange exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version); 4233751fbdSMichael Lange if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename); 4333751fbdSMichael Lange } 4433751fbdSMichael Lange ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr); 4533751fbdSMichael Lange if (!rank) {ierr = ex_close(exoid);CHKERRQ(ierr);} 4633751fbdSMichael Lange #else 4733751fbdSMichael Lange SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 4833751fbdSMichael Lange #endif 4933751fbdSMichael Lange PetscFunctionReturn(0); 5033751fbdSMichael Lange } 5133751fbdSMichael Lange 52552f7358SJed Brown /*@ 5333751fbdSMichael Lange DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID. 54552f7358SJed Brown 55552f7358SJed Brown Collective on comm 56552f7358SJed Brown 57552f7358SJed Brown Input Parameters: 58552f7358SJed Brown + comm - The MPI communicator 59552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 60552f7358SJed Brown - interpolate - Create faces and edges in the mesh 61552f7358SJed Brown 62552f7358SJed Brown Output Parameter: 63552f7358SJed Brown . dm - The DM object representing the mesh 64552f7358SJed Brown 65552f7358SJed Brown Level: beginner 66552f7358SJed Brown 67552f7358SJed Brown .keywords: mesh,ExodusII 68e84e7769SJed Brown .seealso: DMPLEX, DMCreate() 69552f7358SJed Brown @*/ 70552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 71552f7358SJed Brown { 72552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 73552f7358SJed Brown PetscMPIInt num_proc, rank; 74552f7358SJed Brown PetscSection coordSection; 75552f7358SJed Brown Vec coordinates; 76552f7358SJed Brown PetscScalar *coords; 77552f7358SJed Brown PetscInt coordSize, v; 78552f7358SJed Brown PetscErrorCode ierr; 79552f7358SJed Brown /* Read from ex_get_init() */ 80552f7358SJed Brown char title[PETSC_MAX_PATH_LEN+1]; 81552f7358SJed Brown int dim = 0, numVertices = 0, numCells = 0; 82552f7358SJed Brown int num_cs = 0, num_vs = 0, num_fs = 0; 83552f7358SJed Brown #endif 84552f7358SJed Brown 85552f7358SJed Brown PetscFunctionBegin; 86552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 87552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 88552f7358SJed Brown ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 89552f7358SJed Brown ierr = DMCreate(comm, dm);CHKERRQ(ierr); 90552f7358SJed Brown ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 91552f7358SJed Brown /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */ 92552f7358SJed Brown if (!rank) { 9339bba695SBarry Smith ierr = PetscMemzero(title,(PETSC_MAX_PATH_LEN+1)*sizeof(char));CHKERRQ(ierr); 94552f7358SJed Brown ierr = ex_get_init(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs); 95552f7358SJed Brown if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ExodusII ex_get_init() failed with error code %D\n",ierr); 96552f7358SJed Brown if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n"); 97552f7358SJed Brown } 98552f7358SJed Brown ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 99552f7358SJed Brown ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 100552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr); 101c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 102552f7358SJed Brown ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 103552f7358SJed Brown 104552f7358SJed Brown /* Read cell sets information */ 105552f7358SJed Brown if (!rank) { 106552f7358SJed Brown PetscInt *cone; 107552f7358SJed Brown int c, cs, c_loc, v, v_loc; 108552f7358SJed Brown /* Read from ex_get_elem_blk_ids() */ 109552f7358SJed Brown int *cs_id; 110552f7358SJed Brown /* Read from ex_get_elem_block() */ 111552f7358SJed Brown char buffer[PETSC_MAX_PATH_LEN+1]; 112552f7358SJed Brown int num_cell_in_set, num_vertex_per_cell, num_attr; 113552f7358SJed Brown /* Read from ex_get_elem_conn() */ 114552f7358SJed Brown int *cs_connect; 115552f7358SJed Brown 116552f7358SJed Brown /* Get cell sets IDs */ 117785e854fSJed Brown ierr = PetscMalloc1(num_cs, &cs_id);CHKERRQ(ierr); 118f958083aSBlaise Bourdin ierr = ex_get_ids (exoid, EX_ELEM_BLOCK, cs_id);CHKERRQ(ierr); 119552f7358SJed Brown /* Read the cell set connectivity table and build mesh topology 120552f7358SJed Brown EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 121552f7358SJed Brown /* First set sizes */ 122552f7358SJed Brown for (cs = 0, c = 0; cs < num_cs; cs++) { 123f958083aSBlaise Bourdin ierr = ex_get_block(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr);CHKERRQ(ierr); 124552f7358SJed Brown for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 125552f7358SJed Brown ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr); 126552f7358SJed Brown } 127552f7358SJed Brown } 128552f7358SJed Brown ierr = DMSetUp(*dm);CHKERRQ(ierr); 129552f7358SJed Brown for (cs = 0, c = 0; cs < num_cs; cs++) { 130f958083aSBlaise Bourdin ierr = ex_get_block(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr);CHKERRQ(ierr); 131dcca6d9dSJed Brown ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr); 132f958083aSBlaise Bourdin ierr = ex_get_conn (exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL);CHKERRQ(ierr); 133eaf898f9SPatrick Sanan /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 134552f7358SJed Brown for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 135552f7358SJed Brown for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) { 136552f7358SJed Brown cone[v_loc] = cs_connect[v]+numCells-1; 137552f7358SJed Brown } 138731dcdf4SMatthew G. Knepley if (dim == 3) { 1392e1b13c2SMatthew G. Knepley /* Tetrahedra are inverted */ 1402e1b13c2SMatthew G. Knepley if (num_vertex_per_cell == 4) { 1412e1b13c2SMatthew G. Knepley PetscInt tmp = cone[0]; 1422e1b13c2SMatthew G. Knepley cone[0] = cone[1]; 1432e1b13c2SMatthew G. Knepley cone[1] = tmp; 1442e1b13c2SMatthew G. Knepley } 1452e1b13c2SMatthew G. Knepley /* Hexahedra are inverted */ 1462e1b13c2SMatthew G. Knepley if (num_vertex_per_cell == 8) { 1472e1b13c2SMatthew G. Knepley PetscInt tmp = cone[1]; 1482e1b13c2SMatthew G. Knepley cone[1] = cone[3]; 1492e1b13c2SMatthew G. Knepley cone[3] = tmp; 1502e1b13c2SMatthew G. Knepley } 151731dcdf4SMatthew G. Knepley } 152552f7358SJed Brown ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 153c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr); 154552f7358SJed Brown } 155552f7358SJed Brown ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr); 156552f7358SJed Brown } 157552f7358SJed Brown ierr = PetscFree(cs_id);CHKERRQ(ierr); 158552f7358SJed Brown } 159552f7358SJed Brown ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 160552f7358SJed Brown ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 161552f7358SJed Brown if (interpolate) { 1625fd9971aSMatthew G. Knepley DM idm; 163552f7358SJed Brown 1649f074e33SMatthew G Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 165552f7358SJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 166552f7358SJed Brown *dm = idm; 167552f7358SJed Brown } 168552f7358SJed Brown 169552f7358SJed Brown /* Create vertex set label */ 170552f7358SJed Brown if (!rank && (num_vs > 0)) { 171552f7358SJed Brown int vs, v; 172552f7358SJed Brown /* Read from ex_get_node_set_ids() */ 173552f7358SJed Brown int *vs_id; 174552f7358SJed Brown /* Read from ex_get_node_set_param() */ 175f958083aSBlaise Bourdin int num_vertex_in_set; 176552f7358SJed Brown /* Read from ex_get_node_set() */ 177552f7358SJed Brown int *vs_vertex_list; 178552f7358SJed Brown 179552f7358SJed Brown /* Get vertex set ids */ 180785e854fSJed Brown ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr); 181f958083aSBlaise Bourdin ierr = ex_get_ids(exoid, EX_NODE_SET, vs_id);CHKERRQ(ierr); 182552f7358SJed Brown for (vs = 0; vs < num_vs; ++vs) { 183f958083aSBlaise Bourdin ierr = ex_get_set_param(exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL); 184785e854fSJed Brown ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr); 185f958083aSBlaise Bourdin ierr = ex_get_set(exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL);CHKERRQ(ierr); 186552f7358SJed Brown for (v = 0; v < num_vertex_in_set; ++v) { 187c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr); 188552f7358SJed Brown } 189552f7358SJed Brown ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr); 190552f7358SJed Brown } 191552f7358SJed Brown ierr = PetscFree(vs_id);CHKERRQ(ierr); 192552f7358SJed Brown } 193552f7358SJed Brown /* Read coordinates */ 19469d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 195552f7358SJed Brown ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 196552f7358SJed Brown ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 197552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 198552f7358SJed Brown for (v = numCells; v < numCells+numVertices; ++v) { 199552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 200552f7358SJed Brown ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 201552f7358SJed Brown } 202552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 203552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 2048b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 205552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 206552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2074e90ef8eSMatthew G. Knepley ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 2082eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 209552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 2100aba08f6SKarl Rupp if (!rank) { 211*e84f7031SBlaise Bourdin PetscReal *x, *y, *z; 212552f7358SJed Brown 213dcca6d9dSJed Brown ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr); 214552f7358SJed Brown ierr = ex_get_coord(exoid, x, y, z);CHKERRQ(ierr); 2150d644c17SKarl Rupp if (dim > 0) { 2160d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v]; 2170d644c17SKarl Rupp } 2180d644c17SKarl Rupp if (dim > 1) { 2190d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v]; 2200d644c17SKarl Rupp } 2210d644c17SKarl Rupp if (dim > 2) { 2220d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v]; 2230d644c17SKarl Rupp } 224552f7358SJed Brown ierr = PetscFree3(x,y,z);CHKERRQ(ierr); 225552f7358SJed Brown } 226552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 227552f7358SJed Brown ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 228552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 229552f7358SJed Brown 230552f7358SJed Brown /* Create side set label */ 231552f7358SJed Brown if (!rank && interpolate && (num_fs > 0)) { 232552f7358SJed Brown int fs, f, voff; 233552f7358SJed Brown /* Read from ex_get_side_set_ids() */ 234552f7358SJed Brown int *fs_id; 235552f7358SJed Brown /* Read from ex_get_side_set_param() */ 236f958083aSBlaise Bourdin int num_side_in_set; 237552f7358SJed Brown /* Read from ex_get_side_set_node_list() */ 238552f7358SJed Brown int *fs_vertex_count_list, *fs_vertex_list; 239ef073283Smichael_afanasiev /* Read side set labels */ 240ef073283Smichael_afanasiev char fs_name[MAX_STR_LENGTH+1]; 241552f7358SJed Brown 242552f7358SJed Brown /* Get side set ids */ 243785e854fSJed Brown ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr); 2446f1149eeSBlaise Bourdin ierr = ex_get_ids(exoid, EX_SIDE_SET, fs_id);CHKERRQ(ierr); 245552f7358SJed Brown for (fs = 0; fs < num_fs; ++fs) { 246f958083aSBlaise Bourdin ierr = ex_get_set_param(exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL);CHKERRQ(ierr); 247dcca6d9dSJed Brown ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr); 248552f7358SJed Brown ierr = ex_get_side_set_node_list(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);CHKERRQ(ierr); 249ef073283Smichael_afanasiev /* Get the specific name associated with this side set ID. */ 250ef073283Smichael_afanasiev int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name); 251552f7358SJed Brown for (f = 0, voff = 0; f < num_side_in_set; ++f) { 2520298fd71SBarry Smith const PetscInt *faces = NULL; 253552f7358SJed Brown PetscInt faceSize = fs_vertex_count_list[f], numFaces; 254552f7358SJed Brown PetscInt faceVertices[4], v; 255552f7358SJed Brown 256552f7358SJed Brown if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize); 257552f7358SJed Brown for (v = 0; v < faceSize; ++v, ++voff) { 258552f7358SJed Brown faceVertices[v] = fs_vertex_list[voff]+numCells-1; 259552f7358SJed Brown } 260552f7358SJed Brown ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 261552f7358SJed Brown if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces); 262c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr); 263ef073283Smichael_afanasiev /* Only add the label if one has been detected for this side set. */ 264ef073283Smichael_afanasiev if (!fs_name_err) { 265ef073283Smichael_afanasiev ierr = DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]);CHKERRQ(ierr); 266ef073283Smichael_afanasiev } 267552f7358SJed Brown ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 268552f7358SJed Brown } 269552f7358SJed Brown ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr); 270552f7358SJed Brown } 271552f7358SJed Brown ierr = PetscFree(fs_id);CHKERRQ(ierr); 272552f7358SJed Brown } 273552f7358SJed Brown #else 274552f7358SJed Brown SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 275552f7358SJed Brown #endif 276552f7358SJed Brown PetscFunctionReturn(0); 277552f7358SJed Brown } 278