xref: /petsc/src/dm/impls/plex/plexexodusii.c (revision e02669251c099ae2faea4c9cee5397e8f99c2eb2)
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*e0266925SMatthew 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 }
67*e0266925SMatthew G. Knepley #endif
6878569bb4SMatthew G. Knepley 
6978569bb4SMatthew G. Knepley /*
7078569bb4SMatthew G. Knepley   DMPlexView_ExodusII_Internal - Write a DM to disk in exodus format
7178569bb4SMatthew G. Knepley 
7278569bb4SMatthew G. Knepley   Collective on dm
7378569bb4SMatthew G. Knepley 
7478569bb4SMatthew G. Knepley   Input Parameters:
7578569bb4SMatthew G. Knepley + dm  - The dm to be written
7678569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
7778569bb4SMatthew G. Knepley - degree - the degree of the interpolation space
7878569bb4SMatthew G. Knepley 
7978569bb4SMatthew G. Knepley   Notes:
8078569bb4SMatthew G. Knepley   Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels)
8178569bb4SMatthew G. Knepley   consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted.
8278569bb4SMatthew G. Knepley 
8378569bb4SMatthew 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
8478569bb4SMatthew G. Knepley   (including "ghost" cells and vertices) will be written. If all exoid points to the same file, the resulting file will
8578569bb4SMatthew G. Knepley   probably be corrupted.
8678569bb4SMatthew G. Knepley 
8778569bb4SMatthew G. Knepley   DMPlex only represents geometry while most post-processing software expect that a mesh also provides information
8878569bb4SMatthew G. Knepley   on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2.
8978569bb4SMatthew G. Knepley   It should be extended to use PetscFE objects.
9078569bb4SMatthew G. Knepley 
9178569bb4SMatthew G. Knepley   This function will only handle TRI, TET, QUAD and HEX cells.
9278569bb4SMatthew G. Knepley   Level: beginner
9378569bb4SMatthew G. Knepley 
9478569bb4SMatthew G. Knepley .keywords: mesh, ExodusII
9578569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadNodal_PlexEXO(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO()
9678569bb4SMatthew G. Knepley */
9778569bb4SMatthew G. Knepley PetscErrorCode DMPlexView_ExodusII_Internal(DM dm, int exoid, PetscInt degree)
9878569bb4SMatthew G. Knepley {
9978569bb4SMatthew G. Knepley   enum ElemType {TRI, QUAD, TET, HEX};
10078569bb4SMatthew G. Knepley   MPI_Comm        comm;
10178569bb4SMatthew G. Knepley   /* Connectivity Variables */
10278569bb4SMatthew G. Knepley   PetscInt        cellsNotInConnectivity;
10378569bb4SMatthew G. Knepley   /* Cell Sets */
10478569bb4SMatthew G. Knepley   DMLabel         csLabel;
10578569bb4SMatthew G. Knepley   IS              csIS;
10678569bb4SMatthew G. Knepley   const PetscInt *csIdx;
10778569bb4SMatthew G. Knepley   PetscInt        num_cs, cs;
10878569bb4SMatthew G. Knepley   enum ElemType  *type;
10978569bb4SMatthew G. Knepley   /* Coordinate Variables */
11078569bb4SMatthew G. Knepley   DM              cdm;
11178569bb4SMatthew G. Knepley   PetscSection    section;
11278569bb4SMatthew G. Knepley   Vec             coord;
11378569bb4SMatthew G. Knepley   PetscInt      **nodes;
114eae8a387SMatthew G. Knepley   PetscInt        depth, d, dim, skipCells = 0;
11578569bb4SMatthew G. Knepley   PetscInt        pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes;
11678569bb4SMatthew G. Knepley   PetscInt        num_vs, num_fs;
11778569bb4SMatthew G. Knepley   PetscMPIInt     rank, size;
11878569bb4SMatthew G. Knepley   const char     *dmName;
11978569bb4SMatthew G. Knepley   PetscErrorCode  ierr;
12078569bb4SMatthew G. Knepley 
12178569bb4SMatthew G. Knepley   PetscFunctionBegin;
12278569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
12378569bb4SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
12478569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
12578569bb4SMatthew G. Knepley   /* --- Get DM info --- */
12678569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) dm, &dmName);CHKERRQ(ierr);
12778569bb4SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
12878569bb4SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
12978569bb4SMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
13078569bb4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
13178569bb4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
13278569bb4SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);
13378569bb4SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
13478569bb4SMatthew G. Knepley   numCells    = cEnd - cStart;
13578569bb4SMatthew G. Knepley   numEdges    = eEnd - eStart;
13678569bb4SMatthew G. Knepley   numVertices = vEnd - vStart;
13778569bb4SMatthew G. Knepley   if (depth == 3) {numFaces = fEnd - fStart;}
13878569bb4SMatthew G. Knepley   else            {numFaces = 0;}
13978569bb4SMatthew G. Knepley   ierr = DMGetLabelSize(dm, "Cell Sets", &num_cs);CHKERRQ(ierr);
14078569bb4SMatthew G. Knepley   ierr = DMGetLabelSize(dm, "Vertex Sets", &num_vs);CHKERRQ(ierr);
14178569bb4SMatthew G. Knepley   ierr = DMGetLabelSize(dm, "Face Sets", &num_fs);CHKERRQ(ierr);
14278569bb4SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr);
14378569bb4SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
14478569bb4SMatthew G. Knepley   if (num_cs > 0) {
14578569bb4SMatthew G. Knepley     ierr = DMGetLabel(dm, "Cell Sets", &csLabel);CHKERRQ(ierr);
14678569bb4SMatthew G. Knepley     ierr = DMLabelGetValueIS(csLabel, &csIS);CHKERRQ(ierr);
14778569bb4SMatthew G. Knepley     ierr = ISGetIndices(csIS, &csIdx);CHKERRQ(ierr);
14878569bb4SMatthew G. Knepley   }
14978569bb4SMatthew G. Knepley   ierr = PetscMalloc1(num_cs, &nodes);CHKERRQ(ierr);
15078569bb4SMatthew G. Knepley   /* Set element type for each block and compute total number of nodes */
15178569bb4SMatthew G. Knepley   ierr = PetscMalloc1(num_cs, &type);CHKERRQ(ierr);
15278569bb4SMatthew G. Knepley   numNodes = numVertices;
15378569bb4SMatthew G. Knepley   if (degree == 2) {numNodes += numEdges;}
15478569bb4SMatthew G. Knepley   cellsNotInConnectivity = numCells;
15578569bb4SMatthew G. Knepley   for (cs = 0; cs < num_cs; ++cs) {
15678569bb4SMatthew G. Knepley     IS              stratumIS;
15778569bb4SMatthew G. Knepley     const PetscInt *cells;
15878569bb4SMatthew G. Knepley     PetscScalar    *xyz = NULL;
15978569bb4SMatthew G. Knepley     PetscInt        csSize, closureSize;
16078569bb4SMatthew G. Knepley     PetscInt        nodesTriP1[4]  = {3,0,0,0};
16178569bb4SMatthew G. Knepley     PetscInt        nodesTriP2[4]  = {3,3,0,0};
16278569bb4SMatthew G. Knepley     PetscInt        nodesQuadP1[4] = {4,0,0,0};
16378569bb4SMatthew G. Knepley     PetscInt        nodesQuadP2[4] = {4,4,0,1};
16478569bb4SMatthew G. Knepley     PetscInt        nodesTetP1[4]  = {4,0,0,0};
16578569bb4SMatthew G. Knepley     PetscInt        nodesTetP2[4]  = {4,6,0,0};
16678569bb4SMatthew G. Knepley     PetscInt        nodesHexP1[4]  = {8,0,0,0};
16778569bb4SMatthew G. Knepley     PetscInt        nodesHexP2[4]  = {8,12,6,1};
16878569bb4SMatthew G. Knepley 
16978569bb4SMatthew G. Knepley     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
17078569bb4SMatthew G. Knepley     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
17178569bb4SMatthew G. Knepley     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
17278569bb4SMatthew G. Knepley     ierr = DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr);
17378569bb4SMatthew G. Knepley     switch (dim) {
17478569bb4SMatthew G. Knepley     case 2:
17578569bb4SMatthew G. Knepley       if      (closureSize == 3*dim) {type[cs] = TRI;}
17678569bb4SMatthew G. Knepley       else if (closureSize == 4*dim) {type[cs] = QUAD;}
17778569bb4SMatthew 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);
17878569bb4SMatthew G. Knepley       break;
17978569bb4SMatthew G. Knepley     case 3:
18078569bb4SMatthew G. Knepley       if      (closureSize == 4*dim) {type[cs] = TET;}
18178569bb4SMatthew G. Knepley       else if (closureSize == 8*dim) {type[cs] = HEX;}
18278569bb4SMatthew 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);
18378569bb4SMatthew G. Knepley       break;
18478569bb4SMatthew G. Knepley     default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not handled by ExodusII viewer", dim);
18578569bb4SMatthew G. Knepley     }
18678569bb4SMatthew G. Knepley     if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;}
18778569bb4SMatthew G. Knepley     if ((degree == 2) && (type[cs] == HEX))  {numNodes += csSize; numNodes += numFaces;}
18878569bb4SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr);
18978569bb4SMatthew G. Knepley     /* Set nodes and Element type */
19078569bb4SMatthew G. Knepley     if (type[cs] == TRI) {
19178569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesTriP1;
19278569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesTriP2;
19378569bb4SMatthew G. Knepley     } else if (type[cs] == QUAD) {
19478569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesQuadP1;
19578569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesQuadP2;
19678569bb4SMatthew G. Knepley     } else if (type[cs] == TET) {
19778569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesTetP1;
19878569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesTetP2;
19978569bb4SMatthew G. Knepley     } else if (type[cs] == HEX) {
20078569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesHexP1;
20178569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesHexP2;
20278569bb4SMatthew G. Knepley     }
20378569bb4SMatthew G. Knepley     /* Compute the number of cells not in the connectivity table */
20478569bb4SMatthew G. Knepley     cellsNotInConnectivity -= nodes[cs][3]*csSize;
20578569bb4SMatthew G. Knepley 
20678569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
20778569bb4SMatthew G. Knepley     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
20878569bb4SMatthew G. Knepley   }
20978569bb4SMatthew G. Knepley   if (num_cs > 0) {ierr = ex_put_init(exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs);CHKERRQ(ierr);}
21078569bb4SMatthew G. Knepley   /* --- Connectivity --- */
21178569bb4SMatthew G. Knepley   for (cs = 0; cs < num_cs; ++cs) {
21278569bb4SMatthew G. Knepley     IS              stratumIS;
21378569bb4SMatthew G. Knepley     const PetscInt *cells;
21478569bb4SMatthew G. Knepley     PetscInt       *connect;
215eae8a387SMatthew G. Knepley     PetscInt        edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0;
21678569bb4SMatthew G. Knepley     PetscInt        csSize, c, connectSize, closureSize;
21778569bb4SMatthew G. Knepley     char           *elem_type = NULL;
21878569bb4SMatthew G. Knepley     char            elem_type_tri3[]  = "TRI3",  elem_type_quad4[] = "QUAD4";
21978569bb4SMatthew G. Knepley     char            elem_type_tri6[]  = "TRI6",  elem_type_quad9[] = "QUAD9";
22078569bb4SMatthew G. Knepley     char            elem_type_tet4[]  = "TET4",  elem_type_hex8[]  = "HEX8";
22178569bb4SMatthew G. Knepley     char            elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27";
22278569bb4SMatthew G. Knepley 
22378569bb4SMatthew G. Knepley     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
22478569bb4SMatthew G. Knepley     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
22578569bb4SMatthew G. Knepley     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
22678569bb4SMatthew G. Knepley     /* Set Element type */
22778569bb4SMatthew G. Knepley     if (type[cs] == TRI) {
22878569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_tri3;
22978569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_tri6;
23078569bb4SMatthew G. Knepley     } else if (type[cs] == QUAD) {
23178569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_quad4;
23278569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_quad9;
23378569bb4SMatthew G. Knepley     } else if (type[cs] == TET) {
23478569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_tet4;
23578569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_tet10;
23678569bb4SMatthew G. Knepley     } else if (type[cs] == HEX) {
23778569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_hex8;
23878569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_hex27;
23978569bb4SMatthew G. Knepley     }
24078569bb4SMatthew G. Knepley     connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3];
24178569bb4SMatthew G. Knepley     ierr = PetscMalloc1(connectSize, &connect);CHKERRQ(ierr);
24278569bb4SMatthew G. Knepley     ierr = ex_put_block(exoid, EX_ELEM_BLOCK, cs, elem_type, csSize, connectSize, 0, 0, 1);
24378569bb4SMatthew G. Knepley     /* Find number of vertices, edges, and faces in the closure */
24478569bb4SMatthew G. Knepley     verticesInClosure = nodes[cs][0];
24578569bb4SMatthew G. Knepley     if (depth > 1) {
24678569bb4SMatthew G. Knepley       if (dim == 2) {
24778569bb4SMatthew G. Knepley         ierr = DMPlexGetConeSize(dm, cells[0], &edgesInClosure);CHKERRQ(ierr);
24878569bb4SMatthew G. Knepley       } else if (dim == 3) {
24978569bb4SMatthew G. Knepley         PetscInt *closure = NULL;
25078569bb4SMatthew G. Knepley 
25178569bb4SMatthew G. Knepley         ierr = DMPlexGetConeSize(dm, cells[0], &facesInClosure);CHKERRQ(ierr);
25278569bb4SMatthew G. Knepley         ierr = DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
25378569bb4SMatthew G. Knepley         edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure;
25478569bb4SMatthew G. Knepley         ierr = DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
25578569bb4SMatthew G. Knepley       }
25678569bb4SMatthew G. Knepley     }
25778569bb4SMatthew G. Knepley     /* Get connectivity for each cell */
25878569bb4SMatthew G. Knepley     for (c = 0; c < csSize; ++c) {
25978569bb4SMatthew G. Knepley       PetscInt *closure = NULL;
26078569bb4SMatthew G. Knepley       PetscInt  temp, i;
26178569bb4SMatthew G. Knepley 
26278569bb4SMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
26378569bb4SMatthew G. Knepley       for (i = 0; i < connectSize; ++i) {
26478569bb4SMatthew G. Knepley         if (i < nodes[cs][0]) {/* Vertices */
26578569bb4SMatthew G. Knepley           connect[i] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1;
26678569bb4SMatthew G. Knepley           connect[i] -= cellsNotInConnectivity;
26778569bb4SMatthew G. Knepley         } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */
26878569bb4SMatthew G. Knepley           connect[i] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1;
26978569bb4SMatthew G. Knepley           if (nodes[cs][2] == 0) connect[i] -= numFaces;
27078569bb4SMatthew G. Knepley           connect[i] -= cellsNotInConnectivity;
27178569bb4SMatthew G. Knepley         } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]){ /* Cells */
27278569bb4SMatthew G. Knepley           connect[i] = closure[0] + 1;
27378569bb4SMatthew G. Knepley           connect[i] -= skipCells;
27478569bb4SMatthew G. Knepley         } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]){ /* Faces */
27578569bb4SMatthew G. Knepley           connect[i] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1;
27678569bb4SMatthew G. Knepley           connect[i] -= cellsNotInConnectivity;
27778569bb4SMatthew G. Knepley         } else {
27878569bb4SMatthew G. Knepley           connect[i] = -1;
27978569bb4SMatthew G. Knepley         }
28078569bb4SMatthew G. Knepley       }
28178569bb4SMatthew G. Knepley       /* Tetrahedra are inverted */
28278569bb4SMatthew G. Knepley       if (type[cs] == TET) {
28378569bb4SMatthew G. Knepley         temp = connect[0]; connect[0] = connect[1]; connect[1] = temp;
28478569bb4SMatthew G. Knepley         if (degree == 2) {
28578569bb4SMatthew G. Knepley           temp = connect[5]; connect[5] = connect[6]; connect[6] = temp;
28678569bb4SMatthew G. Knepley           temp = connect[7]; connect[7] = connect[8]; connect[8] = temp;
28778569bb4SMatthew G. Knepley         }
28878569bb4SMatthew G. Knepley       }
28978569bb4SMatthew G. Knepley       /* Hexahedra are inverted */
29078569bb4SMatthew G. Knepley       if (type[cs] == HEX) {
29178569bb4SMatthew G. Knepley         temp = connect[1]; connect[1] = connect[3]; connect[3] = temp;
29278569bb4SMatthew G. Knepley         if (degree == 2) {
29378569bb4SMatthew G. Knepley           temp = connect[8];  connect[8]  = connect[11]; connect[11] = temp;
29478569bb4SMatthew G. Knepley           temp = connect[9];  connect[9]  = connect[10]; connect[10] = temp;
29578569bb4SMatthew G. Knepley           temp = connect[16]; connect[16] = connect[17]; connect[17] = temp;
29678569bb4SMatthew G. Knepley           temp = connect[18]; connect[18] = connect[19]; connect[19] = temp;
29778569bb4SMatthew G. Knepley 
29878569bb4SMatthew G. Knepley           temp = connect[12]; connect[12] = connect[16]; connect[16] = temp;
29978569bb4SMatthew G. Knepley           temp = connect[13]; connect[13] = connect[17]; connect[17] = temp;
30078569bb4SMatthew G. Knepley           temp = connect[14]; connect[14] = connect[18]; connect[18] = temp;
30178569bb4SMatthew G. Knepley           temp = connect[15]; connect[15] = connect[19]; connect[19] = temp;
30278569bb4SMatthew G. Knepley 
30378569bb4SMatthew G. Knepley           temp = connect[23]; connect[23] = connect[26]; connect[26] = temp;
30478569bb4SMatthew G. Knepley           temp = connect[24]; connect[24] = connect[25]; connect[25] = temp;
30578569bb4SMatthew G. Knepley           temp = connect[25]; connect[25] = connect[26]; connect[26] = temp;
30678569bb4SMatthew G. Knepley         }
30778569bb4SMatthew G. Knepley       }
30878569bb4SMatthew G. Knepley       ierr = ex_put_partial_conn(exoid, EX_ELEM_BLOCK, cs, c+1, 1, connect, NULL, NULL);
30978569bb4SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
31078569bb4SMatthew G. Knepley     }
31178569bb4SMatthew G. Knepley     skipCells += (nodes[cs][3] == 0)*csSize;
31278569bb4SMatthew G. Knepley     ierr = PetscFree(connect);CHKERRQ(ierr);
31378569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
31478569bb4SMatthew G. Knepley     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
31578569bb4SMatthew G. Knepley   }
31678569bb4SMatthew G. Knepley   ierr = PetscFree(type);CHKERRQ(ierr);
31778569bb4SMatthew G. Knepley   /* --- Coordinates --- */
31878569bb4SMatthew G. Knepley   ierr = PetscSectionCreate(comm, &section);CHKERRQ(ierr);
31978569bb4SMatthew G. Knepley   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
32078569bb4SMatthew G. Knepley   for (d = 0; d < depth; ++d) {
32178569bb4SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
32278569bb4SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
32378569bb4SMatthew G. Knepley       ierr = PetscSectionSetDof(section, p, nodes[0][d] > 0);CHKERRQ(ierr);
32478569bb4SMatthew G. Knepley     }
32578569bb4SMatthew G. Knepley   }
32678569bb4SMatthew G. Knepley   for (cs = 0; cs < num_cs; ++cs) {
32778569bb4SMatthew G. Knepley     IS              stratumIS;
32878569bb4SMatthew G. Knepley     const PetscInt *cells;
32978569bb4SMatthew G. Knepley     PetscInt        csSize, c;
33078569bb4SMatthew G. Knepley 
33178569bb4SMatthew G. Knepley     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
33278569bb4SMatthew G. Knepley     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
33378569bb4SMatthew G. Knepley     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
33478569bb4SMatthew G. Knepley     for (c = 0; c < csSize; ++c) {
33578569bb4SMatthew G. Knepley       ierr = PetscSectionSetDof(section, cells[c], nodes[cs][3] > 0);CHKERRQ(ierr);
33678569bb4SMatthew G. Knepley     }
33778569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
33878569bb4SMatthew G. Knepley     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
33978569bb4SMatthew G. Knepley   }
34078569bb4SMatthew G. Knepley   if (num_cs > 0) {
34178569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(csIS, &csIdx);CHKERRQ(ierr);
34278569bb4SMatthew G. Knepley     ierr = ISDestroy(&csIS);CHKERRQ(ierr);
34378569bb4SMatthew G. Knepley   }
34478569bb4SMatthew G. Knepley   ierr = PetscFree(nodes);CHKERRQ(ierr);
34578569bb4SMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
34678569bb4SMatthew G. Knepley   if (numNodes > 0) {
34778569bb4SMatthew G. Knepley     const char  *coordNames[3] = {"x", "y", "z"};
34878569bb4SMatthew G. Knepley     PetscScalar *coords, *closure;
34978569bb4SMatthew G. Knepley     PetscReal   *cval;
35078569bb4SMatthew G. Knepley     PetscInt     hasDof, n = 0;
35178569bb4SMatthew G. Knepley 
35278569bb4SMatthew G. Knepley     /* There can't be more than 24 values in the closure of a point for the coord section */
35378569bb4SMatthew G. Knepley     ierr = PetscMalloc3(numNodes*3, &coords, dim, &cval, 24, &closure);CHKERRQ(ierr);
35478569bb4SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr);
35578569bb4SMatthew G. Knepley     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
35678569bb4SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
35778569bb4SMatthew G. Knepley       ierr = PetscSectionGetDof(section, p, &hasDof);
35878569bb4SMatthew G. Knepley       if (hasDof) {
35978569bb4SMatthew G. Knepley         PetscInt closureSize = 24, j;
36078569bb4SMatthew G. Knepley 
36178569bb4SMatthew G. Knepley         ierr = DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure);CHKERRQ(ierr);
36278569bb4SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
36378569bb4SMatthew G. Knepley           cval[d] = 0.0;
36478569bb4SMatthew G. Knepley           for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d];
36578569bb4SMatthew G. Knepley           coords[d*numNodes+n] = cval[d] * dim / closureSize;
36678569bb4SMatthew G. Knepley         }
36778569bb4SMatthew G. Knepley         ++n;
36878569bb4SMatthew G. Knepley       }
36978569bb4SMatthew G. Knepley     }
37078569bb4SMatthew G. Knepley     ierr = ex_put_coord(exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes]);CHKERRQ(ierr);
37178569bb4SMatthew G. Knepley     ierr = PetscFree3(coords, cval, closure);CHKERRQ(ierr);
37278569bb4SMatthew G. Knepley     ierr = ex_put_coord_names(exoid, (char **) coordNames);CHKERRQ(ierr);
37378569bb4SMatthew G. Knepley   }
37478569bb4SMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
37578569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
37678569bb4SMatthew G. Knepley }
37778569bb4SMatthew G. Knepley 
37878569bb4SMatthew G. Knepley /*
37978569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file
38078569bb4SMatthew G. Knepley 
38178569bb4SMatthew G. Knepley   Collective on v
38278569bb4SMatthew G. Knepley 
38378569bb4SMatthew G. Knepley   Input Parameters:
38478569bb4SMatthew G. Knepley + v  - The vector to be written
38578569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
38678569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1)
38778569bb4SMatthew G. Knepley 
38878569bb4SMatthew G. Knepley   Notes:
38978569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
39078569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
39178569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
39278569bb4SMatthew G. Knepley   amongst all nodal variables.
39378569bb4SMatthew G. Knepley 
39478569bb4SMatthew G. Knepley   Level: beginner
39578569bb4SMatthew G. Knepley 
39678569bb4SMatthew G. Knepley .keywords: mesh,ExodusII
39778569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecLoadNodal_PlexEXO(),VecViewZonal_PlexEXO(),VecLoadZonal_PlexEXO()
39878569bb4SMatthew G. Knepley @*/
39978569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step)
40078569bb4SMatthew G. Knepley {
40178569bb4SMatthew G. Knepley   MPI_Comm         comm;
40278569bb4SMatthew G. Knepley   PetscMPIInt      size;
40378569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
40478569bb4SMatthew G. Knepley   DM               dm;
40578569bb4SMatthew G. Knepley   Vec              vNatural, vComp;
40678569bb4SMatthew G. Knepley   const PetscReal *varray;
40778569bb4SMatthew G. Knepley   const char      *vecname;
40878569bb4SMatthew G. Knepley   PetscInt         xs, xe, bs;
40978569bb4SMatthew G. Knepley   PetscBool        useNatural;
41078569bb4SMatthew G. Knepley   int              offset;
41178569bb4SMatthew G. Knepley   PetscErrorCode   ierr;
41278569bb4SMatthew G. Knepley #endif
41378569bb4SMatthew G. Knepley 
41478569bb4SMatthew G. Knepley   PetscFunctionBegin;
41578569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
41678569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
41778569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
41878569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
41978569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
42078569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
42178569bb4SMatthew G. Knepley   if (useNatural) {
42278569bb4SMatthew G. Knepley     ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr);
42378569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr);
42478569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr);
42578569bb4SMatthew G. Knepley   } else {
42678569bb4SMatthew G. Knepley     vNatural = v;
42778569bb4SMatthew G. Knepley   }
42878569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
42978569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
43078569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr);
43178569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file\n", vecname);
43278569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
43378569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
43478569bb4SMatthew G. Knepley      We assume that they are stored sequentially */
43578569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
43678569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
43778569bb4SMatthew G. Knepley   if (bs == 1) {
43878569bb4SMatthew G. Knepley     ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr);
43978569bb4SMatthew G. Knepley     ierr = ex_put_partial_var(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);CHKERRQ(ierr);
44078569bb4SMatthew G. Knepley     ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr);
44178569bb4SMatthew G. Knepley   } else {
44278569bb4SMatthew G. Knepley     IS       compIS;
44378569bb4SMatthew G. Knepley     PetscInt c;
44478569bb4SMatthew G. Knepley 
44578569bb4SMatthew G. Knepley     ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);
44678569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
44778569bb4SMatthew G. Knepley       ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
44878569bb4SMatthew G. Knepley       ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
44978569bb4SMatthew G. Knepley       ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr);
45078569bb4SMatthew G. Knepley       ierr = ex_put_partial_var(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);CHKERRQ(ierr);
45178569bb4SMatthew G. Knepley       ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr);
45278569bb4SMatthew G. Knepley       ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
45378569bb4SMatthew G. Knepley     }
45478569bb4SMatthew G. Knepley     ierr = ISDestroy(&compIS);CHKERRQ(ierr);
45578569bb4SMatthew G. Knepley   }
45678569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);}
45778569bb4SMatthew G. Knepley #else
45878569bb4SMatthew G. Knepley   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
45978569bb4SMatthew G. Knepley #endif
46078569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
46178569bb4SMatthew G. Knepley }
46278569bb4SMatthew G. Knepley 
46378569bb4SMatthew G. Knepley /*
46478569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file
46578569bb4SMatthew G. Knepley 
46678569bb4SMatthew G. Knepley   Collective on v
46778569bb4SMatthew G. Knepley 
46878569bb4SMatthew G. Knepley   Input Parameters:
46978569bb4SMatthew G. Knepley + v  - The vector to be written
47078569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
47178569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1)
47278569bb4SMatthew G. Knepley 
47378569bb4SMatthew G. Knepley   Notes:
47478569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
47578569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
47678569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
47778569bb4SMatthew G. Knepley   amongst all nodal variables.
47878569bb4SMatthew G. Knepley 
47978569bb4SMatthew G. Knepley   Level: beginner
48078569bb4SMatthew G. Knepley 
48178569bb4SMatthew G. Knepley .keywords: mesh, ExodusII
48278569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO()
48378569bb4SMatthew G. Knepley */
48478569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step)
48578569bb4SMatthew G. Knepley {
48678569bb4SMatthew G. Knepley   MPI_Comm       comm;
48778569bb4SMatthew G. Knepley   PetscMPIInt    size;
48878569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
48978569bb4SMatthew G. Knepley   DM             dm;
49078569bb4SMatthew G. Knepley   Vec            vNatural, vComp;
49178569bb4SMatthew G. Knepley   PetscReal     *varray;
49278569bb4SMatthew G. Knepley   const char    *vecname;
49378569bb4SMatthew G. Knepley   PetscInt       xs, xe, bs;
49478569bb4SMatthew G. Knepley   PetscBool      useNatural;
49578569bb4SMatthew G. Knepley   int            offset;
49678569bb4SMatthew G. Knepley   PetscErrorCode ierr;
49778569bb4SMatthew G. Knepley #endif
49878569bb4SMatthew G. Knepley 
49978569bb4SMatthew G. Knepley   PetscFunctionBegin;
50078569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
50178569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
50278569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
50378569bb4SMatthew G. Knepley   ierr = VecGetDM(v,&dm);CHKERRQ(ierr);
50478569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
50578569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
50678569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
50778569bb4SMatthew G. Knepley   else            {vNatural = v;}
50878569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
50978569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
51078569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr);
51178569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file", vecname);
51278569bb4SMatthew G. Knepley   /* Read local chunk from the file */
51378569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
51478569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
51578569bb4SMatthew G. Knepley   if (bs == 1) {
51678569bb4SMatthew G. Knepley     ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr);
51778569bb4SMatthew G. Knepley     ierr = ex_get_partial_var(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);CHKERRQ(ierr);
51878569bb4SMatthew G. Knepley     ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr);
51978569bb4SMatthew G. Knepley   } else {
52078569bb4SMatthew G. Knepley     IS       compIS;
52178569bb4SMatthew G. Knepley     PetscInt c;
52278569bb4SMatthew G. Knepley 
52378569bb4SMatthew G. Knepley     ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);
52478569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
52578569bb4SMatthew G. Knepley       ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
52678569bb4SMatthew G. Knepley       ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
52778569bb4SMatthew G. Knepley       ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr);
52878569bb4SMatthew G. Knepley       ierr = ex_get_partial_var(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);CHKERRQ(ierr);
52978569bb4SMatthew G. Knepley       ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr);
53078569bb4SMatthew G. Knepley       ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
53178569bb4SMatthew G. Knepley     }
53278569bb4SMatthew G. Knepley     ierr = ISDestroy(&compIS);CHKERRQ(ierr);
53378569bb4SMatthew G. Knepley   }
53478569bb4SMatthew G. Knepley   if (useNatural) {
53578569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr);
53678569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr);
53778569bb4SMatthew G. Knepley     ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);
53878569bb4SMatthew G. Knepley   }
53978569bb4SMatthew G. Knepley #else
54078569bb4SMatthew G. Knepley   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
54178569bb4SMatthew G. Knepley #endif
54278569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
54378569bb4SMatthew G. Knepley }
54478569bb4SMatthew G. Knepley 
54578569bb4SMatthew G. Knepley /*
54678569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file
54778569bb4SMatthew G. Knepley 
54878569bb4SMatthew G. Knepley   Collective on v
54978569bb4SMatthew G. Knepley 
55078569bb4SMatthew G. Knepley   Input Parameters:
55178569bb4SMatthew G. Knepley + v  - The vector to be written
55278569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
55378569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1)
55478569bb4SMatthew G. Knepley 
55578569bb4SMatthew G. Knepley   Notes:
55678569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
55778569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
55878569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
55978569bb4SMatthew G. Knepley   amongst all zonal variables.
56078569bb4SMatthew G. Knepley 
56178569bb4SMatthew G. Knepley   Level: beginner
56278569bb4SMatthew G. Knepley 
56378569bb4SMatthew G. Knepley .keywords: mesh,ExodusII
56478569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Zonal_Internal()
56578569bb4SMatthew G. Knepley */
56678569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step)
56778569bb4SMatthew G. Knepley {
56878569bb4SMatthew G. Knepley   MPI_Comm          comm;
56978569bb4SMatthew G. Knepley   PetscMPIInt       size;
57078569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
57178569bb4SMatthew G. Knepley   DM                dm;
57278569bb4SMatthew G. Knepley   Vec               vNatural, vComp;
57378569bb4SMatthew G. Knepley   const PetscReal  *varray;
57478569bb4SMatthew G. Knepley   const char       *vecname;
57578569bb4SMatthew G. Knepley   PetscInt          xs, xe, bs;
57678569bb4SMatthew G. Knepley   PetscBool         useNatural;
57778569bb4SMatthew G. Knepley   int               offset;
57878569bb4SMatthew G. Knepley   IS                compIS;
57978569bb4SMatthew G. Knepley   PetscInt         *csSize, *csID;
58078569bb4SMatthew G. Knepley   PetscInt          numCS, set, csxs = 0;
58178569bb4SMatthew G. Knepley   PetscErrorCode    ierr;
58278569bb4SMatthew G. Knepley #endif
58378569bb4SMatthew G. Knepley 
58478569bb4SMatthew G. Knepley   PetscFunctionBegin;
58578569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)v, &comm);CHKERRQ(ierr);
58678569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
58778569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
58878569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
58978569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
59078569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
59178569bb4SMatthew G. Knepley   if (useNatural) {
59278569bb4SMatthew G. Knepley     ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr);
59378569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr);
59478569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr);
59578569bb4SMatthew G. Knepley   } else {
59678569bb4SMatthew G. Knepley     vNatural = v;
59778569bb4SMatthew G. Knepley   }
59878569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
59978569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
60078569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr);
60178569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname);
60278569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
60378569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
60478569bb4SMatthew G. Knepley      We assume that they are stored sequentially
60578569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
60678569bb4SMatthew G. Knepley      but once the vector has been reordered to natural size, we cannot use the label informations
60778569bb4SMatthew G. Knepley      to figure out what to save where. */
60878569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
60978569bb4SMatthew G. Knepley   ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr);
61078569bb4SMatthew G. Knepley   ierr = ex_get_ids(exoid, EX_ELEM_BLOCK, csID);CHKERRQ(ierr);
61178569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
61278569bb4SMatthew G. Knepley     ex_block block;
61378569bb4SMatthew G. Knepley 
61478569bb4SMatthew G. Knepley     block.id   = csID[set];
61578569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
61678569bb4SMatthew G. Knepley     ierr = ex_get_block_param(exoid, &block);CHKERRQ(ierr);
61778569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
61878569bb4SMatthew G. Knepley   }
61978569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
62078569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
62178569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);}
62278569bb4SMatthew G. Knepley   for (set = 0; set < numCS; set++) {
62378569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
62478569bb4SMatthew G. Knepley 
62578569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
62678569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
62778569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
62878569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
62978569bb4SMatthew G. Knepley     if (bs == 1) {
63078569bb4SMatthew G. Knepley       ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr);
63178569bb4SMatthew 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);
63278569bb4SMatthew G. Knepley       ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr);
63378569bb4SMatthew G. Knepley     } else {
63478569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
63578569bb4SMatthew G. Knepley         ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
63678569bb4SMatthew G. Knepley         ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
63778569bb4SMatthew G. Knepley         ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr);
63878569bb4SMatthew 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)]);
63978569bb4SMatthew G. Knepley         ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr);
64078569bb4SMatthew G. Knepley         ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
64178569bb4SMatthew G. Knepley       }
64278569bb4SMatthew G. Knepley     }
64378569bb4SMatthew G. Knepley     csxs += csSize[set];
64478569bb4SMatthew G. Knepley   }
64578569bb4SMatthew G. Knepley   ierr = PetscFree2(csID, csSize);CHKERRQ(ierr);
64678569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);}
64778569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMRestoreGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
64878569bb4SMatthew G. Knepley #else
64978569bb4SMatthew G. Knepley   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
65078569bb4SMatthew G. Knepley #endif
65178569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
65278569bb4SMatthew G. Knepley }
65378569bb4SMatthew G. Knepley 
65478569bb4SMatthew G. Knepley /*
65578569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file
65678569bb4SMatthew G. Knepley 
65778569bb4SMatthew G. Knepley   Collective on v
65878569bb4SMatthew G. Knepley 
65978569bb4SMatthew G. Knepley   Input Parameters:
66078569bb4SMatthew G. Knepley + v  - The vector to be written
66178569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
66278569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1)
66378569bb4SMatthew G. Knepley 
66478569bb4SMatthew G. Knepley   Notes:
66578569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
66678569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
66778569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
66878569bb4SMatthew G. Knepley   amongst all zonal variables.
66978569bb4SMatthew G. Knepley 
67078569bb4SMatthew G. Knepley   Level: beginner
67178569bb4SMatthew G. Knepley 
67278569bb4SMatthew G. Knepley .keywords: mesh,ExodusII
67378569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Zonal_Internal()
67478569bb4SMatthew G. Knepley */
67578569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step)
67678569bb4SMatthew G. Knepley {
67778569bb4SMatthew G. Knepley   MPI_Comm          comm;
67878569bb4SMatthew G. Knepley   PetscMPIInt       size;
67978569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
68078569bb4SMatthew G. Knepley   DM                dm;
68178569bb4SMatthew G. Knepley   Vec               vNatural, vComp;
68278569bb4SMatthew G. Knepley   PetscReal        *varray;
68378569bb4SMatthew G. Knepley   const char       *vecname;
68478569bb4SMatthew G. Knepley   PetscInt          xs, xe, bs;
68578569bb4SMatthew G. Knepley   PetscBool         useNatural;
68678569bb4SMatthew G. Knepley   int               offset;
68778569bb4SMatthew G. Knepley   IS                compIS;
68878569bb4SMatthew G. Knepley   PetscInt         *csSize, *csID;
68978569bb4SMatthew G. Knepley   PetscInt          numCS, set, csxs = 0;
69078569bb4SMatthew G. Knepley   PetscErrorCode    ierr;
69178569bb4SMatthew G. Knepley #endif
69278569bb4SMatthew G. Knepley 
69378569bb4SMatthew G. Knepley   PetscFunctionBegin;
69478569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
69578569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
69678569bb4SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
69778569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
69878569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
69978569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
70078569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
70178569bb4SMatthew G. Knepley   else            {vNatural = v;}
70278569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
70378569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
70478569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr);
70578569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname);
70678569bb4SMatthew G. Knepley   /* Read local chunk of the result in the exodus file
70778569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
70878569bb4SMatthew G. Knepley      We assume that they are stored sequentially
70978569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
71078569bb4SMatthew G. Knepley      but once the vector has been reordered to natural size, we cannot use the label informations
71178569bb4SMatthew G. Knepley      to figure out what to save where. */
71278569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
71378569bb4SMatthew G. Knepley   ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr);
71478569bb4SMatthew G. Knepley   ierr = ex_get_ids(exoid, EX_ELEM_BLOCK, csID);CHKERRQ(ierr);
71578569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
71678569bb4SMatthew G. Knepley     ex_block block;
71778569bb4SMatthew G. Knepley 
71878569bb4SMatthew G. Knepley     block.id   = csID[set];
71978569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
72078569bb4SMatthew G. Knepley     ierr = ex_get_block_param(exoid, &block);CHKERRQ(ierr);
72178569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
72278569bb4SMatthew G. Knepley   }
72378569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
72478569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
72578569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);}
72678569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
72778569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
72878569bb4SMatthew G. Knepley 
72978569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
73078569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
73178569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
73278569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
73378569bb4SMatthew G. Knepley     if (bs == 1) {
73478569bb4SMatthew G. Knepley       ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr);
73578569bb4SMatthew 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);
73678569bb4SMatthew G. Knepley       ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr);
73778569bb4SMatthew G. Knepley     } else {
73878569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
73978569bb4SMatthew G. Knepley         ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
74078569bb4SMatthew G. Knepley         ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
74178569bb4SMatthew G. Knepley         ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr);
74278569bb4SMatthew 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)]);
74378569bb4SMatthew G. Knepley         ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr);
74478569bb4SMatthew G. Knepley         ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
74578569bb4SMatthew G. Knepley       }
74678569bb4SMatthew G. Knepley     }
74778569bb4SMatthew G. Knepley     csxs += csSize[set];
74878569bb4SMatthew G. Knepley   }
74978569bb4SMatthew G. Knepley   ierr = PetscFree2(csID, csSize);CHKERRQ(ierr);
75078569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);}
75178569bb4SMatthew G. Knepley   if (useNatural) {
75278569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr);
75378569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr);
75478569bb4SMatthew G. Knepley     ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);
75578569bb4SMatthew G. Knepley   }
75678569bb4SMatthew G. Knepley #else
75778569bb4SMatthew G. Knepley   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
75878569bb4SMatthew G. Knepley #endif
75978569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
76078569bb4SMatthew G. Knepley }
76178569bb4SMatthew G. Knepley 
76233751fbdSMichael Lange /*@C
763eaf898f9SPatrick Sanan   DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file.
76433751fbdSMichael Lange 
76533751fbdSMichael Lange   Collective on comm
76633751fbdSMichael Lange 
76733751fbdSMichael Lange   Input Parameters:
76833751fbdSMichael Lange + comm  - The MPI communicator
76933751fbdSMichael Lange . filename - The name of the ExodusII file
77033751fbdSMichael Lange - interpolate - Create faces and edges in the mesh
77133751fbdSMichael Lange 
77233751fbdSMichael Lange   Output Parameter:
77333751fbdSMichael Lange . dm  - The DM object representing the mesh
77433751fbdSMichael Lange 
77533751fbdSMichael Lange   Level: beginner
77633751fbdSMichael Lange 
77733751fbdSMichael Lange .keywords: mesh,ExodusII
77833751fbdSMichael Lange .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus()
77933751fbdSMichael Lange @*/
78033751fbdSMichael Lange PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
78133751fbdSMichael Lange {
78233751fbdSMichael Lange   PetscMPIInt    rank;
78333751fbdSMichael Lange   PetscErrorCode ierr;
78433751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
785e84f7031SBlaise Bourdin   int   CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1;
78633751fbdSMichael Lange   float version;
78733751fbdSMichael Lange #endif
78833751fbdSMichael Lange 
78933751fbdSMichael Lange   PetscFunctionBegin;
79033751fbdSMichael Lange   PetscValidCharPointer(filename, 2);
79133751fbdSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
79233751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
79333751fbdSMichael Lange   if (!rank) {
79433751fbdSMichael Lange     exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version);
79533751fbdSMichael Lange     if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename);
79633751fbdSMichael Lange   }
79733751fbdSMichael Lange   ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr);
79833751fbdSMichael Lange   if (!rank) {ierr = ex_close(exoid);CHKERRQ(ierr);}
79933751fbdSMichael Lange #else
80033751fbdSMichael Lange   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
80133751fbdSMichael Lange #endif
80233751fbdSMichael Lange   PetscFunctionReturn(0);
80333751fbdSMichael Lange }
80433751fbdSMichael Lange 
805552f7358SJed Brown /*@
80633751fbdSMichael Lange   DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID.
807552f7358SJed Brown 
808552f7358SJed Brown   Collective on comm
809552f7358SJed Brown 
810552f7358SJed Brown   Input Parameters:
811552f7358SJed Brown + comm  - The MPI communicator
812552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
813552f7358SJed Brown - interpolate - Create faces and edges in the mesh
814552f7358SJed Brown 
815552f7358SJed Brown   Output Parameter:
816552f7358SJed Brown . dm  - The DM object representing the mesh
817552f7358SJed Brown 
818552f7358SJed Brown   Level: beginner
819552f7358SJed Brown 
820552f7358SJed Brown .keywords: mesh,ExodusII
821e84e7769SJed Brown .seealso: DMPLEX, DMCreate()
822552f7358SJed Brown @*/
823552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm)
824552f7358SJed Brown {
825552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
826552f7358SJed Brown   PetscMPIInt    num_proc, rank;
827552f7358SJed Brown   PetscSection   coordSection;
828552f7358SJed Brown   Vec            coordinates;
829552f7358SJed Brown   PetscScalar    *coords;
830552f7358SJed Brown   PetscInt       coordSize, v;
831552f7358SJed Brown   PetscErrorCode ierr;
832552f7358SJed Brown   /* Read from ex_get_init() */
833552f7358SJed Brown   char title[PETSC_MAX_PATH_LEN+1];
834552f7358SJed Brown   int  dim    = 0, numVertices = 0, numCells = 0;
835552f7358SJed Brown   int  num_cs = 0, num_vs = 0, num_fs = 0;
836552f7358SJed Brown #endif
837552f7358SJed Brown 
838552f7358SJed Brown   PetscFunctionBegin;
839552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
840552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
841552f7358SJed Brown   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
842552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
843552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
844552f7358SJed Brown   /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */
845552f7358SJed Brown   if (!rank) {
84639bba695SBarry Smith     ierr = PetscMemzero(title,(PETSC_MAX_PATH_LEN+1)*sizeof(char));CHKERRQ(ierr);
847552f7358SJed Brown     ierr = ex_get_init(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs);
848552f7358SJed Brown     if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ExodusII ex_get_init() failed with error code %D\n",ierr);
849552f7358SJed Brown     if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n");
850552f7358SJed Brown   }
851552f7358SJed Brown   ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
852552f7358SJed Brown   ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
853552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr);
854c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
855552f7358SJed Brown   ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
856552f7358SJed Brown 
857552f7358SJed Brown   /* Read cell sets information */
858552f7358SJed Brown   if (!rank) {
859552f7358SJed Brown     PetscInt *cone;
860552f7358SJed Brown     int      c, cs, c_loc, v, v_loc;
861552f7358SJed Brown     /* Read from ex_get_elem_blk_ids() */
862552f7358SJed Brown     int *cs_id;
863552f7358SJed Brown     /* Read from ex_get_elem_block() */
864552f7358SJed Brown     char buffer[PETSC_MAX_PATH_LEN+1];
865552f7358SJed Brown     int  num_cell_in_set, num_vertex_per_cell, num_attr;
866552f7358SJed Brown     /* Read from ex_get_elem_conn() */
867552f7358SJed Brown     int *cs_connect;
868552f7358SJed Brown 
869552f7358SJed Brown     /* Get cell sets IDs */
870785e854fSJed Brown     ierr = PetscMalloc1(num_cs, &cs_id);CHKERRQ(ierr);
871f958083aSBlaise Bourdin     ierr = ex_get_ids (exoid, EX_ELEM_BLOCK, cs_id);CHKERRQ(ierr);
872552f7358SJed Brown     /* Read the cell set connectivity table and build mesh topology
873552f7358SJed Brown        EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
874552f7358SJed Brown     /* First set sizes */
875552f7358SJed Brown     for (cs = 0, c = 0; cs < num_cs; cs++) {
876f958083aSBlaise 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);
877552f7358SJed Brown       for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
878552f7358SJed Brown         ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr);
879552f7358SJed Brown       }
880552f7358SJed Brown     }
881552f7358SJed Brown     ierr = DMSetUp(*dm);CHKERRQ(ierr);
882552f7358SJed Brown     for (cs = 0, c = 0; cs < num_cs; cs++) {
883f958083aSBlaise 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);
884dcca6d9dSJed Brown       ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr);
885f958083aSBlaise Bourdin       ierr = ex_get_conn (exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL);CHKERRQ(ierr);
886eaf898f9SPatrick Sanan       /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
887552f7358SJed Brown       for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
888552f7358SJed Brown         for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) {
889552f7358SJed Brown           cone[v_loc] = cs_connect[v]+numCells-1;
890552f7358SJed Brown         }
891731dcdf4SMatthew G. Knepley         if (dim == 3) {
8922e1b13c2SMatthew G. Knepley           /* Tetrahedra are inverted */
8932e1b13c2SMatthew G. Knepley           if (num_vertex_per_cell == 4) {
8942e1b13c2SMatthew G. Knepley             PetscInt tmp = cone[0];
8952e1b13c2SMatthew G. Knepley             cone[0] = cone[1];
8962e1b13c2SMatthew G. Knepley             cone[1] = tmp;
8972e1b13c2SMatthew G. Knepley           }
8982e1b13c2SMatthew G. Knepley           /* Hexahedra are inverted */
8992e1b13c2SMatthew G. Knepley           if (num_vertex_per_cell == 8) {
9002e1b13c2SMatthew G. Knepley             PetscInt tmp = cone[1];
9012e1b13c2SMatthew G. Knepley             cone[1] = cone[3];
9022e1b13c2SMatthew G. Knepley             cone[3] = tmp;
9032e1b13c2SMatthew G. Knepley           }
904731dcdf4SMatthew G. Knepley         }
905552f7358SJed Brown         ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
906c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr);
907552f7358SJed Brown       }
908552f7358SJed Brown       ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr);
909552f7358SJed Brown     }
910552f7358SJed Brown     ierr = PetscFree(cs_id);CHKERRQ(ierr);
911552f7358SJed Brown   }
912552f7358SJed Brown   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
913552f7358SJed Brown   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
914552f7358SJed Brown   if (interpolate) {
9155fd9971aSMatthew G. Knepley     DM idm;
916552f7358SJed Brown 
9179f074e33SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
918552f7358SJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
919552f7358SJed Brown     *dm  = idm;
920552f7358SJed Brown   }
921552f7358SJed Brown 
922552f7358SJed Brown   /* Create vertex set label */
923552f7358SJed Brown   if (!rank && (num_vs > 0)) {
924552f7358SJed Brown     int vs, v;
925552f7358SJed Brown     /* Read from ex_get_node_set_ids() */
926552f7358SJed Brown     int *vs_id;
927552f7358SJed Brown     /* Read from ex_get_node_set_param() */
928f958083aSBlaise Bourdin     int num_vertex_in_set;
929552f7358SJed Brown     /* Read from ex_get_node_set() */
930552f7358SJed Brown     int *vs_vertex_list;
931552f7358SJed Brown 
932552f7358SJed Brown     /* Get vertex set ids */
933785e854fSJed Brown     ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr);
934f958083aSBlaise Bourdin     ierr = ex_get_ids(exoid, EX_NODE_SET, vs_id);CHKERRQ(ierr);
935552f7358SJed Brown     for (vs = 0; vs < num_vs; ++vs) {
936f958083aSBlaise Bourdin       ierr = ex_get_set_param(exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL);
937785e854fSJed Brown       ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr);
938f958083aSBlaise Bourdin       ierr = ex_get_set(exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL);CHKERRQ(ierr);
939552f7358SJed Brown       for (v = 0; v < num_vertex_in_set; ++v) {
940c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr);
941552f7358SJed Brown       }
942552f7358SJed Brown       ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr);
943552f7358SJed Brown     }
944552f7358SJed Brown     ierr = PetscFree(vs_id);CHKERRQ(ierr);
945552f7358SJed Brown   }
946552f7358SJed Brown   /* Read coordinates */
94769d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
948552f7358SJed Brown   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
949552f7358SJed Brown   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
950552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
951552f7358SJed Brown   for (v = numCells; v < numCells+numVertices; ++v) {
952552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
953552f7358SJed Brown     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
954552f7358SJed Brown   }
955552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
956552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
9578b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
958552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
959552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
9604e90ef8eSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
9612eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
962552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
9630aba08f6SKarl Rupp   if (!rank) {
964e84f7031SBlaise Bourdin     PetscReal *x, *y, *z;
965552f7358SJed Brown 
966dcca6d9dSJed Brown     ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr);
967552f7358SJed Brown     ierr = ex_get_coord(exoid, x, y, z);CHKERRQ(ierr);
9680d644c17SKarl Rupp     if (dim > 0) {
9690d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v];
9700d644c17SKarl Rupp     }
9710d644c17SKarl Rupp     if (dim > 1) {
9720d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v];
9730d644c17SKarl Rupp     }
9740d644c17SKarl Rupp     if (dim > 2) {
9750d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v];
9760d644c17SKarl Rupp     }
977552f7358SJed Brown     ierr = PetscFree3(x,y,z);CHKERRQ(ierr);
978552f7358SJed Brown   }
979552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
980552f7358SJed Brown   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
981552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
982552f7358SJed Brown 
983552f7358SJed Brown   /* Create side set label */
984552f7358SJed Brown   if (!rank && interpolate && (num_fs > 0)) {
985552f7358SJed Brown     int fs, f, voff;
986552f7358SJed Brown     /* Read from ex_get_side_set_ids() */
987552f7358SJed Brown     int *fs_id;
988552f7358SJed Brown     /* Read from ex_get_side_set_param() */
989f958083aSBlaise Bourdin     int num_side_in_set;
990552f7358SJed Brown     /* Read from ex_get_side_set_node_list() */
991552f7358SJed Brown     int *fs_vertex_count_list, *fs_vertex_list;
992ef073283Smichael_afanasiev     /* Read side set labels */
993ef073283Smichael_afanasiev     char fs_name[MAX_STR_LENGTH+1];
994552f7358SJed Brown 
995552f7358SJed Brown     /* Get side set ids */
996785e854fSJed Brown     ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr);
9976f1149eeSBlaise Bourdin     ierr = ex_get_ids(exoid, EX_SIDE_SET, fs_id);CHKERRQ(ierr);
998552f7358SJed Brown     for (fs = 0; fs < num_fs; ++fs) {
999f958083aSBlaise Bourdin       ierr = ex_get_set_param(exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL);CHKERRQ(ierr);
1000dcca6d9dSJed Brown       ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr);
1001552f7358SJed Brown       ierr = ex_get_side_set_node_list(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);CHKERRQ(ierr);
1002ef073283Smichael_afanasiev       /* Get the specific name associated with this side set ID. */
1003ef073283Smichael_afanasiev       int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name);
1004552f7358SJed Brown       for (f = 0, voff = 0; f < num_side_in_set; ++f) {
10050298fd71SBarry Smith         const PetscInt *faces   = NULL;
1006552f7358SJed Brown         PetscInt       faceSize = fs_vertex_count_list[f], numFaces;
1007552f7358SJed Brown         PetscInt       faceVertices[4], v;
1008552f7358SJed Brown 
1009552f7358SJed Brown         if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize);
1010552f7358SJed Brown         for (v = 0; v < faceSize; ++v, ++voff) {
1011552f7358SJed Brown           faceVertices[v] = fs_vertex_list[voff]+numCells-1;
1012552f7358SJed Brown         }
1013552f7358SJed Brown         ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
1014552f7358SJed Brown         if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces);
1015c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr);
1016ef073283Smichael_afanasiev         /* Only add the label if one has been detected for this side set. */
1017ef073283Smichael_afanasiev         if (!fs_name_err) {
1018ef073283Smichael_afanasiev           ierr = DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]);CHKERRQ(ierr);
1019ef073283Smichael_afanasiev         }
1020552f7358SJed Brown         ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
1021552f7358SJed Brown       }
1022552f7358SJed Brown       ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr);
1023552f7358SJed Brown     }
1024552f7358SJed Brown     ierr = PetscFree(fs_id);CHKERRQ(ierr);
1025552f7358SJed Brown   }
1026552f7358SJed Brown #else
1027552f7358SJed Brown   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
1028552f7358SJed Brown #endif
1029552f7358SJed Brown   PetscFunctionReturn(0);
1030552f7358SJed Brown }
1031