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 9*78569bb4SMatthew G. Knepley /* 10*78569bb4SMatthew G. Knepley EXOGetVarIndex - Locate a result in an exodus file based on its name 11*78569bb4SMatthew G. Knepley 12*78569bb4SMatthew G. Knepley Not collective 13*78569bb4SMatthew G. Knepley 14*78569bb4SMatthew G. Knepley Input Parameters: 15*78569bb4SMatthew G. Knepley + exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 16*78569bb4SMatthew G. Knepley . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK 17*78569bb4SMatthew G. Knepley - name - the name of the result 18*78569bb4SMatthew G. Knepley 19*78569bb4SMatthew G. Knepley Output Parameters: 20*78569bb4SMatthew G. Knepley . varIndex - the location in the exodus file of the result 21*78569bb4SMatthew G. Knepley 22*78569bb4SMatthew G. Knepley Notes: 23*78569bb4SMatthew G. Knepley The exodus variable index is obtained by comparing name and the 24*78569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance if name is "V" 25*78569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 26*78569bb4SMatthew G. Knepley amongst all variables of type obj_type. 27*78569bb4SMatthew G. Knepley 28*78569bb4SMatthew G. Knepley Level: beginner 29*78569bb4SMatthew G. Knepley 30*78569bb4SMatthew G. Knepley .keywords: mesh,ExodusII 31*78569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex(),DMView_Exodus(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadNodal_PlexEXO(),VecLoadZonal_PlexEXO() 32*78569bb4SMatthew G. Knepley */ 33*78569bb4SMatthew G. Knepley static PetscErrorCode EXOGetVarIndex_Private(int exoid, ex_entity_type obj_type, const char name[], int *varIndex) 34*78569bb4SMatthew G. Knepley { 35*78569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 36*78569bb4SMatthew G. Knepley int exoerr, num_vars, i, j; 37*78569bb4SMatthew G. Knepley char ext_name[MAX_STR_LENGTH+1], var_name[MAX_STR_LENGTH+1]; 38*78569bb4SMatthew G. Knepley const int num_suffix = 5; 39*78569bb4SMatthew G. Knepley char *suffix[5]; 40*78569bb4SMatthew G. Knepley PetscBool flg; 41*78569bb4SMatthew G. Knepley PetscErrorCode ierr; 42*78569bb4SMatthew G. Knepley #endif 43*78569bb4SMatthew G. Knepley 44*78569bb4SMatthew G. Knepley PetscFunctionBegin; 45*78569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 46*78569bb4SMatthew G. Knepley suffix[0] = (char *) ""; 47*78569bb4SMatthew G. Knepley suffix[1] = (char *) "_X"; 48*78569bb4SMatthew G. Knepley suffix[2] = (char *) "_XX"; 49*78569bb4SMatthew G. Knepley suffix[3] = (char *) "_1"; 50*78569bb4SMatthew G. Knepley suffix[4] = (char *) "_11"; 51*78569bb4SMatthew G. Knepley 52*78569bb4SMatthew G. Knepley *varIndex = 0; 53*78569bb4SMatthew G. Knepley exoerr = ex_get_variable_param(exoid, obj_type, &num_vars); 54*78569bb4SMatthew G. Knepley for (i = 0; i < num_vars; ++i) { 55*78569bb4SMatthew G. Knepley exoerr = ex_get_variable_name(exoid, obj_type, i+1, var_name); 56*78569bb4SMatthew G. Knepley for (j = 0; j < num_suffix; ++j){ 57*78569bb4SMatthew G. Knepley ierr = PetscStrncpy(ext_name, name, MAX_STR_LENGTH);CHKERRQ(ierr); 58*78569bb4SMatthew G. Knepley ierr = PetscStrncat(ext_name, suffix[j], MAX_STR_LENGTH);CHKERRQ(ierr); 59*78569bb4SMatthew G. Knepley ierr = PetscStrcasecmp(ext_name, var_name, &flg); 60*78569bb4SMatthew G. Knepley if (flg) { 61*78569bb4SMatthew G. Knepley *varIndex = i+1; 62*78569bb4SMatthew G. Knepley PetscFunctionReturn(0); 63*78569bb4SMatthew G. Knepley } 64*78569bb4SMatthew G. Knepley } 65*78569bb4SMatthew G. Knepley } 66*78569bb4SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to locate variable %s in ExodusII file.", name); 67*78569bb4SMatthew G. Knepley #else 68*78569bb4SMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 69*78569bb4SMatthew G. Knepley #endif 70*78569bb4SMatthew G. Knepley PetscFunctionReturn(-1); 71*78569bb4SMatthew G. Knepley } 72*78569bb4SMatthew G. Knepley 73*78569bb4SMatthew G. Knepley /* 74*78569bb4SMatthew G. Knepley DMPlexView_ExodusII_Internal - Write a DM to disk in exodus format 75*78569bb4SMatthew G. Knepley 76*78569bb4SMatthew G. Knepley Collective on dm 77*78569bb4SMatthew G. Knepley 78*78569bb4SMatthew G. Knepley Input Parameters: 79*78569bb4SMatthew G. Knepley + dm - The dm to be written 80*78569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 81*78569bb4SMatthew G. Knepley - degree - the degree of the interpolation space 82*78569bb4SMatthew G. Knepley 83*78569bb4SMatthew G. Knepley Notes: 84*78569bb4SMatthew G. Knepley Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels) 85*78569bb4SMatthew G. Knepley consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted. 86*78569bb4SMatthew G. Knepley 87*78569bb4SMatthew 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 88*78569bb4SMatthew G. Knepley (including "ghost" cells and vertices) will be written. If all exoid points to the same file, the resulting file will 89*78569bb4SMatthew G. Knepley probably be corrupted. 90*78569bb4SMatthew G. Knepley 91*78569bb4SMatthew G. Knepley DMPlex only represents geometry while most post-processing software expect that a mesh also provides information 92*78569bb4SMatthew G. Knepley on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2. 93*78569bb4SMatthew G. Knepley It should be extended to use PetscFE objects. 94*78569bb4SMatthew G. Knepley 95*78569bb4SMatthew G. Knepley This function will only handle TRI, TET, QUAD and HEX cells. 96*78569bb4SMatthew G. Knepley Level: beginner 97*78569bb4SMatthew G. Knepley 98*78569bb4SMatthew G. Knepley .keywords: mesh, ExodusII 99*78569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadNodal_PlexEXO(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO() 100*78569bb4SMatthew G. Knepley */ 101*78569bb4SMatthew G. Knepley PetscErrorCode DMPlexView_ExodusII_Internal(DM dm, int exoid, PetscInt degree) 102*78569bb4SMatthew G. Knepley { 103*78569bb4SMatthew G. Knepley enum ElemType {TRI, QUAD, TET, HEX}; 104*78569bb4SMatthew G. Knepley MPI_Comm comm; 105*78569bb4SMatthew G. Knepley /* Connectivity Variables */ 106*78569bb4SMatthew G. Knepley PetscInt cellsNotInConnectivity; 107*78569bb4SMatthew G. Knepley /* Cell Sets */ 108*78569bb4SMatthew G. Knepley DMLabel csLabel; 109*78569bb4SMatthew G. Knepley IS csIS; 110*78569bb4SMatthew G. Knepley const PetscInt *csIdx; 111*78569bb4SMatthew G. Knepley PetscInt num_cs, cs; 112*78569bb4SMatthew G. Knepley enum ElemType *type; 113*78569bb4SMatthew G. Knepley /* Coordinate Variables */ 114*78569bb4SMatthew G. Knepley DM cdm; 115*78569bb4SMatthew G. Knepley PetscSection section; 116*78569bb4SMatthew G. Knepley Vec coord; 117*78569bb4SMatthew G. Knepley PetscInt **nodes; 118*78569bb4SMatthew G. Knepley PetscInt depth, d, dim; 119*78569bb4SMatthew G. Knepley PetscInt pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes; 120*78569bb4SMatthew G. Knepley PetscInt num_vs, num_fs; 121*78569bb4SMatthew G. Knepley PetscMPIInt rank, size; 122*78569bb4SMatthew G. Knepley const char *dmName; 123*78569bb4SMatthew G. Knepley PetscErrorCode ierr; 124*78569bb4SMatthew G. Knepley 125*78569bb4SMatthew G. Knepley PetscFunctionBegin; 126*78569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 127*78569bb4SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 128*78569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 129*78569bb4SMatthew G. Knepley /* --- Get DM info --- */ 130*78569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm, &dmName);CHKERRQ(ierr); 131*78569bb4SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 132*78569bb4SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 133*78569bb4SMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 134*78569bb4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 135*78569bb4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 136*78569bb4SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr); 137*78569bb4SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 138*78569bb4SMatthew G. Knepley numCells = cEnd - cStart; 139*78569bb4SMatthew G. Knepley numEdges = eEnd - eStart; 140*78569bb4SMatthew G. Knepley numVertices = vEnd - vStart; 141*78569bb4SMatthew G. Knepley if (depth == 3) {numFaces = fEnd - fStart;} 142*78569bb4SMatthew G. Knepley else {numFaces = 0;} 143*78569bb4SMatthew G. Knepley ierr = DMGetLabelSize(dm, "Cell Sets", &num_cs);CHKERRQ(ierr); 144*78569bb4SMatthew G. Knepley ierr = DMGetLabelSize(dm, "Vertex Sets", &num_vs);CHKERRQ(ierr); 145*78569bb4SMatthew G. Knepley ierr = DMGetLabelSize(dm, "Face Sets", &num_fs);CHKERRQ(ierr); 146*78569bb4SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr); 147*78569bb4SMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 148*78569bb4SMatthew G. Knepley if (num_cs > 0) { 149*78569bb4SMatthew G. Knepley ierr = DMGetLabel(dm, "Cell Sets", &csLabel);CHKERRQ(ierr); 150*78569bb4SMatthew G. Knepley ierr = DMLabelGetValueIS(csLabel, &csIS);CHKERRQ(ierr); 151*78569bb4SMatthew G. Knepley ierr = ISGetIndices(csIS, &csIdx);CHKERRQ(ierr); 152*78569bb4SMatthew G. Knepley } 153*78569bb4SMatthew G. Knepley ierr = PetscMalloc1(num_cs, &nodes);CHKERRQ(ierr); 154*78569bb4SMatthew G. Knepley /* Set element type for each block and compute total number of nodes */ 155*78569bb4SMatthew G. Knepley ierr = PetscMalloc1(num_cs, &type);CHKERRQ(ierr); 156*78569bb4SMatthew G. Knepley numNodes = numVertices; 157*78569bb4SMatthew G. Knepley if (degree == 2) {numNodes += numEdges;} 158*78569bb4SMatthew G. Knepley cellsNotInConnectivity = numCells; 159*78569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 160*78569bb4SMatthew G. Knepley IS stratumIS; 161*78569bb4SMatthew G. Knepley const PetscInt *cells; 162*78569bb4SMatthew G. Knepley PetscScalar *xyz = NULL; 163*78569bb4SMatthew G. Knepley PetscInt csSize, closureSize; 164*78569bb4SMatthew G. Knepley PetscInt nodesTriP1[4] = {3,0,0,0}; 165*78569bb4SMatthew G. Knepley PetscInt nodesTriP2[4] = {3,3,0,0}; 166*78569bb4SMatthew G. Knepley PetscInt nodesQuadP1[4] = {4,0,0,0}; 167*78569bb4SMatthew G. Knepley PetscInt nodesQuadP2[4] = {4,4,0,1}; 168*78569bb4SMatthew G. Knepley PetscInt nodesTetP1[4] = {4,0,0,0}; 169*78569bb4SMatthew G. Knepley PetscInt nodesTetP2[4] = {4,6,0,0}; 170*78569bb4SMatthew G. Knepley PetscInt nodesHexP1[4] = {8,0,0,0}; 171*78569bb4SMatthew G. Knepley PetscInt nodesHexP2[4] = {8,12,6,1}; 172*78569bb4SMatthew G. Knepley 173*78569bb4SMatthew G. Knepley ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 174*78569bb4SMatthew G. Knepley ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 175*78569bb4SMatthew G. Knepley ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 176*78569bb4SMatthew G. Knepley ierr = DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr); 177*78569bb4SMatthew G. Knepley switch (dim) { 178*78569bb4SMatthew G. Knepley case 2: 179*78569bb4SMatthew G. Knepley if (closureSize == 3*dim) {type[cs] = TRI;} 180*78569bb4SMatthew G. Knepley else if (closureSize == 4*dim) {type[cs] = QUAD;} 181*78569bb4SMatthew 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); 182*78569bb4SMatthew G. Knepley break; 183*78569bb4SMatthew G. Knepley case 3: 184*78569bb4SMatthew G. Knepley if (closureSize == 4*dim) {type[cs] = TET;} 185*78569bb4SMatthew G. Knepley else if (closureSize == 8*dim) {type[cs] = HEX;} 186*78569bb4SMatthew 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); 187*78569bb4SMatthew G. Knepley break; 188*78569bb4SMatthew G. Knepley default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not handled by ExodusII viewer", dim); 189*78569bb4SMatthew G. Knepley } 190*78569bb4SMatthew G. Knepley if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;} 191*78569bb4SMatthew G. Knepley if ((degree == 2) && (type[cs] == HEX)) {numNodes += csSize; numNodes += numFaces;} 192*78569bb4SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr); 193*78569bb4SMatthew G. Knepley /* Set nodes and Element type */ 194*78569bb4SMatthew G. Knepley if (type[cs] == TRI) { 195*78569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTriP1; 196*78569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTriP2; 197*78569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 198*78569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesQuadP1; 199*78569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesQuadP2; 200*78569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 201*78569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTetP1; 202*78569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTetP2; 203*78569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 204*78569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesHexP1; 205*78569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesHexP2; 206*78569bb4SMatthew G. Knepley } 207*78569bb4SMatthew G. Knepley /* Compute the number of cells not in the connectivity table */ 208*78569bb4SMatthew G. Knepley cellsNotInConnectivity -= nodes[cs][3]*csSize; 209*78569bb4SMatthew G. Knepley 210*78569bb4SMatthew G. Knepley ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 211*78569bb4SMatthew G. Knepley ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 212*78569bb4SMatthew G. Knepley } 213*78569bb4SMatthew G. Knepley if (num_cs > 0) {ierr = ex_put_init(exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs);CHKERRQ(ierr);} 214*78569bb4SMatthew G. Knepley /* --- Connectivity --- */ 215*78569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 216*78569bb4SMatthew G. Knepley IS stratumIS; 217*78569bb4SMatthew G. Knepley const PetscInt *cells; 218*78569bb4SMatthew G. Knepley PetscInt *connect; 219*78569bb4SMatthew G. Knepley PetscInt edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0, skipCells = 0; 220*78569bb4SMatthew G. Knepley PetscInt csSize, c, connectSize, closureSize; 221*78569bb4SMatthew G. Knepley char *elem_type = NULL; 222*78569bb4SMatthew G. Knepley char elem_type_tri3[] = "TRI3", elem_type_quad4[] = "QUAD4"; 223*78569bb4SMatthew G. Knepley char elem_type_tri6[] = "TRI6", elem_type_quad9[] = "QUAD9"; 224*78569bb4SMatthew G. Knepley char elem_type_tet4[] = "TET4", elem_type_hex8[] = "HEX8"; 225*78569bb4SMatthew G. Knepley char elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27"; 226*78569bb4SMatthew G. Knepley 227*78569bb4SMatthew G. Knepley ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 228*78569bb4SMatthew G. Knepley ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 229*78569bb4SMatthew G. Knepley ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 230*78569bb4SMatthew G. Knepley /* Set Element type */ 231*78569bb4SMatthew G. Knepley if (type[cs] == TRI) { 232*78569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tri3; 233*78569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tri6; 234*78569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 235*78569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_quad4; 236*78569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_quad9; 237*78569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 238*78569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tet4; 239*78569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tet10; 240*78569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 241*78569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_hex8; 242*78569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_hex27; 243*78569bb4SMatthew G. Knepley } 244*78569bb4SMatthew G. Knepley connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3]; 245*78569bb4SMatthew G. Knepley ierr = PetscMalloc1(connectSize, &connect);CHKERRQ(ierr); 246*78569bb4SMatthew G. Knepley ierr = ex_put_block(exoid, EX_ELEM_BLOCK, cs, elem_type, csSize, connectSize, 0, 0, 1); 247*78569bb4SMatthew G. Knepley /* Find number of vertices, edges, and faces in the closure */ 248*78569bb4SMatthew G. Knepley verticesInClosure = nodes[cs][0]; 249*78569bb4SMatthew G. Knepley if (depth > 1) { 250*78569bb4SMatthew G. Knepley if (dim == 2) { 251*78569bb4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cells[0], &edgesInClosure);CHKERRQ(ierr); 252*78569bb4SMatthew G. Knepley } else if (dim == 3) { 253*78569bb4SMatthew G. Knepley PetscInt *closure = NULL; 254*78569bb4SMatthew G. Knepley 255*78569bb4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cells[0], &facesInClosure);CHKERRQ(ierr); 256*78569bb4SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 257*78569bb4SMatthew G. Knepley edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure; 258*78569bb4SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 259*78569bb4SMatthew G. Knepley } 260*78569bb4SMatthew G. Knepley } 261*78569bb4SMatthew G. Knepley /* Get connectivity for each cell */ 262*78569bb4SMatthew G. Knepley for (c = 0; c < csSize; ++c) { 263*78569bb4SMatthew G. Knepley PetscInt *closure = NULL; 264*78569bb4SMatthew G. Knepley PetscInt temp, i; 265*78569bb4SMatthew G. Knepley 266*78569bb4SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 267*78569bb4SMatthew G. Knepley for (i = 0; i < connectSize; ++i) { 268*78569bb4SMatthew G. Knepley if (i < nodes[cs][0]) {/* Vertices */ 269*78569bb4SMatthew G. Knepley connect[i] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1; 270*78569bb4SMatthew G. Knepley connect[i] -= cellsNotInConnectivity; 271*78569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */ 272*78569bb4SMatthew G. Knepley connect[i] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1; 273*78569bb4SMatthew G. Knepley if (nodes[cs][2] == 0) connect[i] -= numFaces; 274*78569bb4SMatthew G. Knepley connect[i] -= cellsNotInConnectivity; 275*78569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]){ /* Cells */ 276*78569bb4SMatthew G. Knepley connect[i] = closure[0] + 1; 277*78569bb4SMatthew G. Knepley connect[i] -= skipCells; 278*78569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]){ /* Faces */ 279*78569bb4SMatthew G. Knepley connect[i] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1; 280*78569bb4SMatthew G. Knepley connect[i] -= cellsNotInConnectivity; 281*78569bb4SMatthew G. Knepley } else { 282*78569bb4SMatthew G. Knepley connect[i] = -1; 283*78569bb4SMatthew G. Knepley } 284*78569bb4SMatthew G. Knepley } 285*78569bb4SMatthew G. Knepley /* Tetrahedra are inverted */ 286*78569bb4SMatthew G. Knepley if (type[cs] == TET) { 287*78569bb4SMatthew G. Knepley temp = connect[0]; connect[0] = connect[1]; connect[1] = temp; 288*78569bb4SMatthew G. Knepley if (degree == 2) { 289*78569bb4SMatthew G. Knepley temp = connect[5]; connect[5] = connect[6]; connect[6] = temp; 290*78569bb4SMatthew G. Knepley temp = connect[7]; connect[7] = connect[8]; connect[8] = temp; 291*78569bb4SMatthew G. Knepley } 292*78569bb4SMatthew G. Knepley } 293*78569bb4SMatthew G. Knepley /* Hexahedra are inverted */ 294*78569bb4SMatthew G. Knepley if (type[cs] == HEX) { 295*78569bb4SMatthew G. Knepley temp = connect[1]; connect[1] = connect[3]; connect[3] = temp; 296*78569bb4SMatthew G. Knepley if (degree == 2) { 297*78569bb4SMatthew G. Knepley temp = connect[8]; connect[8] = connect[11]; connect[11] = temp; 298*78569bb4SMatthew G. Knepley temp = connect[9]; connect[9] = connect[10]; connect[10] = temp; 299*78569bb4SMatthew G. Knepley temp = connect[16]; connect[16] = connect[17]; connect[17] = temp; 300*78569bb4SMatthew G. Knepley temp = connect[18]; connect[18] = connect[19]; connect[19] = temp; 301*78569bb4SMatthew G. Knepley 302*78569bb4SMatthew G. Knepley temp = connect[12]; connect[12] = connect[16]; connect[16] = temp; 303*78569bb4SMatthew G. Knepley temp = connect[13]; connect[13] = connect[17]; connect[17] = temp; 304*78569bb4SMatthew G. Knepley temp = connect[14]; connect[14] = connect[18]; connect[18] = temp; 305*78569bb4SMatthew G. Knepley temp = connect[15]; connect[15] = connect[19]; connect[19] = temp; 306*78569bb4SMatthew G. Knepley 307*78569bb4SMatthew G. Knepley temp = connect[23]; connect[23] = connect[26]; connect[26] = temp; 308*78569bb4SMatthew G. Knepley temp = connect[24]; connect[24] = connect[25]; connect[25] = temp; 309*78569bb4SMatthew G. Knepley temp = connect[25]; connect[25] = connect[26]; connect[26] = temp; 310*78569bb4SMatthew G. Knepley } 311*78569bb4SMatthew G. Knepley } 312*78569bb4SMatthew G. Knepley ierr = ex_put_partial_conn(exoid, EX_ELEM_BLOCK, cs, c+1, 1, connect, NULL, NULL); 313*78569bb4SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 314*78569bb4SMatthew G. Knepley } 315*78569bb4SMatthew G. Knepley skipCells += (nodes[cs][3] == 0)*csSize; 316*78569bb4SMatthew G. Knepley ierr = PetscFree(connect);CHKERRQ(ierr); 317*78569bb4SMatthew G. Knepley ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 318*78569bb4SMatthew G. Knepley ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 319*78569bb4SMatthew G. Knepley } 320*78569bb4SMatthew G. Knepley ierr = PetscFree(type);CHKERRQ(ierr); 321*78569bb4SMatthew G. Knepley /* --- Coordinates --- */ 322*78569bb4SMatthew G. Knepley ierr = PetscSectionCreate(comm, §ion);CHKERRQ(ierr); 323*78569bb4SMatthew G. Knepley ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 324*78569bb4SMatthew G. Knepley for (d = 0; d < depth; ++d) { 325*78569bb4SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 326*78569bb4SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 327*78569bb4SMatthew G. Knepley ierr = PetscSectionSetDof(section, p, nodes[0][d] > 0);CHKERRQ(ierr); 328*78569bb4SMatthew G. Knepley } 329*78569bb4SMatthew G. Knepley } 330*78569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 331*78569bb4SMatthew G. Knepley IS stratumIS; 332*78569bb4SMatthew G. Knepley const PetscInt *cells; 333*78569bb4SMatthew G. Knepley PetscInt csSize, c; 334*78569bb4SMatthew G. Knepley 335*78569bb4SMatthew G. Knepley ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 336*78569bb4SMatthew G. Knepley ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 337*78569bb4SMatthew G. Knepley ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 338*78569bb4SMatthew G. Knepley for (c = 0; c < csSize; ++c) { 339*78569bb4SMatthew G. Knepley ierr = PetscSectionSetDof(section, cells[c], nodes[cs][3] > 0);CHKERRQ(ierr); 340*78569bb4SMatthew G. Knepley } 341*78569bb4SMatthew G. Knepley ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 342*78569bb4SMatthew G. Knepley ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 343*78569bb4SMatthew G. Knepley } 344*78569bb4SMatthew G. Knepley if (num_cs > 0) { 345*78569bb4SMatthew G. Knepley ierr = ISRestoreIndices(csIS, &csIdx);CHKERRQ(ierr); 346*78569bb4SMatthew G. Knepley ierr = ISDestroy(&csIS);CHKERRQ(ierr); 347*78569bb4SMatthew G. Knepley } 348*78569bb4SMatthew G. Knepley ierr = PetscFree(nodes);CHKERRQ(ierr); 349*78569bb4SMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 350*78569bb4SMatthew G. Knepley if (numNodes > 0) { 351*78569bb4SMatthew G. Knepley const char *coordNames[3] = {"x", "y", "z"}; 352*78569bb4SMatthew G. Knepley PetscScalar *coords, *closure; 353*78569bb4SMatthew G. Knepley PetscReal *cval; 354*78569bb4SMatthew G. Knepley PetscInt hasDof, n = 0; 355*78569bb4SMatthew G. Knepley 356*78569bb4SMatthew G. Knepley /* There can't be more than 24 values in the closure of a point for the coord section */ 357*78569bb4SMatthew G. Knepley ierr = PetscMalloc3(numNodes*3, &coords, dim, &cval, 24, &closure);CHKERRQ(ierr); 358*78569bb4SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr); 359*78569bb4SMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 360*78569bb4SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 361*78569bb4SMatthew G. Knepley ierr = PetscSectionGetDof(section, p, &hasDof); 362*78569bb4SMatthew G. Knepley if (hasDof) { 363*78569bb4SMatthew G. Knepley PetscInt closureSize = 24, j; 364*78569bb4SMatthew G. Knepley 365*78569bb4SMatthew G. Knepley ierr = DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure);CHKERRQ(ierr); 366*78569bb4SMatthew G. Knepley for (d = 0; d < dim; ++d) { 367*78569bb4SMatthew G. Knepley cval[d] = 0.0; 368*78569bb4SMatthew G. Knepley for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d]; 369*78569bb4SMatthew G. Knepley coords[d*numNodes+n] = cval[d] * dim / closureSize; 370*78569bb4SMatthew G. Knepley } 371*78569bb4SMatthew G. Knepley ++n; 372*78569bb4SMatthew G. Knepley } 373*78569bb4SMatthew G. Knepley } 374*78569bb4SMatthew G. Knepley ierr = ex_put_coord(exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes]);CHKERRQ(ierr); 375*78569bb4SMatthew G. Knepley ierr = PetscFree3(coords, cval, closure);CHKERRQ(ierr); 376*78569bb4SMatthew G. Knepley ierr = ex_put_coord_names(exoid, (char **) coordNames);CHKERRQ(ierr); 377*78569bb4SMatthew G. Knepley } 378*78569bb4SMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 379*78569bb4SMatthew G. Knepley PetscFunctionReturn(0); 380*78569bb4SMatthew G. Knepley } 381*78569bb4SMatthew G. Knepley 382*78569bb4SMatthew G. Knepley /* 383*78569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file 384*78569bb4SMatthew G. Knepley 385*78569bb4SMatthew G. Knepley Collective on v 386*78569bb4SMatthew G. Knepley 387*78569bb4SMatthew G. Knepley Input Parameters: 388*78569bb4SMatthew G. Knepley + v - The vector to be written 389*78569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 390*78569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1) 391*78569bb4SMatthew G. Knepley 392*78569bb4SMatthew G. Knepley Notes: 393*78569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 394*78569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 395*78569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 396*78569bb4SMatthew G. Knepley amongst all nodal variables. 397*78569bb4SMatthew G. Knepley 398*78569bb4SMatthew G. Knepley Level: beginner 399*78569bb4SMatthew G. Knepley 400*78569bb4SMatthew G. Knepley .keywords: mesh,ExodusII 401*78569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecLoadNodal_PlexEXO(),VecViewZonal_PlexEXO(),VecLoadZonal_PlexEXO() 402*78569bb4SMatthew G. Knepley @*/ 403*78569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step) 404*78569bb4SMatthew G. Knepley { 405*78569bb4SMatthew G. Knepley MPI_Comm comm; 406*78569bb4SMatthew G. Knepley PetscMPIInt size; 407*78569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 408*78569bb4SMatthew G. Knepley DM dm; 409*78569bb4SMatthew G. Knepley Vec vNatural, vComp; 410*78569bb4SMatthew G. Knepley const PetscReal *varray; 411*78569bb4SMatthew G. Knepley const char *vecname; 412*78569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 413*78569bb4SMatthew G. Knepley PetscBool useNatural; 414*78569bb4SMatthew G. Knepley int offset; 415*78569bb4SMatthew G. Knepley PetscErrorCode ierr; 416*78569bb4SMatthew G. Knepley #endif 417*78569bb4SMatthew G. Knepley 418*78569bb4SMatthew G. Knepley PetscFunctionBegin; 419*78569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 420*78569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 421*78569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 422*78569bb4SMatthew G. Knepley ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 423*78569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 424*78569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 425*78569bb4SMatthew G. Knepley if (useNatural) { 426*78569bb4SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr); 427*78569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr); 428*78569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr); 429*78569bb4SMatthew G. Knepley } else { 430*78569bb4SMatthew G. Knepley vNatural = v; 431*78569bb4SMatthew G. Knepley } 432*78569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 433*78569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 434*78569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr); 435*78569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file\n", vecname); 436*78569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 437*78569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 438*78569bb4SMatthew G. Knepley We assume that they are stored sequentially */ 439*78569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 440*78569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 441*78569bb4SMatthew G. Knepley if (bs == 1) { 442*78569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr); 443*78569bb4SMatthew G. Knepley ierr = ex_put_partial_var(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);CHKERRQ(ierr); 444*78569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr); 445*78569bb4SMatthew G. Knepley } else { 446*78569bb4SMatthew G. Knepley IS compIS; 447*78569bb4SMatthew G. Knepley PetscInt c; 448*78569bb4SMatthew G. Knepley 449*78569bb4SMatthew G. Knepley ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr); 450*78569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 451*78569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 452*78569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 453*78569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr); 454*78569bb4SMatthew G. Knepley ierr = ex_put_partial_var(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);CHKERRQ(ierr); 455*78569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr); 456*78569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 457*78569bb4SMatthew G. Knepley } 458*78569bb4SMatthew G. Knepley ierr = ISDestroy(&compIS);CHKERRQ(ierr); 459*78569bb4SMatthew G. Knepley } 460*78569bb4SMatthew G. Knepley if (useNatural) {ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);} 461*78569bb4SMatthew G. Knepley #else 462*78569bb4SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 463*78569bb4SMatthew G. Knepley #endif 464*78569bb4SMatthew G. Knepley PetscFunctionReturn(0); 465*78569bb4SMatthew G. Knepley } 466*78569bb4SMatthew G. Knepley 467*78569bb4SMatthew G. Knepley /* 468*78569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file 469*78569bb4SMatthew G. Knepley 470*78569bb4SMatthew G. Knepley Collective on v 471*78569bb4SMatthew G. Knepley 472*78569bb4SMatthew G. Knepley Input Parameters: 473*78569bb4SMatthew G. Knepley + v - The vector to be written 474*78569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 475*78569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1) 476*78569bb4SMatthew G. Knepley 477*78569bb4SMatthew G. Knepley Notes: 478*78569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 479*78569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 480*78569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 481*78569bb4SMatthew G. Knepley amongst all nodal variables. 482*78569bb4SMatthew G. Knepley 483*78569bb4SMatthew G. Knepley Level: beginner 484*78569bb4SMatthew G. Knepley 485*78569bb4SMatthew G. Knepley .keywords: mesh, ExodusII 486*78569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO() 487*78569bb4SMatthew G. Knepley */ 488*78569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step) 489*78569bb4SMatthew G. Knepley { 490*78569bb4SMatthew G. Knepley MPI_Comm comm; 491*78569bb4SMatthew G. Knepley PetscMPIInt size; 492*78569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 493*78569bb4SMatthew G. Knepley DM dm; 494*78569bb4SMatthew G. Knepley Vec vNatural, vComp; 495*78569bb4SMatthew G. Knepley PetscReal *varray; 496*78569bb4SMatthew G. Knepley const char *vecname; 497*78569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 498*78569bb4SMatthew G. Knepley PetscBool useNatural; 499*78569bb4SMatthew G. Knepley int offset; 500*78569bb4SMatthew G. Knepley PetscErrorCode ierr; 501*78569bb4SMatthew G. Knepley #endif 502*78569bb4SMatthew G. Knepley 503*78569bb4SMatthew G. Knepley PetscFunctionBegin; 504*78569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 505*78569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 506*78569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 507*78569bb4SMatthew G. Knepley ierr = VecGetDM(v,&dm);CHKERRQ(ierr); 508*78569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 509*78569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 510*78569bb4SMatthew G. Knepley if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 511*78569bb4SMatthew G. Knepley else {vNatural = v;} 512*78569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 513*78569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 514*78569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr); 515*78569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file", vecname); 516*78569bb4SMatthew G. Knepley /* Read local chunk from the file */ 517*78569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 518*78569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 519*78569bb4SMatthew G. Knepley if (bs == 1) { 520*78569bb4SMatthew G. Knepley ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr); 521*78569bb4SMatthew G. Knepley ierr = ex_get_partial_var(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);CHKERRQ(ierr); 522*78569bb4SMatthew G. Knepley ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr); 523*78569bb4SMatthew G. Knepley } else { 524*78569bb4SMatthew G. Knepley IS compIS; 525*78569bb4SMatthew G. Knepley PetscInt c; 526*78569bb4SMatthew G. Knepley 527*78569bb4SMatthew G. Knepley ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr); 528*78569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 529*78569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 530*78569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 531*78569bb4SMatthew G. Knepley ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr); 532*78569bb4SMatthew G. Knepley ierr = ex_get_partial_var(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);CHKERRQ(ierr); 533*78569bb4SMatthew G. Knepley ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr); 534*78569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 535*78569bb4SMatthew G. Knepley } 536*78569bb4SMatthew G. Knepley ierr = ISDestroy(&compIS);CHKERRQ(ierr); 537*78569bb4SMatthew G. Knepley } 538*78569bb4SMatthew G. Knepley if (useNatural) { 539*78569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr); 540*78569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr); 541*78569bb4SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr); 542*78569bb4SMatthew G. Knepley } 543*78569bb4SMatthew G. Knepley #else 544*78569bb4SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 545*78569bb4SMatthew G. Knepley #endif 546*78569bb4SMatthew G. Knepley PetscFunctionReturn(0); 547*78569bb4SMatthew G. Knepley } 548*78569bb4SMatthew G. Knepley 549*78569bb4SMatthew G. Knepley /* 550*78569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file 551*78569bb4SMatthew G. Knepley 552*78569bb4SMatthew G. Knepley Collective on v 553*78569bb4SMatthew G. Knepley 554*78569bb4SMatthew G. Knepley Input Parameters: 555*78569bb4SMatthew G. Knepley + v - The vector to be written 556*78569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 557*78569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1) 558*78569bb4SMatthew G. Knepley 559*78569bb4SMatthew G. Knepley Notes: 560*78569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 561*78569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 562*78569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 563*78569bb4SMatthew G. Knepley amongst all zonal variables. 564*78569bb4SMatthew G. Knepley 565*78569bb4SMatthew G. Knepley Level: beginner 566*78569bb4SMatthew G. Knepley 567*78569bb4SMatthew G. Knepley .keywords: mesh,ExodusII 568*78569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Zonal_Internal() 569*78569bb4SMatthew G. Knepley */ 570*78569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step) 571*78569bb4SMatthew G. Knepley { 572*78569bb4SMatthew G. Knepley MPI_Comm comm; 573*78569bb4SMatthew G. Knepley PetscMPIInt size; 574*78569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 575*78569bb4SMatthew G. Knepley DM dm; 576*78569bb4SMatthew G. Knepley Vec vNatural, vComp; 577*78569bb4SMatthew G. Knepley const PetscReal *varray; 578*78569bb4SMatthew G. Knepley const char *vecname; 579*78569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 580*78569bb4SMatthew G. Knepley PetscBool useNatural; 581*78569bb4SMatthew G. Knepley int offset; 582*78569bb4SMatthew G. Knepley IS compIS; 583*78569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 584*78569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 585*78569bb4SMatthew G. Knepley PetscErrorCode ierr; 586*78569bb4SMatthew G. Knepley #endif 587*78569bb4SMatthew G. Knepley 588*78569bb4SMatthew G. Knepley PetscFunctionBegin; 589*78569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)v, &comm);CHKERRQ(ierr); 590*78569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 591*78569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 592*78569bb4SMatthew G. Knepley ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 593*78569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 594*78569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 595*78569bb4SMatthew G. Knepley if (useNatural) { 596*78569bb4SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr); 597*78569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr); 598*78569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr); 599*78569bb4SMatthew G. Knepley } else { 600*78569bb4SMatthew G. Knepley vNatural = v; 601*78569bb4SMatthew G. Knepley } 602*78569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 603*78569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 604*78569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr); 605*78569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname); 606*78569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 607*78569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 608*78569bb4SMatthew G. Knepley We assume that they are stored sequentially 609*78569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 610*78569bb4SMatthew G. Knepley but once the vector has been reordered to natural size, we cannot use the label informations 611*78569bb4SMatthew G. Knepley to figure out what to save where. */ 612*78569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 613*78569bb4SMatthew G. Knepley ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr); 614*78569bb4SMatthew G. Knepley ierr = ex_get_ids(exoid, EX_ELEM_BLOCK, csID);CHKERRQ(ierr); 615*78569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 616*78569bb4SMatthew G. Knepley ex_block block; 617*78569bb4SMatthew G. Knepley 618*78569bb4SMatthew G. Knepley block.id = csID[set]; 619*78569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 620*78569bb4SMatthew G. Knepley ierr = ex_get_block_param(exoid, &block);CHKERRQ(ierr); 621*78569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 622*78569bb4SMatthew G. Knepley } 623*78569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 624*78569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 625*78569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);} 626*78569bb4SMatthew G. Knepley for (set = 0; set < numCS; set++) { 627*78569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 628*78569bb4SMatthew G. Knepley 629*78569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 630*78569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 631*78569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 632*78569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 633*78569bb4SMatthew G. Knepley if (bs == 1) { 634*78569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr); 635*78569bb4SMatthew 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); 636*78569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr); 637*78569bb4SMatthew G. Knepley } else { 638*78569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 639*78569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 640*78569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 641*78569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr); 642*78569bb4SMatthew 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)]); 643*78569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr); 644*78569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 645*78569bb4SMatthew G. Knepley } 646*78569bb4SMatthew G. Knepley } 647*78569bb4SMatthew G. Knepley csxs += csSize[set]; 648*78569bb4SMatthew G. Knepley } 649*78569bb4SMatthew G. Knepley ierr = PetscFree2(csID, csSize);CHKERRQ(ierr); 650*78569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);} 651*78569bb4SMatthew G. Knepley if (useNatural) {ierr = DMRestoreGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 652*78569bb4SMatthew G. Knepley #else 653*78569bb4SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 654*78569bb4SMatthew G. Knepley #endif 655*78569bb4SMatthew G. Knepley PetscFunctionReturn(0); 656*78569bb4SMatthew G. Knepley } 657*78569bb4SMatthew G. Knepley 658*78569bb4SMatthew G. Knepley /* 659*78569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file 660*78569bb4SMatthew G. Knepley 661*78569bb4SMatthew G. Knepley Collective on v 662*78569bb4SMatthew G. Knepley 663*78569bb4SMatthew G. Knepley Input Parameters: 664*78569bb4SMatthew G. Knepley + v - The vector to be written 665*78569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 666*78569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1) 667*78569bb4SMatthew G. Knepley 668*78569bb4SMatthew G. Knepley Notes: 669*78569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 670*78569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 671*78569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 672*78569bb4SMatthew G. Knepley amongst all zonal variables. 673*78569bb4SMatthew G. Knepley 674*78569bb4SMatthew G. Knepley Level: beginner 675*78569bb4SMatthew G. Knepley 676*78569bb4SMatthew G. Knepley .keywords: mesh,ExodusII 677*78569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Zonal_Internal() 678*78569bb4SMatthew G. Knepley */ 679*78569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step) 680*78569bb4SMatthew G. Knepley { 681*78569bb4SMatthew G. Knepley MPI_Comm comm; 682*78569bb4SMatthew G. Knepley PetscMPIInt size; 683*78569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 684*78569bb4SMatthew G. Knepley DM dm; 685*78569bb4SMatthew G. Knepley Vec vNatural, vComp; 686*78569bb4SMatthew G. Knepley PetscReal *varray; 687*78569bb4SMatthew G. Knepley const char *vecname; 688*78569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 689*78569bb4SMatthew G. Knepley PetscBool useNatural; 690*78569bb4SMatthew G. Knepley int offset; 691*78569bb4SMatthew G. Knepley IS compIS; 692*78569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 693*78569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 694*78569bb4SMatthew G. Knepley PetscErrorCode ierr; 695*78569bb4SMatthew G. Knepley #endif 696*78569bb4SMatthew G. Knepley 697*78569bb4SMatthew G. Knepley PetscFunctionBegin; 698*78569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr); 699*78569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 700*78569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 701*78569bb4SMatthew G. Knepley ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 702*78569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 703*78569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 704*78569bb4SMatthew G. Knepley if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 705*78569bb4SMatthew G. Knepley else {vNatural = v;} 706*78569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 707*78569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 708*78569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr); 709*78569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname); 710*78569bb4SMatthew G. Knepley /* Read local chunk of the result in the exodus file 711*78569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 712*78569bb4SMatthew G. Knepley We assume that they are stored sequentially 713*78569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 714*78569bb4SMatthew G. Knepley but once the vector has been reordered to natural size, we cannot use the label informations 715*78569bb4SMatthew G. Knepley to figure out what to save where. */ 716*78569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 717*78569bb4SMatthew G. Knepley ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr); 718*78569bb4SMatthew G. Knepley ierr = ex_get_ids(exoid, EX_ELEM_BLOCK, csID);CHKERRQ(ierr); 719*78569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 720*78569bb4SMatthew G. Knepley ex_block block; 721*78569bb4SMatthew G. Knepley 722*78569bb4SMatthew G. Knepley block.id = csID[set]; 723*78569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 724*78569bb4SMatthew G. Knepley ierr = ex_get_block_param(exoid, &block);CHKERRQ(ierr); 725*78569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 726*78569bb4SMatthew G. Knepley } 727*78569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 728*78569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 729*78569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);} 730*78569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 731*78569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 732*78569bb4SMatthew G. Knepley 733*78569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 734*78569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 735*78569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 736*78569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 737*78569bb4SMatthew G. Knepley if (bs == 1) { 738*78569bb4SMatthew G. Knepley ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr); 739*78569bb4SMatthew 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); 740*78569bb4SMatthew G. Knepley ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr); 741*78569bb4SMatthew G. Knepley } else { 742*78569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 743*78569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 744*78569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 745*78569bb4SMatthew G. Knepley ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr); 746*78569bb4SMatthew 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)]); 747*78569bb4SMatthew G. Knepley ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr); 748*78569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 749*78569bb4SMatthew G. Knepley } 750*78569bb4SMatthew G. Knepley } 751*78569bb4SMatthew G. Knepley csxs += csSize[set]; 752*78569bb4SMatthew G. Knepley } 753*78569bb4SMatthew G. Knepley ierr = PetscFree2(csID, csSize);CHKERRQ(ierr); 754*78569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);} 755*78569bb4SMatthew G. Knepley if (useNatural) { 756*78569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr); 757*78569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr); 758*78569bb4SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr); 759*78569bb4SMatthew G. Knepley } 760*78569bb4SMatthew G. Knepley #else 761*78569bb4SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 762*78569bb4SMatthew G. Knepley #endif 763*78569bb4SMatthew G. Knepley PetscFunctionReturn(0); 764*78569bb4SMatthew G. Knepley } 765*78569bb4SMatthew G. Knepley 76633751fbdSMichael Lange /*@C 767eaf898f9SPatrick Sanan DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file. 76833751fbdSMichael Lange 76933751fbdSMichael Lange Collective on comm 77033751fbdSMichael Lange 77133751fbdSMichael Lange Input Parameters: 77233751fbdSMichael Lange + comm - The MPI communicator 77333751fbdSMichael Lange . filename - The name of the ExodusII file 77433751fbdSMichael Lange - interpolate - Create faces and edges in the mesh 77533751fbdSMichael Lange 77633751fbdSMichael Lange Output Parameter: 77733751fbdSMichael Lange . dm - The DM object representing the mesh 77833751fbdSMichael Lange 77933751fbdSMichael Lange Level: beginner 78033751fbdSMichael Lange 78133751fbdSMichael Lange .keywords: mesh,ExodusII 78233751fbdSMichael Lange .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus() 78333751fbdSMichael Lange @*/ 78433751fbdSMichael Lange PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 78533751fbdSMichael Lange { 78633751fbdSMichael Lange PetscMPIInt rank; 78733751fbdSMichael Lange PetscErrorCode ierr; 78833751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 789e84f7031SBlaise Bourdin int CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1; 79033751fbdSMichael Lange float version; 79133751fbdSMichael Lange #endif 79233751fbdSMichael Lange 79333751fbdSMichael Lange PetscFunctionBegin; 79433751fbdSMichael Lange PetscValidCharPointer(filename, 2); 79533751fbdSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 79633751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 79733751fbdSMichael Lange if (!rank) { 79833751fbdSMichael Lange exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version); 79933751fbdSMichael Lange if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename); 80033751fbdSMichael Lange } 80133751fbdSMichael Lange ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr); 80233751fbdSMichael Lange if (!rank) {ierr = ex_close(exoid);CHKERRQ(ierr);} 80333751fbdSMichael Lange #else 80433751fbdSMichael Lange SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 80533751fbdSMichael Lange #endif 80633751fbdSMichael Lange PetscFunctionReturn(0); 80733751fbdSMichael Lange } 80833751fbdSMichael Lange 809552f7358SJed Brown /*@ 81033751fbdSMichael Lange DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID. 811552f7358SJed Brown 812552f7358SJed Brown Collective on comm 813552f7358SJed Brown 814552f7358SJed Brown Input Parameters: 815552f7358SJed Brown + comm - The MPI communicator 816552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 817552f7358SJed Brown - interpolate - Create faces and edges in the mesh 818552f7358SJed Brown 819552f7358SJed Brown Output Parameter: 820552f7358SJed Brown . dm - The DM object representing the mesh 821552f7358SJed Brown 822552f7358SJed Brown Level: beginner 823552f7358SJed Brown 824552f7358SJed Brown .keywords: mesh,ExodusII 825e84e7769SJed Brown .seealso: DMPLEX, DMCreate() 826552f7358SJed Brown @*/ 827552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 828552f7358SJed Brown { 829552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 830552f7358SJed Brown PetscMPIInt num_proc, rank; 831552f7358SJed Brown PetscSection coordSection; 832552f7358SJed Brown Vec coordinates; 833552f7358SJed Brown PetscScalar *coords; 834552f7358SJed Brown PetscInt coordSize, v; 835552f7358SJed Brown PetscErrorCode ierr; 836552f7358SJed Brown /* Read from ex_get_init() */ 837552f7358SJed Brown char title[PETSC_MAX_PATH_LEN+1]; 838552f7358SJed Brown int dim = 0, numVertices = 0, numCells = 0; 839552f7358SJed Brown int num_cs = 0, num_vs = 0, num_fs = 0; 840552f7358SJed Brown #endif 841552f7358SJed Brown 842552f7358SJed Brown PetscFunctionBegin; 843552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 844552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 845552f7358SJed Brown ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 846552f7358SJed Brown ierr = DMCreate(comm, dm);CHKERRQ(ierr); 847552f7358SJed Brown ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 848552f7358SJed Brown /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */ 849552f7358SJed Brown if (!rank) { 85039bba695SBarry Smith ierr = PetscMemzero(title,(PETSC_MAX_PATH_LEN+1)*sizeof(char));CHKERRQ(ierr); 851552f7358SJed Brown ierr = ex_get_init(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs); 852552f7358SJed Brown if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ExodusII ex_get_init() failed with error code %D\n",ierr); 853552f7358SJed Brown if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n"); 854552f7358SJed Brown } 855552f7358SJed Brown ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 856552f7358SJed Brown ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 857552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr); 858c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 859552f7358SJed Brown ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 860552f7358SJed Brown 861552f7358SJed Brown /* Read cell sets information */ 862552f7358SJed Brown if (!rank) { 863552f7358SJed Brown PetscInt *cone; 864552f7358SJed Brown int c, cs, c_loc, v, v_loc; 865552f7358SJed Brown /* Read from ex_get_elem_blk_ids() */ 866552f7358SJed Brown int *cs_id; 867552f7358SJed Brown /* Read from ex_get_elem_block() */ 868552f7358SJed Brown char buffer[PETSC_MAX_PATH_LEN+1]; 869552f7358SJed Brown int num_cell_in_set, num_vertex_per_cell, num_attr; 870552f7358SJed Brown /* Read from ex_get_elem_conn() */ 871552f7358SJed Brown int *cs_connect; 872552f7358SJed Brown 873552f7358SJed Brown /* Get cell sets IDs */ 874785e854fSJed Brown ierr = PetscMalloc1(num_cs, &cs_id);CHKERRQ(ierr); 875f958083aSBlaise Bourdin ierr = ex_get_ids (exoid, EX_ELEM_BLOCK, cs_id);CHKERRQ(ierr); 876552f7358SJed Brown /* Read the cell set connectivity table and build mesh topology 877552f7358SJed Brown EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 878552f7358SJed Brown /* First set sizes */ 879552f7358SJed Brown for (cs = 0, c = 0; cs < num_cs; cs++) { 880f958083aSBlaise 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); 881552f7358SJed Brown for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 882552f7358SJed Brown ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr); 883552f7358SJed Brown } 884552f7358SJed Brown } 885552f7358SJed Brown ierr = DMSetUp(*dm);CHKERRQ(ierr); 886552f7358SJed Brown for (cs = 0, c = 0; cs < num_cs; cs++) { 887f958083aSBlaise 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); 888dcca6d9dSJed Brown ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr); 889f958083aSBlaise Bourdin ierr = ex_get_conn (exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL);CHKERRQ(ierr); 890eaf898f9SPatrick Sanan /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 891552f7358SJed Brown for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 892552f7358SJed Brown for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) { 893552f7358SJed Brown cone[v_loc] = cs_connect[v]+numCells-1; 894552f7358SJed Brown } 895731dcdf4SMatthew G. Knepley if (dim == 3) { 8962e1b13c2SMatthew G. Knepley /* Tetrahedra are inverted */ 8972e1b13c2SMatthew G. Knepley if (num_vertex_per_cell == 4) { 8982e1b13c2SMatthew G. Knepley PetscInt tmp = cone[0]; 8992e1b13c2SMatthew G. Knepley cone[0] = cone[1]; 9002e1b13c2SMatthew G. Knepley cone[1] = tmp; 9012e1b13c2SMatthew G. Knepley } 9022e1b13c2SMatthew G. Knepley /* Hexahedra are inverted */ 9032e1b13c2SMatthew G. Knepley if (num_vertex_per_cell == 8) { 9042e1b13c2SMatthew G. Knepley PetscInt tmp = cone[1]; 9052e1b13c2SMatthew G. Knepley cone[1] = cone[3]; 9062e1b13c2SMatthew G. Knepley cone[3] = tmp; 9072e1b13c2SMatthew G. Knepley } 908731dcdf4SMatthew G. Knepley } 909552f7358SJed Brown ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 910c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr); 911552f7358SJed Brown } 912552f7358SJed Brown ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr); 913552f7358SJed Brown } 914552f7358SJed Brown ierr = PetscFree(cs_id);CHKERRQ(ierr); 915552f7358SJed Brown } 916552f7358SJed Brown ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 917552f7358SJed Brown ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 918552f7358SJed Brown if (interpolate) { 9195fd9971aSMatthew G. Knepley DM idm; 920552f7358SJed Brown 9219f074e33SMatthew G Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 922552f7358SJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 923552f7358SJed Brown *dm = idm; 924552f7358SJed Brown } 925552f7358SJed Brown 926552f7358SJed Brown /* Create vertex set label */ 927552f7358SJed Brown if (!rank && (num_vs > 0)) { 928552f7358SJed Brown int vs, v; 929552f7358SJed Brown /* Read from ex_get_node_set_ids() */ 930552f7358SJed Brown int *vs_id; 931552f7358SJed Brown /* Read from ex_get_node_set_param() */ 932f958083aSBlaise Bourdin int num_vertex_in_set; 933552f7358SJed Brown /* Read from ex_get_node_set() */ 934552f7358SJed Brown int *vs_vertex_list; 935552f7358SJed Brown 936552f7358SJed Brown /* Get vertex set ids */ 937785e854fSJed Brown ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr); 938f958083aSBlaise Bourdin ierr = ex_get_ids(exoid, EX_NODE_SET, vs_id);CHKERRQ(ierr); 939552f7358SJed Brown for (vs = 0; vs < num_vs; ++vs) { 940f958083aSBlaise Bourdin ierr = ex_get_set_param(exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL); 941785e854fSJed Brown ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr); 942f958083aSBlaise Bourdin ierr = ex_get_set(exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL);CHKERRQ(ierr); 943552f7358SJed Brown for (v = 0; v < num_vertex_in_set; ++v) { 944c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr); 945552f7358SJed Brown } 946552f7358SJed Brown ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr); 947552f7358SJed Brown } 948552f7358SJed Brown ierr = PetscFree(vs_id);CHKERRQ(ierr); 949552f7358SJed Brown } 950552f7358SJed Brown /* Read coordinates */ 95169d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 952552f7358SJed Brown ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 953552f7358SJed Brown ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 954552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 955552f7358SJed Brown for (v = numCells; v < numCells+numVertices; ++v) { 956552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 957552f7358SJed Brown ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 958552f7358SJed Brown } 959552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 960552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 9618b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 962552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 963552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 9644e90ef8eSMatthew G. Knepley ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 9652eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 966552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 9670aba08f6SKarl Rupp if (!rank) { 968e84f7031SBlaise Bourdin PetscReal *x, *y, *z; 969552f7358SJed Brown 970dcca6d9dSJed Brown ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr); 971552f7358SJed Brown ierr = ex_get_coord(exoid, x, y, z);CHKERRQ(ierr); 9720d644c17SKarl Rupp if (dim > 0) { 9730d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v]; 9740d644c17SKarl Rupp } 9750d644c17SKarl Rupp if (dim > 1) { 9760d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v]; 9770d644c17SKarl Rupp } 9780d644c17SKarl Rupp if (dim > 2) { 9790d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v]; 9800d644c17SKarl Rupp } 981552f7358SJed Brown ierr = PetscFree3(x,y,z);CHKERRQ(ierr); 982552f7358SJed Brown } 983552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 984552f7358SJed Brown ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 985552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 986552f7358SJed Brown 987552f7358SJed Brown /* Create side set label */ 988552f7358SJed Brown if (!rank && interpolate && (num_fs > 0)) { 989552f7358SJed Brown int fs, f, voff; 990552f7358SJed Brown /* Read from ex_get_side_set_ids() */ 991552f7358SJed Brown int *fs_id; 992552f7358SJed Brown /* Read from ex_get_side_set_param() */ 993f958083aSBlaise Bourdin int num_side_in_set; 994552f7358SJed Brown /* Read from ex_get_side_set_node_list() */ 995552f7358SJed Brown int *fs_vertex_count_list, *fs_vertex_list; 996ef073283Smichael_afanasiev /* Read side set labels */ 997ef073283Smichael_afanasiev char fs_name[MAX_STR_LENGTH+1]; 998552f7358SJed Brown 999552f7358SJed Brown /* Get side set ids */ 1000785e854fSJed Brown ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr); 10016f1149eeSBlaise Bourdin ierr = ex_get_ids(exoid, EX_SIDE_SET, fs_id);CHKERRQ(ierr); 1002552f7358SJed Brown for (fs = 0; fs < num_fs; ++fs) { 1003f958083aSBlaise Bourdin ierr = ex_get_set_param(exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL);CHKERRQ(ierr); 1004dcca6d9dSJed Brown ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr); 1005552f7358SJed Brown ierr = ex_get_side_set_node_list(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);CHKERRQ(ierr); 1006ef073283Smichael_afanasiev /* Get the specific name associated with this side set ID. */ 1007ef073283Smichael_afanasiev int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name); 1008552f7358SJed Brown for (f = 0, voff = 0; f < num_side_in_set; ++f) { 10090298fd71SBarry Smith const PetscInt *faces = NULL; 1010552f7358SJed Brown PetscInt faceSize = fs_vertex_count_list[f], numFaces; 1011552f7358SJed Brown PetscInt faceVertices[4], v; 1012552f7358SJed Brown 1013552f7358SJed Brown if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize); 1014552f7358SJed Brown for (v = 0; v < faceSize; ++v, ++voff) { 1015552f7358SJed Brown faceVertices[v] = fs_vertex_list[voff]+numCells-1; 1016552f7358SJed Brown } 1017552f7358SJed Brown ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 1018552f7358SJed Brown if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces); 1019c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr); 1020ef073283Smichael_afanasiev /* Only add the label if one has been detected for this side set. */ 1021ef073283Smichael_afanasiev if (!fs_name_err) { 1022ef073283Smichael_afanasiev ierr = DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]);CHKERRQ(ierr); 1023ef073283Smichael_afanasiev } 1024552f7358SJed Brown ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 1025552f7358SJed Brown } 1026552f7358SJed Brown ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr); 1027552f7358SJed Brown } 1028552f7358SJed Brown ierr = PetscFree(fs_id);CHKERRQ(ierr); 1029552f7358SJed Brown } 1030552f7358SJed Brown #else 1031552f7358SJed Brown SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 1032552f7358SJed Brown #endif 1033552f7358SJed Brown PetscFunctionReturn(0); 1034552f7358SJed Brown } 1035