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 9e0266925SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 1078569bb4SMatthew G. Knepley /* 1178569bb4SMatthew G. Knepley EXOGetVarIndex - Locate a result in an exodus file based on its name 1278569bb4SMatthew G. Knepley 1378569bb4SMatthew G. Knepley Not collective 1478569bb4SMatthew G. Knepley 1578569bb4SMatthew G. Knepley Input Parameters: 1678569bb4SMatthew G. Knepley + exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 1778569bb4SMatthew G. Knepley . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK 1878569bb4SMatthew G. Knepley - name - the name of the result 1978569bb4SMatthew G. Knepley 2078569bb4SMatthew G. Knepley Output Parameters: 2178569bb4SMatthew G. Knepley . varIndex - the location in the exodus file of the result 2278569bb4SMatthew G. Knepley 2378569bb4SMatthew G. Knepley Notes: 2478569bb4SMatthew G. Knepley The exodus variable index is obtained by comparing name and the 2578569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance if name is "V" 2678569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 2778569bb4SMatthew G. Knepley amongst all variables of type obj_type. 2878569bb4SMatthew G. Knepley 2978569bb4SMatthew G. Knepley Level: beginner 3078569bb4SMatthew G. Knepley 3178569bb4SMatthew G. Knepley .keywords: mesh,ExodusII 3278569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex(),DMView_Exodus(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadNodal_PlexEXO(),VecLoadZonal_PlexEXO() 3378569bb4SMatthew G. Knepley */ 3478569bb4SMatthew G. Knepley static PetscErrorCode EXOGetVarIndex_Private(int exoid, ex_entity_type obj_type, const char name[], int *varIndex) 3578569bb4SMatthew G. Knepley { 3678569bb4SMatthew G. Knepley int exoerr, num_vars, i, j; 3778569bb4SMatthew G. Knepley char ext_name[MAX_STR_LENGTH+1], var_name[MAX_STR_LENGTH+1]; 3878569bb4SMatthew G. Knepley const int num_suffix = 5; 3978569bb4SMatthew G. Knepley char *suffix[5]; 4078569bb4SMatthew G. Knepley PetscBool flg; 4178569bb4SMatthew G. Knepley PetscErrorCode ierr; 4278569bb4SMatthew G. Knepley 4378569bb4SMatthew G. Knepley PetscFunctionBegin; 4478569bb4SMatthew G. Knepley suffix[0] = (char *) ""; 4578569bb4SMatthew G. Knepley suffix[1] = (char *) "_X"; 4678569bb4SMatthew G. Knepley suffix[2] = (char *) "_XX"; 4778569bb4SMatthew G. Knepley suffix[3] = (char *) "_1"; 4878569bb4SMatthew G. Knepley suffix[4] = (char *) "_11"; 4978569bb4SMatthew G. Knepley 5078569bb4SMatthew G. Knepley *varIndex = 0; 5178569bb4SMatthew G. Knepley exoerr = ex_get_variable_param(exoid, obj_type, &num_vars); 5278569bb4SMatthew G. Knepley for (i = 0; i < num_vars; ++i) { 5378569bb4SMatthew G. Knepley exoerr = ex_get_variable_name(exoid, obj_type, i+1, var_name); 5478569bb4SMatthew G. Knepley for (j = 0; j < num_suffix; ++j){ 5578569bb4SMatthew G. Knepley ierr = PetscStrncpy(ext_name, name, MAX_STR_LENGTH);CHKERRQ(ierr); 5678569bb4SMatthew G. Knepley ierr = PetscStrncat(ext_name, suffix[j], MAX_STR_LENGTH);CHKERRQ(ierr); 5778569bb4SMatthew G. Knepley ierr = PetscStrcasecmp(ext_name, var_name, &flg); 5878569bb4SMatthew G. Knepley if (flg) { 5978569bb4SMatthew G. Knepley *varIndex = i+1; 6078569bb4SMatthew G. Knepley PetscFunctionReturn(0); 6178569bb4SMatthew G. Knepley } 6278569bb4SMatthew G. Knepley } 6378569bb4SMatthew G. Knepley } 6478569bb4SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to locate variable %s in ExodusII file.", name); 6578569bb4SMatthew G. Knepley PetscFunctionReturn(-1); 6678569bb4SMatthew G. Knepley } 6778569bb4SMatthew G. Knepley 6878569bb4SMatthew G. Knepley /* 6978569bb4SMatthew G. Knepley DMPlexView_ExodusII_Internal - Write a DM to disk in exodus format 7078569bb4SMatthew G. Knepley 7178569bb4SMatthew G. Knepley Collective on dm 7278569bb4SMatthew G. Knepley 7378569bb4SMatthew G. Knepley Input Parameters: 7478569bb4SMatthew G. Knepley + dm - The dm to be written 7578569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 7678569bb4SMatthew G. Knepley - degree - the degree of the interpolation space 7778569bb4SMatthew G. Knepley 7878569bb4SMatthew G. Knepley Notes: 7978569bb4SMatthew G. Knepley Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels) 8078569bb4SMatthew G. Knepley consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted. 8178569bb4SMatthew G. Knepley 8278569bb4SMatthew G. Knepley If the dm has been distributed and exoid points to different files on each MPI rank, only the "local" part of the DM 8378569bb4SMatthew G. Knepley (including "ghost" cells and vertices) will be written. If all exoid points to the same file, the resulting file will 8478569bb4SMatthew G. Knepley probably be corrupted. 8578569bb4SMatthew G. Knepley 8678569bb4SMatthew G. Knepley DMPlex only represents geometry while most post-processing software expect that a mesh also provides information 8778569bb4SMatthew G. Knepley on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2. 8878569bb4SMatthew G. Knepley It should be extended to use PetscFE objects. 8978569bb4SMatthew G. Knepley 9078569bb4SMatthew G. Knepley This function will only handle TRI, TET, QUAD and HEX cells. 9178569bb4SMatthew G. Knepley Level: beginner 9278569bb4SMatthew G. Knepley 9378569bb4SMatthew G. Knepley .keywords: mesh, ExodusII 9478569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadNodal_PlexEXO(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO() 9578569bb4SMatthew G. Knepley */ 9678569bb4SMatthew G. Knepley PetscErrorCode DMPlexView_ExodusII_Internal(DM dm, int exoid, PetscInt degree) 9778569bb4SMatthew G. Knepley { 9878569bb4SMatthew G. Knepley enum ElemType {TRI, QUAD, TET, HEX}; 9978569bb4SMatthew G. Knepley MPI_Comm comm; 10078569bb4SMatthew G. Knepley /* Connectivity Variables */ 10178569bb4SMatthew G. Knepley PetscInt cellsNotInConnectivity; 10278569bb4SMatthew G. Knepley /* Cell Sets */ 10378569bb4SMatthew G. Knepley DMLabel csLabel; 10478569bb4SMatthew G. Knepley IS csIS; 10578569bb4SMatthew G. Knepley const PetscInt *csIdx; 10678569bb4SMatthew G. Knepley PetscInt num_cs, cs; 10778569bb4SMatthew G. Knepley enum ElemType *type; 10878569bb4SMatthew G. Knepley /* Coordinate Variables */ 10978569bb4SMatthew G. Knepley DM cdm; 11078569bb4SMatthew G. Knepley PetscSection section; 11178569bb4SMatthew G. Knepley Vec coord; 11278569bb4SMatthew G. Knepley PetscInt **nodes; 113eae8a387SMatthew G. Knepley PetscInt depth, d, dim, skipCells = 0; 11478569bb4SMatthew G. Knepley PetscInt pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes; 11578569bb4SMatthew G. Knepley PetscInt num_vs, num_fs; 11678569bb4SMatthew G. Knepley PetscMPIInt rank, size; 11778569bb4SMatthew G. Knepley const char *dmName; 11878569bb4SMatthew G. Knepley PetscErrorCode ierr; 11978569bb4SMatthew G. Knepley 12078569bb4SMatthew G. Knepley PetscFunctionBegin; 12178569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 12278569bb4SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 12378569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 12478569bb4SMatthew G. Knepley /* --- Get DM info --- */ 12578569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm, &dmName);CHKERRQ(ierr); 12678569bb4SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 12778569bb4SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 12878569bb4SMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 12978569bb4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 13078569bb4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 13178569bb4SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr); 13278569bb4SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 13378569bb4SMatthew G. Knepley numCells = cEnd - cStart; 13478569bb4SMatthew G. Knepley numEdges = eEnd - eStart; 13578569bb4SMatthew G. Knepley numVertices = vEnd - vStart; 13678569bb4SMatthew G. Knepley if (depth == 3) {numFaces = fEnd - fStart;} 13778569bb4SMatthew G. Knepley else {numFaces = 0;} 13878569bb4SMatthew G. Knepley ierr = DMGetLabelSize(dm, "Cell Sets", &num_cs);CHKERRQ(ierr); 13978569bb4SMatthew G. Knepley ierr = DMGetLabelSize(dm, "Vertex Sets", &num_vs);CHKERRQ(ierr); 14078569bb4SMatthew G. Knepley ierr = DMGetLabelSize(dm, "Face Sets", &num_fs);CHKERRQ(ierr); 14178569bb4SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr); 14278569bb4SMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 14378569bb4SMatthew G. Knepley if (num_cs > 0) { 14478569bb4SMatthew G. Knepley ierr = DMGetLabel(dm, "Cell Sets", &csLabel);CHKERRQ(ierr); 14578569bb4SMatthew G. Knepley ierr = DMLabelGetValueIS(csLabel, &csIS);CHKERRQ(ierr); 14678569bb4SMatthew G. Knepley ierr = ISGetIndices(csIS, &csIdx);CHKERRQ(ierr); 14778569bb4SMatthew G. Knepley } 14878569bb4SMatthew G. Knepley ierr = PetscMalloc1(num_cs, &nodes);CHKERRQ(ierr); 14978569bb4SMatthew G. Knepley /* Set element type for each block and compute total number of nodes */ 15078569bb4SMatthew G. Knepley ierr = PetscMalloc1(num_cs, &type);CHKERRQ(ierr); 15178569bb4SMatthew G. Knepley numNodes = numVertices; 15278569bb4SMatthew G. Knepley if (degree == 2) {numNodes += numEdges;} 15378569bb4SMatthew G. Knepley cellsNotInConnectivity = numCells; 15478569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 15578569bb4SMatthew G. Knepley IS stratumIS; 15678569bb4SMatthew G. Knepley const PetscInt *cells; 15778569bb4SMatthew G. Knepley PetscScalar *xyz = NULL; 15878569bb4SMatthew G. Knepley PetscInt csSize, closureSize; 15978569bb4SMatthew G. Knepley PetscInt nodesTriP1[4] = {3,0,0,0}; 16078569bb4SMatthew G. Knepley PetscInt nodesTriP2[4] = {3,3,0,0}; 16178569bb4SMatthew G. Knepley PetscInt nodesQuadP1[4] = {4,0,0,0}; 16278569bb4SMatthew G. Knepley PetscInt nodesQuadP2[4] = {4,4,0,1}; 16378569bb4SMatthew G. Knepley PetscInt nodesTetP1[4] = {4,0,0,0}; 16478569bb4SMatthew G. Knepley PetscInt nodesTetP2[4] = {4,6,0,0}; 16578569bb4SMatthew G. Knepley PetscInt nodesHexP1[4] = {8,0,0,0}; 16678569bb4SMatthew G. Knepley PetscInt nodesHexP2[4] = {8,12,6,1}; 16778569bb4SMatthew G. Knepley 16878569bb4SMatthew G. Knepley ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 16978569bb4SMatthew G. Knepley ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 17078569bb4SMatthew G. Knepley ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 17178569bb4SMatthew G. Knepley ierr = DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr); 17278569bb4SMatthew G. Knepley switch (dim) { 17378569bb4SMatthew G. Knepley case 2: 17478569bb4SMatthew G. Knepley if (closureSize == 3*dim) {type[cs] = TRI;} 17578569bb4SMatthew G. Knepley else if (closureSize == 4*dim) {type[cs] = QUAD;} 17678569bb4SMatthew G. Knepley else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim); 17778569bb4SMatthew G. Knepley break; 17878569bb4SMatthew G. Knepley case 3: 17978569bb4SMatthew G. Knepley if (closureSize == 4*dim) {type[cs] = TET;} 18078569bb4SMatthew G. Knepley else if (closureSize == 8*dim) {type[cs] = HEX;} 18178569bb4SMatthew G. Knepley else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim); 18278569bb4SMatthew G. Knepley break; 18378569bb4SMatthew G. Knepley default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not handled by ExodusII viewer", dim); 18478569bb4SMatthew G. Knepley } 18578569bb4SMatthew G. Knepley if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;} 18678569bb4SMatthew G. Knepley if ((degree == 2) && (type[cs] == HEX)) {numNodes += csSize; numNodes += numFaces;} 18778569bb4SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr); 18878569bb4SMatthew G. Knepley /* Set nodes and Element type */ 18978569bb4SMatthew G. Knepley if (type[cs] == TRI) { 19078569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTriP1; 19178569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTriP2; 19278569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 19378569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesQuadP1; 19478569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesQuadP2; 19578569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 19678569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTetP1; 19778569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTetP2; 19878569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 19978569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesHexP1; 20078569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesHexP2; 20178569bb4SMatthew G. Knepley } 20278569bb4SMatthew G. Knepley /* Compute the number of cells not in the connectivity table */ 20378569bb4SMatthew G. Knepley cellsNotInConnectivity -= nodes[cs][3]*csSize; 20478569bb4SMatthew G. Knepley 20578569bb4SMatthew G. Knepley ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 20678569bb4SMatthew G. Knepley ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 20778569bb4SMatthew G. Knepley } 20878569bb4SMatthew G. Knepley if (num_cs > 0) {ierr = ex_put_init(exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs);CHKERRQ(ierr);} 20978569bb4SMatthew G. Knepley /* --- Connectivity --- */ 21078569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 21178569bb4SMatthew G. Knepley IS stratumIS; 21278569bb4SMatthew G. Knepley const PetscInt *cells; 21378569bb4SMatthew G. Knepley PetscInt *connect; 214eae8a387SMatthew G. Knepley PetscInt edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0; 21578569bb4SMatthew G. Knepley PetscInt csSize, c, connectSize, closureSize; 21678569bb4SMatthew G. Knepley char *elem_type = NULL; 21778569bb4SMatthew G. Knepley char elem_type_tri3[] = "TRI3", elem_type_quad4[] = "QUAD4"; 21878569bb4SMatthew G. Knepley char elem_type_tri6[] = "TRI6", elem_type_quad9[] = "QUAD9"; 21978569bb4SMatthew G. Knepley char elem_type_tet4[] = "TET4", elem_type_hex8[] = "HEX8"; 22078569bb4SMatthew G. Knepley char elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27"; 22178569bb4SMatthew G. Knepley 22278569bb4SMatthew G. Knepley ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 22378569bb4SMatthew G. Knepley ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 22478569bb4SMatthew G. Knepley ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 22578569bb4SMatthew G. Knepley /* Set Element type */ 22678569bb4SMatthew G. Knepley if (type[cs] == TRI) { 22778569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tri3; 22878569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tri6; 22978569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 23078569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_quad4; 23178569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_quad9; 23278569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 23378569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tet4; 23478569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tet10; 23578569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 23678569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_hex8; 23778569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_hex27; 23878569bb4SMatthew G. Knepley } 23978569bb4SMatthew G. Knepley connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3]; 24078569bb4SMatthew G. Knepley ierr = PetscMalloc1(connectSize, &connect);CHKERRQ(ierr); 24178569bb4SMatthew G. Knepley ierr = ex_put_block(exoid, EX_ELEM_BLOCK, cs, elem_type, csSize, connectSize, 0, 0, 1); 24278569bb4SMatthew G. Knepley /* Find number of vertices, edges, and faces in the closure */ 24378569bb4SMatthew G. Knepley verticesInClosure = nodes[cs][0]; 24478569bb4SMatthew G. Knepley if (depth > 1) { 24578569bb4SMatthew G. Knepley if (dim == 2) { 24678569bb4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cells[0], &edgesInClosure);CHKERRQ(ierr); 24778569bb4SMatthew G. Knepley } else if (dim == 3) { 24878569bb4SMatthew G. Knepley PetscInt *closure = NULL; 24978569bb4SMatthew G. Knepley 25078569bb4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cells[0], &facesInClosure);CHKERRQ(ierr); 25178569bb4SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 25278569bb4SMatthew G. Knepley edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure; 25378569bb4SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 25478569bb4SMatthew G. Knepley } 25578569bb4SMatthew G. Knepley } 25678569bb4SMatthew G. Knepley /* Get connectivity for each cell */ 25778569bb4SMatthew G. Knepley for (c = 0; c < csSize; ++c) { 25878569bb4SMatthew G. Knepley PetscInt *closure = NULL; 25978569bb4SMatthew G. Knepley PetscInt temp, i; 26078569bb4SMatthew G. Knepley 26178569bb4SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 26278569bb4SMatthew G. Knepley for (i = 0; i < connectSize; ++i) { 26378569bb4SMatthew G. Knepley if (i < nodes[cs][0]) {/* Vertices */ 26478569bb4SMatthew G. Knepley connect[i] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1; 26578569bb4SMatthew G. Knepley connect[i] -= cellsNotInConnectivity; 26678569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */ 26778569bb4SMatthew G. Knepley connect[i] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1; 26878569bb4SMatthew G. Knepley if (nodes[cs][2] == 0) connect[i] -= numFaces; 26978569bb4SMatthew G. Knepley connect[i] -= cellsNotInConnectivity; 27078569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]){ /* Cells */ 27178569bb4SMatthew G. Knepley connect[i] = closure[0] + 1; 27278569bb4SMatthew G. Knepley connect[i] -= skipCells; 27378569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]){ /* Faces */ 27478569bb4SMatthew G. Knepley connect[i] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1; 27578569bb4SMatthew G. Knepley connect[i] -= cellsNotInConnectivity; 27678569bb4SMatthew G. Knepley } else { 27778569bb4SMatthew G. Knepley connect[i] = -1; 27878569bb4SMatthew G. Knepley } 27978569bb4SMatthew G. Knepley } 28078569bb4SMatthew G. Knepley /* Tetrahedra are inverted */ 28178569bb4SMatthew G. Knepley if (type[cs] == TET) { 28278569bb4SMatthew G. Knepley temp = connect[0]; connect[0] = connect[1]; connect[1] = temp; 28378569bb4SMatthew G. Knepley if (degree == 2) { 28478569bb4SMatthew G. Knepley temp = connect[5]; connect[5] = connect[6]; connect[6] = temp; 28578569bb4SMatthew G. Knepley temp = connect[7]; connect[7] = connect[8]; connect[8] = temp; 28678569bb4SMatthew G. Knepley } 28778569bb4SMatthew G. Knepley } 28878569bb4SMatthew G. Knepley /* Hexahedra are inverted */ 28978569bb4SMatthew G. Knepley if (type[cs] == HEX) { 29078569bb4SMatthew G. Knepley temp = connect[1]; connect[1] = connect[3]; connect[3] = temp; 29178569bb4SMatthew G. Knepley if (degree == 2) { 29278569bb4SMatthew G. Knepley temp = connect[8]; connect[8] = connect[11]; connect[11] = temp; 29378569bb4SMatthew G. Knepley temp = connect[9]; connect[9] = connect[10]; connect[10] = temp; 29478569bb4SMatthew G. Knepley temp = connect[16]; connect[16] = connect[17]; connect[17] = temp; 29578569bb4SMatthew G. Knepley temp = connect[18]; connect[18] = connect[19]; connect[19] = temp; 29678569bb4SMatthew G. Knepley 29778569bb4SMatthew G. Knepley temp = connect[12]; connect[12] = connect[16]; connect[16] = temp; 29878569bb4SMatthew G. Knepley temp = connect[13]; connect[13] = connect[17]; connect[17] = temp; 29978569bb4SMatthew G. Knepley temp = connect[14]; connect[14] = connect[18]; connect[18] = temp; 30078569bb4SMatthew G. Knepley temp = connect[15]; connect[15] = connect[19]; connect[19] = temp; 30178569bb4SMatthew G. Knepley 30278569bb4SMatthew G. Knepley temp = connect[23]; connect[23] = connect[26]; connect[26] = temp; 30378569bb4SMatthew G. Knepley temp = connect[24]; connect[24] = connect[25]; connect[25] = temp; 30478569bb4SMatthew G. Knepley temp = connect[25]; connect[25] = connect[26]; connect[26] = temp; 30578569bb4SMatthew G. Knepley } 30678569bb4SMatthew G. Knepley } 30778569bb4SMatthew G. Knepley ierr = ex_put_partial_conn(exoid, EX_ELEM_BLOCK, cs, c+1, 1, connect, NULL, NULL); 30878569bb4SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 30978569bb4SMatthew G. Knepley } 31078569bb4SMatthew G. Knepley skipCells += (nodes[cs][3] == 0)*csSize; 31178569bb4SMatthew G. Knepley ierr = PetscFree(connect);CHKERRQ(ierr); 31278569bb4SMatthew G. Knepley ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 31378569bb4SMatthew G. Knepley ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 31478569bb4SMatthew G. Knepley } 31578569bb4SMatthew G. Knepley ierr = PetscFree(type);CHKERRQ(ierr); 31678569bb4SMatthew G. Knepley /* --- Coordinates --- */ 31778569bb4SMatthew G. Knepley ierr = PetscSectionCreate(comm, §ion);CHKERRQ(ierr); 31878569bb4SMatthew G. Knepley ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 31978569bb4SMatthew G. Knepley for (d = 0; d < depth; ++d) { 32078569bb4SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 32178569bb4SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 32278569bb4SMatthew G. Knepley ierr = PetscSectionSetDof(section, p, nodes[0][d] > 0);CHKERRQ(ierr); 32378569bb4SMatthew G. Knepley } 32478569bb4SMatthew G. Knepley } 32578569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 32678569bb4SMatthew G. Knepley IS stratumIS; 32778569bb4SMatthew G. Knepley const PetscInt *cells; 32878569bb4SMatthew G. Knepley PetscInt csSize, c; 32978569bb4SMatthew G. Knepley 33078569bb4SMatthew G. Knepley ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 33178569bb4SMatthew G. Knepley ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 33278569bb4SMatthew G. Knepley ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 33378569bb4SMatthew G. Knepley for (c = 0; c < csSize; ++c) { 33478569bb4SMatthew G. Knepley ierr = PetscSectionSetDof(section, cells[c], nodes[cs][3] > 0);CHKERRQ(ierr); 33578569bb4SMatthew G. Knepley } 33678569bb4SMatthew G. Knepley ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 33778569bb4SMatthew G. Knepley ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 33878569bb4SMatthew G. Knepley } 33978569bb4SMatthew G. Knepley if (num_cs > 0) { 34078569bb4SMatthew G. Knepley ierr = ISRestoreIndices(csIS, &csIdx);CHKERRQ(ierr); 34178569bb4SMatthew G. Knepley ierr = ISDestroy(&csIS);CHKERRQ(ierr); 34278569bb4SMatthew G. Knepley } 34378569bb4SMatthew G. Knepley ierr = PetscFree(nodes);CHKERRQ(ierr); 34478569bb4SMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 34578569bb4SMatthew G. Knepley if (numNodes > 0) { 34678569bb4SMatthew G. Knepley const char *coordNames[3] = {"x", "y", "z"}; 34778569bb4SMatthew G. Knepley PetscScalar *coords, *closure; 34878569bb4SMatthew G. Knepley PetscReal *cval; 34978569bb4SMatthew G. Knepley PetscInt hasDof, n = 0; 35078569bb4SMatthew G. Knepley 35178569bb4SMatthew G. Knepley /* There can't be more than 24 values in the closure of a point for the coord section */ 35278569bb4SMatthew G. Knepley ierr = PetscMalloc3(numNodes*3, &coords, dim, &cval, 24, &closure);CHKERRQ(ierr); 35378569bb4SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr); 35478569bb4SMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 35578569bb4SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 35678569bb4SMatthew G. Knepley ierr = PetscSectionGetDof(section, p, &hasDof); 35778569bb4SMatthew G. Knepley if (hasDof) { 35878569bb4SMatthew G. Knepley PetscInt closureSize = 24, j; 35978569bb4SMatthew G. Knepley 36078569bb4SMatthew G. Knepley ierr = DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure);CHKERRQ(ierr); 36178569bb4SMatthew G. Knepley for (d = 0; d < dim; ++d) { 36278569bb4SMatthew G. Knepley cval[d] = 0.0; 36378569bb4SMatthew G. Knepley for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d]; 36478569bb4SMatthew G. Knepley coords[d*numNodes+n] = cval[d] * dim / closureSize; 36578569bb4SMatthew G. Knepley } 36678569bb4SMatthew G. Knepley ++n; 36778569bb4SMatthew G. Knepley } 36878569bb4SMatthew G. Knepley } 36978569bb4SMatthew G. Knepley ierr = ex_put_coord(exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes]);CHKERRQ(ierr); 37078569bb4SMatthew G. Knepley ierr = PetscFree3(coords, cval, closure);CHKERRQ(ierr); 37178569bb4SMatthew G. Knepley ierr = ex_put_coord_names(exoid, (char **) coordNames);CHKERRQ(ierr); 37278569bb4SMatthew G. Knepley } 37378569bb4SMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 37478569bb4SMatthew G. Knepley PetscFunctionReturn(0); 37578569bb4SMatthew G. Knepley } 37678569bb4SMatthew G. Knepley 37778569bb4SMatthew G. Knepley /* 37878569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file 37978569bb4SMatthew G. Knepley 38078569bb4SMatthew G. Knepley Collective on v 38178569bb4SMatthew G. Knepley 38278569bb4SMatthew G. Knepley Input Parameters: 38378569bb4SMatthew G. Knepley + v - The vector to be written 38478569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 38578569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1) 38678569bb4SMatthew G. Knepley 38778569bb4SMatthew G. Knepley Notes: 38878569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 38978569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 39078569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 39178569bb4SMatthew G. Knepley amongst all nodal variables. 39278569bb4SMatthew G. Knepley 39378569bb4SMatthew G. Knepley Level: beginner 39478569bb4SMatthew G. Knepley 39578569bb4SMatthew G. Knepley .keywords: mesh,ExodusII 39678569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecLoadNodal_PlexEXO(),VecViewZonal_PlexEXO(),VecLoadZonal_PlexEXO() 39778569bb4SMatthew G. Knepley @*/ 39878569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step) 39978569bb4SMatthew G. Knepley { 40078569bb4SMatthew G. Knepley MPI_Comm comm; 40178569bb4SMatthew G. Knepley PetscMPIInt size; 40278569bb4SMatthew G. Knepley DM dm; 40378569bb4SMatthew G. Knepley Vec vNatural, vComp; 40478569bb4SMatthew G. Knepley const PetscReal *varray; 40578569bb4SMatthew G. Knepley const char *vecname; 40678569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 40778569bb4SMatthew G. Knepley PetscBool useNatural; 40878569bb4SMatthew G. Knepley int offset; 40978569bb4SMatthew G. Knepley PetscErrorCode ierr; 41078569bb4SMatthew G. Knepley 41178569bb4SMatthew G. Knepley PetscFunctionBegin; 41278569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 41378569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 41478569bb4SMatthew G. Knepley ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 41578569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 41678569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 41778569bb4SMatthew G. Knepley if (useNatural) { 41878569bb4SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr); 41978569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr); 42078569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr); 42178569bb4SMatthew G. Knepley } else { 42278569bb4SMatthew G. Knepley vNatural = v; 42378569bb4SMatthew G. Knepley } 42478569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 42578569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 42678569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr); 42778569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file\n", vecname); 42878569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 42978569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 43078569bb4SMatthew G. Knepley We assume that they are stored sequentially */ 43178569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 43278569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 43378569bb4SMatthew G. Knepley if (bs == 1) { 43478569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr); 43578569bb4SMatthew G. Knepley ierr = ex_put_partial_var(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);CHKERRQ(ierr); 43678569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr); 43778569bb4SMatthew G. Knepley } else { 43878569bb4SMatthew G. Knepley IS compIS; 43978569bb4SMatthew G. Knepley PetscInt c; 44078569bb4SMatthew G. Knepley 44178569bb4SMatthew G. Knepley ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr); 44278569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 44378569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 44478569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 44578569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr); 44678569bb4SMatthew G. Knepley ierr = ex_put_partial_var(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);CHKERRQ(ierr); 44778569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr); 44878569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 44978569bb4SMatthew G. Knepley } 45078569bb4SMatthew G. Knepley ierr = ISDestroy(&compIS);CHKERRQ(ierr); 45178569bb4SMatthew G. Knepley } 45278569bb4SMatthew G. Knepley if (useNatural) {ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);} 45378569bb4SMatthew G. Knepley PetscFunctionReturn(0); 45478569bb4SMatthew G. Knepley } 45578569bb4SMatthew G. Knepley 45678569bb4SMatthew G. Knepley /* 45778569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file 45878569bb4SMatthew G. Knepley 45978569bb4SMatthew G. Knepley Collective on v 46078569bb4SMatthew G. Knepley 46178569bb4SMatthew G. Knepley Input Parameters: 46278569bb4SMatthew G. Knepley + v - The vector to be written 46378569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 46478569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1) 46578569bb4SMatthew G. Knepley 46678569bb4SMatthew G. Knepley Notes: 46778569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 46878569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 46978569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 47078569bb4SMatthew G. Knepley amongst all nodal variables. 47178569bb4SMatthew G. Knepley 47278569bb4SMatthew G. Knepley Level: beginner 47378569bb4SMatthew G. Knepley 47478569bb4SMatthew G. Knepley .keywords: mesh, ExodusII 47578569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO() 47678569bb4SMatthew G. Knepley */ 47778569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step) 47878569bb4SMatthew G. Knepley { 47978569bb4SMatthew G. Knepley MPI_Comm comm; 48078569bb4SMatthew G. Knepley PetscMPIInt size; 48178569bb4SMatthew G. Knepley DM dm; 48278569bb4SMatthew G. Knepley Vec vNatural, vComp; 48378569bb4SMatthew G. Knepley PetscReal *varray; 48478569bb4SMatthew G. Knepley const char *vecname; 48578569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 48678569bb4SMatthew G. Knepley PetscBool useNatural; 48778569bb4SMatthew G. Knepley int offset; 48878569bb4SMatthew G. Knepley PetscErrorCode ierr; 48978569bb4SMatthew G. Knepley 49078569bb4SMatthew G. Knepley PetscFunctionBegin; 49178569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 49278569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 49378569bb4SMatthew G. Knepley ierr = VecGetDM(v,&dm);CHKERRQ(ierr); 49478569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 49578569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 49678569bb4SMatthew G. Knepley if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 49778569bb4SMatthew G. Knepley else {vNatural = v;} 49878569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 49978569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 50078569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr); 50178569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file", vecname); 50278569bb4SMatthew G. Knepley /* Read local chunk from the file */ 50378569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 50478569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 50578569bb4SMatthew G. Knepley if (bs == 1) { 50678569bb4SMatthew G. Knepley ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr); 50778569bb4SMatthew G. Knepley ierr = ex_get_partial_var(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);CHKERRQ(ierr); 50878569bb4SMatthew G. Knepley ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr); 50978569bb4SMatthew G. Knepley } else { 51078569bb4SMatthew G. Knepley IS compIS; 51178569bb4SMatthew G. Knepley PetscInt c; 51278569bb4SMatthew G. Knepley 51378569bb4SMatthew G. Knepley ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr); 51478569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 51578569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 51678569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 51778569bb4SMatthew G. Knepley ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr); 51878569bb4SMatthew G. Knepley ierr = ex_get_partial_var(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);CHKERRQ(ierr); 51978569bb4SMatthew G. Knepley ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr); 52078569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 52178569bb4SMatthew G. Knepley } 52278569bb4SMatthew G. Knepley ierr = ISDestroy(&compIS);CHKERRQ(ierr); 52378569bb4SMatthew G. Knepley } 52478569bb4SMatthew G. Knepley if (useNatural) { 52578569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr); 52678569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr); 52778569bb4SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr); 52878569bb4SMatthew G. Knepley } 52978569bb4SMatthew G. Knepley PetscFunctionReturn(0); 53078569bb4SMatthew G. Knepley } 53178569bb4SMatthew G. Knepley 53278569bb4SMatthew G. Knepley /* 53378569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file 53478569bb4SMatthew G. Knepley 53578569bb4SMatthew G. Knepley Collective on v 53678569bb4SMatthew G. Knepley 53778569bb4SMatthew G. Knepley Input Parameters: 53878569bb4SMatthew G. Knepley + v - The vector to be written 53978569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 54078569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1) 54178569bb4SMatthew G. Knepley 54278569bb4SMatthew G. Knepley Notes: 54378569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 54478569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 54578569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 54678569bb4SMatthew G. Knepley amongst all zonal variables. 54778569bb4SMatthew G. Knepley 54878569bb4SMatthew G. Knepley Level: beginner 54978569bb4SMatthew G. Knepley 55078569bb4SMatthew G. Knepley .keywords: mesh,ExodusII 55178569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Zonal_Internal() 55278569bb4SMatthew G. Knepley */ 55378569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step) 55478569bb4SMatthew G. Knepley { 55578569bb4SMatthew G. Knepley MPI_Comm comm; 55678569bb4SMatthew G. Knepley PetscMPIInt size; 55778569bb4SMatthew G. Knepley DM dm; 55878569bb4SMatthew G. Knepley Vec vNatural, vComp; 55978569bb4SMatthew G. Knepley const PetscReal *varray; 56078569bb4SMatthew G. Knepley const char *vecname; 56178569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 56278569bb4SMatthew G. Knepley PetscBool useNatural; 56378569bb4SMatthew G. Knepley int offset; 56478569bb4SMatthew G. Knepley IS compIS; 56578569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 56678569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 56778569bb4SMatthew G. Knepley PetscErrorCode ierr; 56878569bb4SMatthew G. Knepley 56978569bb4SMatthew G. Knepley PetscFunctionBegin; 57078569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)v, &comm);CHKERRQ(ierr); 57178569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 57278569bb4SMatthew G. Knepley ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 57378569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 57478569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 57578569bb4SMatthew G. Knepley if (useNatural) { 57678569bb4SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr); 57778569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr); 57878569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr); 57978569bb4SMatthew G. Knepley } else { 58078569bb4SMatthew G. Knepley vNatural = v; 58178569bb4SMatthew G. Knepley } 58278569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 58378569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 58478569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr); 58578569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname); 58678569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 58778569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 58878569bb4SMatthew G. Knepley We assume that they are stored sequentially 58978569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 59078569bb4SMatthew G. Knepley but once the vector has been reordered to natural size, we cannot use the label informations 59178569bb4SMatthew G. Knepley to figure out what to save where. */ 59278569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 59378569bb4SMatthew G. Knepley ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr); 59478569bb4SMatthew G. Knepley ierr = ex_get_ids(exoid, EX_ELEM_BLOCK, csID);CHKERRQ(ierr); 59578569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 59678569bb4SMatthew G. Knepley ex_block block; 59778569bb4SMatthew G. Knepley 59878569bb4SMatthew G. Knepley block.id = csID[set]; 59978569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 60078569bb4SMatthew G. Knepley ierr = ex_get_block_param(exoid, &block);CHKERRQ(ierr); 60178569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 60278569bb4SMatthew G. Knepley } 60378569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 60478569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 60578569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);} 60678569bb4SMatthew G. Knepley for (set = 0; set < numCS; set++) { 60778569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 60878569bb4SMatthew G. Knepley 60978569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 61078569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 61178569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 61278569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 61378569bb4SMatthew G. Knepley if (bs == 1) { 61478569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr); 61578569bb4SMatthew G. Knepley ierr = ex_put_partial_var(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]);CHKERRQ(ierr); 61678569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr); 61778569bb4SMatthew G. Knepley } else { 61878569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 61978569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 62078569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 62178569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr); 62278569bb4SMatthew G. Knepley ierr = ex_put_partial_var(exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]); 62378569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr); 62478569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 62578569bb4SMatthew G. Knepley } 62678569bb4SMatthew G. Knepley } 62778569bb4SMatthew G. Knepley csxs += csSize[set]; 62878569bb4SMatthew G. Knepley } 62978569bb4SMatthew G. Knepley ierr = PetscFree2(csID, csSize);CHKERRQ(ierr); 63078569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);} 63178569bb4SMatthew G. Knepley if (useNatural) {ierr = DMRestoreGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 63278569bb4SMatthew G. Knepley PetscFunctionReturn(0); 63378569bb4SMatthew G. Knepley } 63478569bb4SMatthew G. Knepley 63578569bb4SMatthew G. Knepley /* 63678569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file 63778569bb4SMatthew G. Knepley 63878569bb4SMatthew G. Knepley Collective on v 63978569bb4SMatthew G. Knepley 64078569bb4SMatthew G. Knepley Input Parameters: 64178569bb4SMatthew G. Knepley + v - The vector to be written 64278569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 64378569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1) 64478569bb4SMatthew G. Knepley 64578569bb4SMatthew G. Knepley Notes: 64678569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 64778569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 64878569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 64978569bb4SMatthew G. Knepley amongst all zonal variables. 65078569bb4SMatthew G. Knepley 65178569bb4SMatthew G. Knepley Level: beginner 65278569bb4SMatthew G. Knepley 65378569bb4SMatthew G. Knepley .keywords: mesh,ExodusII 65478569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Zonal_Internal() 65578569bb4SMatthew G. Knepley */ 65678569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step) 65778569bb4SMatthew G. Knepley { 65878569bb4SMatthew G. Knepley MPI_Comm comm; 65978569bb4SMatthew G. Knepley PetscMPIInt size; 66078569bb4SMatthew G. Knepley DM dm; 66178569bb4SMatthew G. Knepley Vec vNatural, vComp; 66278569bb4SMatthew G. Knepley PetscReal *varray; 66378569bb4SMatthew G. Knepley const char *vecname; 66478569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 66578569bb4SMatthew G. Knepley PetscBool useNatural; 66678569bb4SMatthew G. Knepley int offset; 66778569bb4SMatthew G. Knepley IS compIS; 66878569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 66978569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 67078569bb4SMatthew G. Knepley PetscErrorCode ierr; 67178569bb4SMatthew G. Knepley 67278569bb4SMatthew G. Knepley PetscFunctionBegin; 67378569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr); 67478569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 67578569bb4SMatthew G. Knepley ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 67678569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 67778569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 67878569bb4SMatthew G. Knepley if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 67978569bb4SMatthew G. Knepley else {vNatural = v;} 68078569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 68178569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 68278569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr); 68378569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname); 68478569bb4SMatthew G. Knepley /* Read local chunk of the result in the exodus file 68578569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 68678569bb4SMatthew G. Knepley We assume that they are stored sequentially 68778569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 68878569bb4SMatthew G. Knepley but once the vector has been reordered to natural size, we cannot use the label informations 68978569bb4SMatthew G. Knepley to figure out what to save where. */ 69078569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 69178569bb4SMatthew G. Knepley ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr); 69278569bb4SMatthew G. Knepley ierr = ex_get_ids(exoid, EX_ELEM_BLOCK, csID);CHKERRQ(ierr); 69378569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 69478569bb4SMatthew G. Knepley ex_block block; 69578569bb4SMatthew G. Knepley 69678569bb4SMatthew G. Knepley block.id = csID[set]; 69778569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 69878569bb4SMatthew G. Knepley ierr = ex_get_block_param(exoid, &block);CHKERRQ(ierr); 69978569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 70078569bb4SMatthew G. Knepley } 70178569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 70278569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 70378569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);} 70478569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 70578569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 70678569bb4SMatthew G. Knepley 70778569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 70878569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 70978569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 71078569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 71178569bb4SMatthew G. Knepley if (bs == 1) { 71278569bb4SMatthew G. Knepley ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr); 71378569bb4SMatthew G. Knepley ierr = ex_get_partial_var(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]);CHKERRQ(ierr); 71478569bb4SMatthew G. Knepley ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr); 71578569bb4SMatthew G. Knepley } else { 71678569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 71778569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 71878569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 71978569bb4SMatthew G. Knepley ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr); 72078569bb4SMatthew G. Knepley ierr = ex_get_partial_var(exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]); 72178569bb4SMatthew G. Knepley ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr); 72278569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 72378569bb4SMatthew G. Knepley } 72478569bb4SMatthew G. Knepley } 72578569bb4SMatthew G. Knepley csxs += csSize[set]; 72678569bb4SMatthew G. Knepley } 72778569bb4SMatthew G. Knepley ierr = PetscFree2(csID, csSize);CHKERRQ(ierr); 72878569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);} 72978569bb4SMatthew G. Knepley if (useNatural) { 73078569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr); 73178569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr); 73278569bb4SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr); 73378569bb4SMatthew G. Knepley } 73478569bb4SMatthew G. Knepley PetscFunctionReturn(0); 73578569bb4SMatthew G. Knepley } 736*b53c8456SSatish Balay #endif 73778569bb4SMatthew G. Knepley 73833751fbdSMichael Lange /*@C 739eaf898f9SPatrick Sanan DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file. 74033751fbdSMichael Lange 74133751fbdSMichael Lange Collective on comm 74233751fbdSMichael Lange 74333751fbdSMichael Lange Input Parameters: 74433751fbdSMichael Lange + comm - The MPI communicator 74533751fbdSMichael Lange . filename - The name of the ExodusII file 74633751fbdSMichael Lange - interpolate - Create faces and edges in the mesh 74733751fbdSMichael Lange 74833751fbdSMichael Lange Output Parameter: 74933751fbdSMichael Lange . dm - The DM object representing the mesh 75033751fbdSMichael Lange 75133751fbdSMichael Lange Level: beginner 75233751fbdSMichael Lange 75333751fbdSMichael Lange .keywords: mesh,ExodusII 75433751fbdSMichael Lange .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus() 75533751fbdSMichael Lange @*/ 75633751fbdSMichael Lange PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 75733751fbdSMichael Lange { 75833751fbdSMichael Lange PetscMPIInt rank; 75933751fbdSMichael Lange PetscErrorCode ierr; 76033751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 761e84f7031SBlaise Bourdin int CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1; 76233751fbdSMichael Lange float version; 76333751fbdSMichael Lange #endif 76433751fbdSMichael Lange 76533751fbdSMichael Lange PetscFunctionBegin; 76633751fbdSMichael Lange PetscValidCharPointer(filename, 2); 76733751fbdSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 76833751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 76933751fbdSMichael Lange if (!rank) { 77033751fbdSMichael Lange exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version); 77133751fbdSMichael Lange if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename); 77233751fbdSMichael Lange } 77333751fbdSMichael Lange ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr); 77433751fbdSMichael Lange if (!rank) {ierr = ex_close(exoid);CHKERRQ(ierr);} 77533751fbdSMichael Lange #else 77633751fbdSMichael Lange SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 77733751fbdSMichael Lange #endif 77833751fbdSMichael Lange PetscFunctionReturn(0); 77933751fbdSMichael Lange } 78033751fbdSMichael Lange 781552f7358SJed Brown /*@ 78233751fbdSMichael Lange DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID. 783552f7358SJed Brown 784552f7358SJed Brown Collective on comm 785552f7358SJed Brown 786552f7358SJed Brown Input Parameters: 787552f7358SJed Brown + comm - The MPI communicator 788552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 789552f7358SJed Brown - interpolate - Create faces and edges in the mesh 790552f7358SJed Brown 791552f7358SJed Brown Output Parameter: 792552f7358SJed Brown . dm - The DM object representing the mesh 793552f7358SJed Brown 794552f7358SJed Brown Level: beginner 795552f7358SJed Brown 796552f7358SJed Brown .keywords: mesh,ExodusII 797e84e7769SJed Brown .seealso: DMPLEX, DMCreate() 798552f7358SJed Brown @*/ 799552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 800552f7358SJed Brown { 801552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 802552f7358SJed Brown PetscMPIInt num_proc, rank; 803552f7358SJed Brown PetscSection coordSection; 804552f7358SJed Brown Vec coordinates; 805552f7358SJed Brown PetscScalar *coords; 806552f7358SJed Brown PetscInt coordSize, v; 807552f7358SJed Brown PetscErrorCode ierr; 808552f7358SJed Brown /* Read from ex_get_init() */ 809552f7358SJed Brown char title[PETSC_MAX_PATH_LEN+1]; 810552f7358SJed Brown int dim = 0, numVertices = 0, numCells = 0; 811552f7358SJed Brown int num_cs = 0, num_vs = 0, num_fs = 0; 812552f7358SJed Brown #endif 813552f7358SJed Brown 814552f7358SJed Brown PetscFunctionBegin; 815552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 816552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 817552f7358SJed Brown ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 818552f7358SJed Brown ierr = DMCreate(comm, dm);CHKERRQ(ierr); 819552f7358SJed Brown ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 820552f7358SJed Brown /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */ 821552f7358SJed Brown if (!rank) { 82239bba695SBarry Smith ierr = PetscMemzero(title,(PETSC_MAX_PATH_LEN+1)*sizeof(char));CHKERRQ(ierr); 823552f7358SJed Brown ierr = ex_get_init(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs); 824552f7358SJed Brown if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ExodusII ex_get_init() failed with error code %D\n",ierr); 825552f7358SJed Brown if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n"); 826552f7358SJed Brown } 827552f7358SJed Brown ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 828552f7358SJed Brown ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 829552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr); 830c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 831552f7358SJed Brown ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 832552f7358SJed Brown 833552f7358SJed Brown /* Read cell sets information */ 834552f7358SJed Brown if (!rank) { 835552f7358SJed Brown PetscInt *cone; 836552f7358SJed Brown int c, cs, c_loc, v, v_loc; 837552f7358SJed Brown /* Read from ex_get_elem_blk_ids() */ 838552f7358SJed Brown int *cs_id; 839552f7358SJed Brown /* Read from ex_get_elem_block() */ 840552f7358SJed Brown char buffer[PETSC_MAX_PATH_LEN+1]; 841552f7358SJed Brown int num_cell_in_set, num_vertex_per_cell, num_attr; 842552f7358SJed Brown /* Read from ex_get_elem_conn() */ 843552f7358SJed Brown int *cs_connect; 844552f7358SJed Brown 845552f7358SJed Brown /* Get cell sets IDs */ 846785e854fSJed Brown ierr = PetscMalloc1(num_cs, &cs_id);CHKERRQ(ierr); 847f958083aSBlaise Bourdin ierr = ex_get_ids (exoid, EX_ELEM_BLOCK, cs_id);CHKERRQ(ierr); 848552f7358SJed Brown /* Read the cell set connectivity table and build mesh topology 849552f7358SJed Brown EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 850552f7358SJed Brown /* First set sizes */ 851552f7358SJed Brown for (cs = 0, c = 0; cs < num_cs; cs++) { 852f958083aSBlaise 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); 853552f7358SJed Brown for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 854552f7358SJed Brown ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr); 855552f7358SJed Brown } 856552f7358SJed Brown } 857552f7358SJed Brown ierr = DMSetUp(*dm);CHKERRQ(ierr); 858552f7358SJed Brown for (cs = 0, c = 0; cs < num_cs; cs++) { 859f958083aSBlaise 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); 860dcca6d9dSJed Brown ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr); 861f958083aSBlaise Bourdin ierr = ex_get_conn (exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL);CHKERRQ(ierr); 862eaf898f9SPatrick Sanan /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 863552f7358SJed Brown for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 864552f7358SJed Brown for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) { 865552f7358SJed Brown cone[v_loc] = cs_connect[v]+numCells-1; 866552f7358SJed Brown } 867731dcdf4SMatthew G. Knepley if (dim == 3) { 8682e1b13c2SMatthew G. Knepley /* Tetrahedra are inverted */ 8692e1b13c2SMatthew G. Knepley if (num_vertex_per_cell == 4) { 8702e1b13c2SMatthew G. Knepley PetscInt tmp = cone[0]; 8712e1b13c2SMatthew G. Knepley cone[0] = cone[1]; 8722e1b13c2SMatthew G. Knepley cone[1] = tmp; 8732e1b13c2SMatthew G. Knepley } 8742e1b13c2SMatthew G. Knepley /* Hexahedra are inverted */ 8752e1b13c2SMatthew G. Knepley if (num_vertex_per_cell == 8) { 8762e1b13c2SMatthew G. Knepley PetscInt tmp = cone[1]; 8772e1b13c2SMatthew G. Knepley cone[1] = cone[3]; 8782e1b13c2SMatthew G. Knepley cone[3] = tmp; 8792e1b13c2SMatthew G. Knepley } 880731dcdf4SMatthew G. Knepley } 881552f7358SJed Brown ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 882c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr); 883552f7358SJed Brown } 884552f7358SJed Brown ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr); 885552f7358SJed Brown } 886552f7358SJed Brown ierr = PetscFree(cs_id);CHKERRQ(ierr); 887552f7358SJed Brown } 888552f7358SJed Brown ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 889552f7358SJed Brown ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 890552f7358SJed Brown if (interpolate) { 8915fd9971aSMatthew G. Knepley DM idm; 892552f7358SJed Brown 8939f074e33SMatthew G Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 894552f7358SJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 895552f7358SJed Brown *dm = idm; 896552f7358SJed Brown } 897552f7358SJed Brown 898552f7358SJed Brown /* Create vertex set label */ 899552f7358SJed Brown if (!rank && (num_vs > 0)) { 900552f7358SJed Brown int vs, v; 901552f7358SJed Brown /* Read from ex_get_node_set_ids() */ 902552f7358SJed Brown int *vs_id; 903552f7358SJed Brown /* Read from ex_get_node_set_param() */ 904f958083aSBlaise Bourdin int num_vertex_in_set; 905552f7358SJed Brown /* Read from ex_get_node_set() */ 906552f7358SJed Brown int *vs_vertex_list; 907552f7358SJed Brown 908552f7358SJed Brown /* Get vertex set ids */ 909785e854fSJed Brown ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr); 910f958083aSBlaise Bourdin ierr = ex_get_ids(exoid, EX_NODE_SET, vs_id);CHKERRQ(ierr); 911552f7358SJed Brown for (vs = 0; vs < num_vs; ++vs) { 912f958083aSBlaise Bourdin ierr = ex_get_set_param(exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL); 913785e854fSJed Brown ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr); 914f958083aSBlaise Bourdin ierr = ex_get_set(exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL);CHKERRQ(ierr); 915552f7358SJed Brown for (v = 0; v < num_vertex_in_set; ++v) { 916c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr); 917552f7358SJed Brown } 918552f7358SJed Brown ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr); 919552f7358SJed Brown } 920552f7358SJed Brown ierr = PetscFree(vs_id);CHKERRQ(ierr); 921552f7358SJed Brown } 922552f7358SJed Brown /* Read coordinates */ 92369d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 924552f7358SJed Brown ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 925552f7358SJed Brown ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 926552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 927552f7358SJed Brown for (v = numCells; v < numCells+numVertices; ++v) { 928552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 929552f7358SJed Brown ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 930552f7358SJed Brown } 931552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 932552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 9338b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 934552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 935552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 9364e90ef8eSMatthew G. Knepley ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 9372eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 938552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 9390aba08f6SKarl Rupp if (!rank) { 940e84f7031SBlaise Bourdin PetscReal *x, *y, *z; 941552f7358SJed Brown 942dcca6d9dSJed Brown ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr); 943552f7358SJed Brown ierr = ex_get_coord(exoid, x, y, z);CHKERRQ(ierr); 9440d644c17SKarl Rupp if (dim > 0) { 9450d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v]; 9460d644c17SKarl Rupp } 9470d644c17SKarl Rupp if (dim > 1) { 9480d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v]; 9490d644c17SKarl Rupp } 9500d644c17SKarl Rupp if (dim > 2) { 9510d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v]; 9520d644c17SKarl Rupp } 953552f7358SJed Brown ierr = PetscFree3(x,y,z);CHKERRQ(ierr); 954552f7358SJed Brown } 955552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 956552f7358SJed Brown ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 957552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 958552f7358SJed Brown 959552f7358SJed Brown /* Create side set label */ 960552f7358SJed Brown if (!rank && interpolate && (num_fs > 0)) { 961552f7358SJed Brown int fs, f, voff; 962552f7358SJed Brown /* Read from ex_get_side_set_ids() */ 963552f7358SJed Brown int *fs_id; 964552f7358SJed Brown /* Read from ex_get_side_set_param() */ 965f958083aSBlaise Bourdin int num_side_in_set; 966552f7358SJed Brown /* Read from ex_get_side_set_node_list() */ 967552f7358SJed Brown int *fs_vertex_count_list, *fs_vertex_list; 968ef073283Smichael_afanasiev /* Read side set labels */ 969ef073283Smichael_afanasiev char fs_name[MAX_STR_LENGTH+1]; 970552f7358SJed Brown 971552f7358SJed Brown /* Get side set ids */ 972785e854fSJed Brown ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr); 9736f1149eeSBlaise Bourdin ierr = ex_get_ids(exoid, EX_SIDE_SET, fs_id);CHKERRQ(ierr); 974552f7358SJed Brown for (fs = 0; fs < num_fs; ++fs) { 975f958083aSBlaise Bourdin ierr = ex_get_set_param(exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL);CHKERRQ(ierr); 976dcca6d9dSJed Brown ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr); 977552f7358SJed Brown ierr = ex_get_side_set_node_list(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);CHKERRQ(ierr); 978ef073283Smichael_afanasiev /* Get the specific name associated with this side set ID. */ 979ef073283Smichael_afanasiev int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name); 980552f7358SJed Brown for (f = 0, voff = 0; f < num_side_in_set; ++f) { 9810298fd71SBarry Smith const PetscInt *faces = NULL; 982552f7358SJed Brown PetscInt faceSize = fs_vertex_count_list[f], numFaces; 983552f7358SJed Brown PetscInt faceVertices[4], v; 984552f7358SJed Brown 985552f7358SJed Brown if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize); 986552f7358SJed Brown for (v = 0; v < faceSize; ++v, ++voff) { 987552f7358SJed Brown faceVertices[v] = fs_vertex_list[voff]+numCells-1; 988552f7358SJed Brown } 989552f7358SJed Brown ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 990552f7358SJed Brown if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces); 991c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr); 992ef073283Smichael_afanasiev /* Only add the label if one has been detected for this side set. */ 993ef073283Smichael_afanasiev if (!fs_name_err) { 994ef073283Smichael_afanasiev ierr = DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]);CHKERRQ(ierr); 995ef073283Smichael_afanasiev } 996552f7358SJed Brown ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 997552f7358SJed Brown } 998552f7358SJed Brown ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr); 999552f7358SJed Brown } 1000552f7358SJed Brown ierr = PetscFree(fs_id);CHKERRQ(ierr); 1001552f7358SJed Brown } 1002552f7358SJed Brown #else 1003552f7358SJed Brown SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 1004552f7358SJed Brown #endif 1005552f7358SJed Brown PetscFunctionReturn(0); 1006552f7358SJed Brown } 1007