xref: /petsc/src/dm/impls/plex/plexexodusii.c (revision e8f6893f55b15bd22a434b43a473a8aa4cfe7fdc)
1552f7358SJed Brown #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown 
4519f805aSKarl Rupp #if defined(PETSC_HAVE_EXODUSII)
5552f7358SJed Brown #include <netcdf.h>
6552f7358SJed Brown #include <exodusII.h>
7552f7358SJed Brown #endif
8552f7358SJed Brown 
9e0266925SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
1078569bb4SMatthew G. Knepley /*
1178569bb4SMatthew G. Knepley   EXOGetVarIndex - Locate a result in an exodus file based on its name
1278569bb4SMatthew G. Knepley 
1378569bb4SMatthew G. Knepley   Not collective
1478569bb4SMatthew G. Knepley 
1578569bb4SMatthew G. Knepley   Input Parameters:
1678569bb4SMatthew G. Knepley + exoid    - the exodus id of a file (obtained from ex_open or ex_create for instance)
1778569bb4SMatthew G. Knepley . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK
1878569bb4SMatthew G. Knepley - name     - the name of the result
1978569bb4SMatthew G. Knepley 
2078569bb4SMatthew G. Knepley   Output Parameters:
2178569bb4SMatthew G. Knepley . varIndex - the location in the exodus file of the result
2278569bb4SMatthew G. Knepley 
2378569bb4SMatthew G. Knepley   Notes:
2478569bb4SMatthew G. Knepley   The exodus variable index is obtained by comparing name and the
2578569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance if name is "V"
2678569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
2778569bb4SMatthew G. Knepley   amongst all variables of type obj_type.
2878569bb4SMatthew G. Knepley 
2978569bb4SMatthew G. Knepley   Level: beginner
3078569bb4SMatthew G. Knepley 
31cd921712SStefano Zampini .seealso: EXOGetVarIndex(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadNodal_PlexEXO(),VecLoadZonal_PlexEXO()
3278569bb4SMatthew G. Knepley */
3378569bb4SMatthew G. Knepley static PetscErrorCode EXOGetVarIndex_Private(int exoid, ex_entity_type obj_type, const char name[], int *varIndex)
3478569bb4SMatthew G. Knepley {
356de834b4SMatthew G. Knepley   int            num_vars, i, j;
3678569bb4SMatthew G. Knepley   char           ext_name[MAX_STR_LENGTH+1], var_name[MAX_STR_LENGTH+1];
3778569bb4SMatthew G. Knepley   const int      num_suffix = 5;
3878569bb4SMatthew G. Knepley   char          *suffix[5];
3978569bb4SMatthew G. Knepley   PetscBool      flg;
4078569bb4SMatthew G. Knepley   PetscErrorCode ierr;
4178569bb4SMatthew G. Knepley 
4278569bb4SMatthew G. Knepley   PetscFunctionBegin;
4378569bb4SMatthew G. Knepley   suffix[0] = (char *) "";
4478569bb4SMatthew G. Knepley   suffix[1] = (char *) "_X";
4578569bb4SMatthew G. Knepley   suffix[2] = (char *) "_XX";
4678569bb4SMatthew G. Knepley   suffix[3] = (char *) "_1";
4778569bb4SMatthew G. Knepley   suffix[4] = (char *) "_11";
4878569bb4SMatthew G. Knepley 
4978569bb4SMatthew G. Knepley   *varIndex = 0;
506de834b4SMatthew G. Knepley   PetscStackCallStandard(ex_get_variable_param,(exoid, obj_type, &num_vars));
5178569bb4SMatthew G. Knepley   for (i = 0; i < num_vars; ++i) {
526de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_variable_name,(exoid, obj_type, i+1, var_name));
5378569bb4SMatthew G. Knepley     for (j = 0; j < num_suffix; ++j){
5478569bb4SMatthew G. Knepley       ierr = PetscStrncpy(ext_name, name, MAX_STR_LENGTH);CHKERRQ(ierr);
55a126751eSBarry Smith       ierr = PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH);CHKERRQ(ierr);
5678569bb4SMatthew G. Knepley       ierr = PetscStrcasecmp(ext_name, var_name, &flg);
5778569bb4SMatthew G. Knepley       if (flg) {
5878569bb4SMatthew G. Knepley         *varIndex = i+1;
5978569bb4SMatthew G. Knepley         PetscFunctionReturn(0);
6078569bb4SMatthew G. Knepley       }
6178569bb4SMatthew G. Knepley     }
6278569bb4SMatthew G. Knepley   }
6378569bb4SMatthew G. Knepley   SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to locate variable %s in ExodusII file.", name);
6478569bb4SMatthew G. Knepley  PetscFunctionReturn(-1);
6578569bb4SMatthew G. Knepley }
6678569bb4SMatthew G. Knepley 
6778569bb4SMatthew G. Knepley /*
6878569bb4SMatthew G. Knepley   DMPlexView_ExodusII_Internal - Write a DM to disk in exodus format
6978569bb4SMatthew G. Knepley 
7078569bb4SMatthew G. Knepley   Collective on dm
7178569bb4SMatthew G. Knepley 
7278569bb4SMatthew G. Knepley   Input Parameters:
7378569bb4SMatthew G. Knepley + dm  - The dm to be written
7478569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
7578569bb4SMatthew G. Knepley - degree - the degree of the interpolation space
7678569bb4SMatthew G. Knepley 
7778569bb4SMatthew G. Knepley   Notes:
7878569bb4SMatthew G. Knepley   Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels)
7978569bb4SMatthew G. Knepley   consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted.
8078569bb4SMatthew G. Knepley 
8178569bb4SMatthew 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
8278569bb4SMatthew G. Knepley   (including "ghost" cells and vertices) will be written. If all exoid points to the same file, the resulting file will
8378569bb4SMatthew G. Knepley   probably be corrupted.
8478569bb4SMatthew G. Knepley 
8578569bb4SMatthew G. Knepley   DMPlex only represents geometry while most post-processing software expect that a mesh also provides information
8678569bb4SMatthew G. Knepley   on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2.
8778569bb4SMatthew G. Knepley   It should be extended to use PetscFE objects.
8878569bb4SMatthew G. Knepley 
8978569bb4SMatthew G. Knepley   This function will only handle TRI, TET, QUAD and HEX cells.
9078569bb4SMatthew G. Knepley   Level: beginner
9178569bb4SMatthew G. Knepley 
9278569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadNodal_PlexEXO(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO()
9378569bb4SMatthew G. Knepley */
9478569bb4SMatthew G. Knepley PetscErrorCode DMPlexView_ExodusII_Internal(DM dm, int exoid, PetscInt degree)
9578569bb4SMatthew G. Knepley {
9678569bb4SMatthew G. Knepley   enum ElemType {TRI, QUAD, TET, HEX};
9778569bb4SMatthew G. Knepley   MPI_Comm        comm;
9878569bb4SMatthew G. Knepley   /* Connectivity Variables */
9978569bb4SMatthew G. Knepley   PetscInt        cellsNotInConnectivity;
10078569bb4SMatthew G. Knepley   /* Cell Sets */
10178569bb4SMatthew G. Knepley   DMLabel         csLabel;
10278569bb4SMatthew G. Knepley   IS              csIS;
10378569bb4SMatthew G. Knepley   const PetscInt *csIdx;
10478569bb4SMatthew G. Knepley   PetscInt        num_cs, cs;
10578569bb4SMatthew G. Knepley   enum ElemType  *type;
106fe209ef9SBlaise Bourdin   PetscBool       hasLabel;
10778569bb4SMatthew G. Knepley   /* Coordinate Variables */
10878569bb4SMatthew G. Knepley   DM              cdm;
10978569bb4SMatthew G. Knepley   PetscSection    section;
11078569bb4SMatthew G. Knepley   Vec             coord;
11178569bb4SMatthew G. Knepley   PetscInt      **nodes;
112eae8a387SMatthew G. Knepley   PetscInt        depth, d, dim, skipCells = 0;
11378569bb4SMatthew G. Knepley   PetscInt        pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes;
11478569bb4SMatthew G. Knepley   PetscInt        num_vs, num_fs;
11578569bb4SMatthew G. Knepley   PetscMPIInt     rank, size;
11678569bb4SMatthew G. Knepley   const char     *dmName;
117fe209ef9SBlaise Bourdin   PetscInt        nodesTriP1[4]  = {3,0,0,0};
118fe209ef9SBlaise Bourdin   PetscInt        nodesTriP2[4]  = {3,3,0,0};
119fe209ef9SBlaise Bourdin   PetscInt        nodesQuadP1[4] = {4,0,0,0};
120fe209ef9SBlaise Bourdin   PetscInt        nodesQuadP2[4] = {4,4,0,1};
121fe209ef9SBlaise Bourdin   PetscInt        nodesTetP1[4]  = {4,0,0,0};
122fe209ef9SBlaise Bourdin   PetscInt        nodesTetP2[4]  = {4,6,0,0};
123fe209ef9SBlaise Bourdin   PetscInt        nodesHexP1[4]  = {8,0,0,0};
124fe209ef9SBlaise Bourdin   PetscInt        nodesHexP2[4]  = {8,12,6,1};
12578569bb4SMatthew G. Knepley   PetscErrorCode  ierr;
12678569bb4SMatthew G. Knepley 
12778569bb4SMatthew G. Knepley   PetscFunctionBegin;
12878569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
12978569bb4SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
13078569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
13178569bb4SMatthew G. Knepley   /* --- Get DM info --- */
13278569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) dm, &dmName);CHKERRQ(ierr);
13378569bb4SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
13478569bb4SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
13578569bb4SMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
13678569bb4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
13778569bb4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
13878569bb4SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);
13978569bb4SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
14078569bb4SMatthew G. Knepley   numCells    = cEnd - cStart;
14178569bb4SMatthew G. Knepley   numEdges    = eEnd - eStart;
14278569bb4SMatthew G. Knepley   numVertices = vEnd - vStart;
14378569bb4SMatthew G. Knepley   if (depth == 3) {numFaces = fEnd - fStart;}
14478569bb4SMatthew G. Knepley   else            {numFaces = 0;}
14578569bb4SMatthew G. Knepley   ierr = DMGetLabelSize(dm, "Cell Sets", &num_cs);CHKERRQ(ierr);
14678569bb4SMatthew G. Knepley   ierr = DMGetLabelSize(dm, "Vertex Sets", &num_vs);CHKERRQ(ierr);
14778569bb4SMatthew G. Knepley   ierr = DMGetLabelSize(dm, "Face Sets", &num_fs);CHKERRQ(ierr);
14878569bb4SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr);
14978569bb4SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
15078569bb4SMatthew G. Knepley   if (num_cs > 0) {
15178569bb4SMatthew G. Knepley     ierr = DMGetLabel(dm, "Cell Sets", &csLabel);CHKERRQ(ierr);
15278569bb4SMatthew G. Knepley     ierr = DMLabelGetValueIS(csLabel, &csIS);CHKERRQ(ierr);
15378569bb4SMatthew G. Knepley     ierr = ISGetIndices(csIS, &csIdx);CHKERRQ(ierr);
15478569bb4SMatthew G. Knepley   }
15578569bb4SMatthew G. Knepley   ierr = PetscMalloc1(num_cs, &nodes);CHKERRQ(ierr);
15678569bb4SMatthew G. Knepley   /* Set element type for each block and compute total number of nodes */
15778569bb4SMatthew G. Knepley   ierr = PetscMalloc1(num_cs, &type);CHKERRQ(ierr);
15878569bb4SMatthew G. Knepley   numNodes = numVertices;
15978569bb4SMatthew G. Knepley   if (degree == 2) {numNodes += numEdges;}
16078569bb4SMatthew G. Knepley   cellsNotInConnectivity = numCells;
16178569bb4SMatthew G. Knepley   for (cs = 0; cs < num_cs; ++cs) {
16278569bb4SMatthew G. Knepley     IS              stratumIS;
16378569bb4SMatthew G. Knepley     const PetscInt *cells;
16478569bb4SMatthew G. Knepley     PetscScalar    *xyz = NULL;
16578569bb4SMatthew G. Knepley     PetscInt        csSize, closureSize;
16678569bb4SMatthew G. Knepley 
16778569bb4SMatthew G. Knepley     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
16878569bb4SMatthew G. Knepley     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
16978569bb4SMatthew G. Knepley     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
17078569bb4SMatthew G. Knepley     ierr = DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr);
17178569bb4SMatthew G. Knepley     switch (dim) {
17278569bb4SMatthew G. Knepley       case 2:
17378569bb4SMatthew G. Knepley         if      (closureSize == 3*dim) {type[cs] = TRI;}
17478569bb4SMatthew G. Knepley         else if (closureSize == 4*dim) {type[cs] = QUAD;}
17578569bb4SMatthew 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);
17678569bb4SMatthew G. Knepley         break;
17778569bb4SMatthew G. Knepley       case 3:
17878569bb4SMatthew G. Knepley         if      (closureSize == 4*dim) {type[cs] = TET;}
17978569bb4SMatthew G. Knepley         else if (closureSize == 8*dim) {type[cs] = HEX;}
18078569bb4SMatthew 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);
18178569bb4SMatthew G. Knepley         break;
18278569bb4SMatthew G. Knepley       default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not handled by ExodusII viewer", dim);
18378569bb4SMatthew G. Knepley     }
18478569bb4SMatthew G. Knepley     if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;}
18578569bb4SMatthew G. Knepley     if ((degree == 2) && (type[cs] == HEX))  {numNodes += csSize; numNodes += numFaces;}
18678569bb4SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr);
18778569bb4SMatthew G. Knepley     /* Set nodes and Element type */
18878569bb4SMatthew G. Knepley     if (type[cs] == TRI) {
18978569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesTriP1;
19078569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesTriP2;
19178569bb4SMatthew G. Knepley     } else if (type[cs] == QUAD) {
19278569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesQuadP1;
19378569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesQuadP2;
19478569bb4SMatthew G. Knepley     } else if (type[cs] == TET) {
19578569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesTetP1;
19678569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesTetP2;
19778569bb4SMatthew G. Knepley     } else if (type[cs] == HEX) {
19878569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesHexP1;
19978569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesHexP2;
20078569bb4SMatthew G. Knepley     }
20178569bb4SMatthew G. Knepley     /* Compute the number of cells not in the connectivity table */
20278569bb4SMatthew G. Knepley     cellsNotInConnectivity -= nodes[cs][3]*csSize;
20378569bb4SMatthew G. Knepley 
20478569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
20578569bb4SMatthew G. Knepley     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
20678569bb4SMatthew G. Knepley   }
2076de834b4SMatthew G. Knepley   if (num_cs > 0) {PetscStackCallStandard(ex_put_init,(exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs));}
20878569bb4SMatthew G. Knepley   /* --- Connectivity --- */
20978569bb4SMatthew G. Knepley   for (cs = 0; cs < num_cs; ++cs) {
21078569bb4SMatthew G. Knepley     IS              stratumIS;
21178569bb4SMatthew G. Knepley     const PetscInt *cells;
212fe209ef9SBlaise Bourdin     PetscInt       *connect, off = 0;
213eae8a387SMatthew G. Knepley     PetscInt        edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0;
21478569bb4SMatthew G. Knepley     PetscInt        csSize, c, connectSize, closureSize;
21578569bb4SMatthew G. Knepley     char           *elem_type = NULL;
21678569bb4SMatthew G. Knepley     char            elem_type_tri3[]  = "TRI3",  elem_type_quad4[] = "QUAD4";
21778569bb4SMatthew G. Knepley     char            elem_type_tri6[]  = "TRI6",  elem_type_quad9[] = "QUAD9";
21878569bb4SMatthew G. Knepley     char            elem_type_tet4[]  = "TET4",  elem_type_hex8[]  = "HEX8";
21978569bb4SMatthew G. Knepley     char            elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27";
22078569bb4SMatthew G. Knepley 
22178569bb4SMatthew G. Knepley     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
22278569bb4SMatthew G. Knepley     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
22378569bb4SMatthew G. Knepley     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
22478569bb4SMatthew G. Knepley     /* Set Element type */
22578569bb4SMatthew G. Knepley     if (type[cs] == TRI) {
22678569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_tri3;
22778569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_tri6;
22878569bb4SMatthew G. Knepley     } else if (type[cs] == QUAD) {
22978569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_quad4;
23078569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_quad9;
23178569bb4SMatthew G. Knepley     } else if (type[cs] == TET) {
23278569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_tet4;
23378569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_tet10;
23478569bb4SMatthew G. Knepley     } else if (type[cs] == HEX) {
23578569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_hex8;
23678569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_hex27;
23778569bb4SMatthew G. Knepley     }
23878569bb4SMatthew G. Knepley     connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3];
239fe209ef9SBlaise Bourdin     ierr = PetscMalloc1(PetscMax(27,connectSize)*csSize, &connect);CHKERRQ(ierr);
240fe209ef9SBlaise Bourdin     PetscStackCallStandard(ex_put_block,(exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1));
24178569bb4SMatthew G. Knepley     /* Find number of vertices, edges, and faces in the closure */
24278569bb4SMatthew G. Knepley     verticesInClosure = nodes[cs][0];
24378569bb4SMatthew G. Knepley     if (depth > 1) {
24478569bb4SMatthew G. Knepley       if (dim == 2) {
24578569bb4SMatthew G. Knepley         ierr = DMPlexGetConeSize(dm, cells[0], &edgesInClosure);CHKERRQ(ierr);
24678569bb4SMatthew G. Knepley       } else if (dim == 3) {
24778569bb4SMatthew G. Knepley         PetscInt *closure = NULL;
24878569bb4SMatthew G. Knepley 
24978569bb4SMatthew G. Knepley         ierr = DMPlexGetConeSize(dm, cells[0], &facesInClosure);CHKERRQ(ierr);
25078569bb4SMatthew G. Knepley         ierr = DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
25178569bb4SMatthew G. Knepley         edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure;
25278569bb4SMatthew G. Knepley         ierr = DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
25378569bb4SMatthew G. Knepley       }
25478569bb4SMatthew G. Knepley     }
25578569bb4SMatthew G. Knepley     /* Get connectivity for each cell */
25678569bb4SMatthew G. Knepley     for (c = 0; c < csSize; ++c) {
25778569bb4SMatthew G. Knepley       PetscInt *closure = NULL;
25878569bb4SMatthew G. Knepley       PetscInt  temp, i;
25978569bb4SMatthew G. Knepley 
26078569bb4SMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
26178569bb4SMatthew G. Knepley       for (i = 0; i < connectSize; ++i) {
26278569bb4SMatthew G. Knepley         if (i < nodes[cs][0]) {/* Vertices */
263fe209ef9SBlaise Bourdin           connect[i+off] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1;
264fe209ef9SBlaise Bourdin           connect[i+off] -= cellsNotInConnectivity;
26578569bb4SMatthew G. Knepley         } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */
266fe209ef9SBlaise Bourdin           connect[i+off] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1;
267fe209ef9SBlaise Bourdin           if (nodes[cs][2] == 0) connect[i+off] -= numFaces;
268fe209ef9SBlaise Bourdin           connect[i+off] -= cellsNotInConnectivity;
26978569bb4SMatthew G. Knepley         } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]){ /* Cells */
270fe209ef9SBlaise Bourdin           connect[i+off] = closure[0] + 1;
271fe209ef9SBlaise Bourdin           connect[i+off] -= skipCells;
27278569bb4SMatthew G. Knepley         } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]){ /* Faces */
273fe209ef9SBlaise Bourdin           connect[i+off] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1;
274fe209ef9SBlaise Bourdin           connect[i+off] -= cellsNotInConnectivity;
27578569bb4SMatthew G. Knepley         } else {
276fe209ef9SBlaise Bourdin           connect[i+off] = -1;
27778569bb4SMatthew G. Knepley         }
27878569bb4SMatthew G. Knepley       }
27978569bb4SMatthew G. Knepley       /* Tetrahedra are inverted */
28078569bb4SMatthew G. Knepley       if (type[cs] == TET) {
281fe209ef9SBlaise Bourdin         temp = connect[0+off]; connect[0+off] = connect[1+off]; connect[1+off] = temp;
28278569bb4SMatthew G. Knepley         if (degree == 2) {
283e2c9586dSBlaise Bourdin           temp = connect[5+off]; connect[5+off] = connect[6+off]; connect[6+off] = temp;
284e2c9586dSBlaise Bourdin           temp = connect[7+off]; connect[7+off] = connect[8+off]; connect[8+off] = temp;
28578569bb4SMatthew G. Knepley         }
28678569bb4SMatthew G. Knepley       }
28778569bb4SMatthew G. Knepley       /* Hexahedra are inverted */
28878569bb4SMatthew G. Knepley       if (type[cs] == HEX) {
289fe209ef9SBlaise Bourdin         temp = connect[1+off]; connect[1+off] = connect[3+off]; connect[3+off] = temp;
29078569bb4SMatthew G. Knepley         if (degree == 2) {
291fe209ef9SBlaise Bourdin           temp = connect[8+off];  connect[8+off]  = connect[11+off]; connect[11+off] = temp;
292fe209ef9SBlaise Bourdin           temp = connect[9+off];  connect[9+off]  = connect[10+off]; connect[10+off] = temp;
293fe209ef9SBlaise Bourdin           temp = connect[16+off]; connect[16+off] = connect[17+off]; connect[17+off] = temp;
294fe209ef9SBlaise Bourdin           temp = connect[18+off]; connect[18+off] = connect[19+off]; connect[19+off] = temp;
29578569bb4SMatthew G. Knepley 
296fe209ef9SBlaise Bourdin           temp = connect[12+off]; connect[12+off] = connect[16+off]; connect[16+off] = temp;
297fe209ef9SBlaise Bourdin           temp = connect[13+off]; connect[13+off] = connect[17+off]; connect[17+off] = temp;
298fe209ef9SBlaise Bourdin           temp = connect[14+off]; connect[14+off] = connect[18+off]; connect[18+off] = temp;
299fe209ef9SBlaise Bourdin           temp = connect[15+off]; connect[15+off] = connect[19+off]; connect[19+off] = temp;
30078569bb4SMatthew G. Knepley 
301fe209ef9SBlaise Bourdin           temp = connect[23+off]; connect[23+off] = connect[26+off]; connect[26+off] = temp;
302fe209ef9SBlaise Bourdin           temp = connect[24+off]; connect[24+off] = connect[25+off]; connect[25+off] = temp;
303fe209ef9SBlaise Bourdin           temp = connect[25+off]; connect[25+off] = connect[26+off]; connect[26+off] = temp;
30478569bb4SMatthew G. Knepley         }
30578569bb4SMatthew G. Knepley       }
306fe209ef9SBlaise Bourdin       off += connectSize;
30778569bb4SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
30878569bb4SMatthew G. Knepley     }
309fe209ef9SBlaise Bourdin     PetscStackCallStandard(ex_put_conn,(exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0));
31078569bb4SMatthew G. Knepley     skipCells += (nodes[cs][3] == 0)*csSize;
31178569bb4SMatthew G. Knepley     ierr = PetscFree(connect);CHKERRQ(ierr);
31278569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
31378569bb4SMatthew G. Knepley     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
31478569bb4SMatthew G. Knepley   }
31578569bb4SMatthew G. Knepley   ierr = PetscFree(type);CHKERRQ(ierr);
31678569bb4SMatthew G. Knepley   /* --- Coordinates --- */
31778569bb4SMatthew G. Knepley   ierr = PetscSectionCreate(comm, &section);CHKERRQ(ierr);
31878569bb4SMatthew G. Knepley   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
31978569bb4SMatthew G. Knepley   for (d = 0; d < depth; ++d) {
32078569bb4SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
32178569bb4SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
32278569bb4SMatthew G. Knepley       ierr = PetscSectionSetDof(section, p, nodes[0][d] > 0);CHKERRQ(ierr);
32378569bb4SMatthew G. Knepley     }
32478569bb4SMatthew G. Knepley   }
32578569bb4SMatthew G. Knepley   for (cs = 0; cs < num_cs; ++cs) {
32678569bb4SMatthew G. Knepley     IS              stratumIS;
32778569bb4SMatthew G. Knepley     const PetscInt *cells;
32878569bb4SMatthew G. Knepley     PetscInt        csSize, c;
32978569bb4SMatthew G. Knepley 
33078569bb4SMatthew G. Knepley     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
33178569bb4SMatthew G. Knepley     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
33278569bb4SMatthew G. Knepley     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
33378569bb4SMatthew G. Knepley     for (c = 0; c < csSize; ++c) {
33478569bb4SMatthew G. Knepley       ierr = PetscSectionSetDof(section, cells[c], nodes[cs][3] > 0);CHKERRQ(ierr);
33578569bb4SMatthew G. Knepley     }
33678569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
33778569bb4SMatthew G. Knepley     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
33878569bb4SMatthew G. Knepley   }
33978569bb4SMatthew G. Knepley   if (num_cs > 0) {
34078569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(csIS, &csIdx);CHKERRQ(ierr);
34178569bb4SMatthew G. Knepley     ierr = ISDestroy(&csIS);CHKERRQ(ierr);
34278569bb4SMatthew G. Knepley   }
34378569bb4SMatthew G. Knepley   ierr = PetscFree(nodes);CHKERRQ(ierr);
34478569bb4SMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
34578569bb4SMatthew G. Knepley   if (numNodes > 0) {
34678569bb4SMatthew G. Knepley     const char  *coordNames[3] = {"x", "y", "z"};
34778569bb4SMatthew G. Knepley     PetscScalar *coords, *closure;
34878569bb4SMatthew G. Knepley     PetscReal   *cval;
34978569bb4SMatthew G. Knepley     PetscInt     hasDof, n = 0;
35078569bb4SMatthew G. Knepley 
35178569bb4SMatthew G. Knepley     /* There can't be more than 24 values in the closure of a point for the coord section */
3526de834b4SMatthew G. Knepley     ierr = PetscCalloc3(numNodes*3, &coords, dim, &cval, 24, &closure);CHKERRQ(ierr);
35378569bb4SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr);
35478569bb4SMatthew G. Knepley     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
35578569bb4SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
35678569bb4SMatthew G. Knepley       ierr = PetscSectionGetDof(section, p, &hasDof);
35778569bb4SMatthew G. Knepley       if (hasDof) {
35878569bb4SMatthew G. Knepley         PetscInt closureSize = 24, j;
35978569bb4SMatthew G. Knepley 
36078569bb4SMatthew G. Knepley         ierr = DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure);CHKERRQ(ierr);
36178569bb4SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
36278569bb4SMatthew G. Knepley           cval[d] = 0.0;
36378569bb4SMatthew G. Knepley           for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d];
36478569bb4SMatthew G. Knepley           coords[d*numNodes+n] = cval[d] * dim / closureSize;
36578569bb4SMatthew G. Knepley         }
36678569bb4SMatthew G. Knepley         ++n;
36778569bb4SMatthew G. Knepley       }
36878569bb4SMatthew G. Knepley     }
3696de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_put_coord,(exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes]));
37078569bb4SMatthew G. Knepley     ierr = PetscFree3(coords, cval, closure);CHKERRQ(ierr);
3716de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_put_coord_names,(exoid, (char **) coordNames));
37278569bb4SMatthew G. Knepley   }
37378569bb4SMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
374fe209ef9SBlaise Bourdin   /* --- Node Sets/Vertex Sets --- */
375fe209ef9SBlaise Bourdin   DMHasLabel(dm, "Vertex Sets", &hasLabel);
376fe209ef9SBlaise Bourdin   if (hasLabel) {
377fe209ef9SBlaise Bourdin     PetscInt        i, vs, vsSize;
378fe209ef9SBlaise Bourdin     const PetscInt *vsIdx, *vertices;
379fe209ef9SBlaise Bourdin     PetscInt       *nodeList;
380fe209ef9SBlaise Bourdin     IS              vsIS, stratumIS;
381fe209ef9SBlaise Bourdin     DMLabel         vsLabel;
382fe209ef9SBlaise Bourdin     ierr = DMGetLabel(dm, "Vertex Sets", &vsLabel);CHKERRQ(ierr);
383fe209ef9SBlaise Bourdin     ierr = DMLabelGetValueIS(vsLabel, &vsIS);CHKERRQ(ierr);
384fe209ef9SBlaise Bourdin     ierr = ISGetIndices(vsIS, &vsIdx);CHKERRQ(ierr);
385fe209ef9SBlaise Bourdin     for (vs=0; vs<num_vs; ++vs) {
386fe209ef9SBlaise Bourdin       ierr = DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS);CHKERRQ(ierr);
387fe209ef9SBlaise Bourdin       ierr = ISGetIndices(stratumIS, &vertices);CHKERRQ(ierr);
388fe209ef9SBlaise Bourdin       ierr = ISGetSize(stratumIS, &vsSize);CHKERRQ(ierr);
389fe209ef9SBlaise Bourdin       ierr = PetscMalloc1(vsSize, &nodeList);
390fe209ef9SBlaise Bourdin       for (i=0; i<vsSize; ++i) {
391fe209ef9SBlaise Bourdin         nodeList[i] = vertices[i] - skipCells + 1;
392fe209ef9SBlaise Bourdin       }
393fe209ef9SBlaise Bourdin       PetscStackCallStandard(ex_put_set_param,(exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0));
394fe209ef9SBlaise Bourdin       PetscStackCallStandard(ex_put_set,(exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL));
395fe209ef9SBlaise Bourdin       ierr = ISRestoreIndices(stratumIS, &vertices);CHKERRQ(ierr);
396fe209ef9SBlaise Bourdin       ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
397fe209ef9SBlaise Bourdin       ierr = PetscFree(nodeList);CHKERRQ(ierr);
398fe209ef9SBlaise Bourdin     }
399fe209ef9SBlaise Bourdin     ierr = ISRestoreIndices(vsIS, &vsIdx);CHKERRQ(ierr);
400fe209ef9SBlaise Bourdin     ierr = ISDestroy(&vsIS);CHKERRQ(ierr);
401fe209ef9SBlaise Bourdin   }
402fe209ef9SBlaise Bourdin   /* --- Side Sets/Face Sets --- */
403fe209ef9SBlaise Bourdin   ierr = DMHasLabel(dm, "Face Sets", &hasLabel);CHKERRQ(ierr);
404fe209ef9SBlaise Bourdin   if (hasLabel) {
405fe209ef9SBlaise Bourdin     PetscInt        i, j, fs, fsSize;
406fe209ef9SBlaise Bourdin     const PetscInt *fsIdx, *faces;
407fe209ef9SBlaise Bourdin     IS              fsIS, stratumIS;
408fe209ef9SBlaise Bourdin     DMLabel         fsLabel;
409fe209ef9SBlaise Bourdin     PetscInt        numPoints, *points;
410fe209ef9SBlaise Bourdin     PetscInt        elem_list_size = 0;
411fe209ef9SBlaise Bourdin     PetscInt       *elem_list, *elem_ind, *side_list;
412fe209ef9SBlaise Bourdin 
413fe209ef9SBlaise Bourdin     ierr = DMGetLabel(dm, "Face Sets", &fsLabel);CHKERRQ(ierr);
414fe209ef9SBlaise Bourdin     /* Compute size of Node List and Element List */
415fe209ef9SBlaise Bourdin     ierr = DMLabelGetValueIS(fsLabel, &fsIS);CHKERRQ(ierr);
416fe209ef9SBlaise Bourdin     ierr = ISGetIndices(fsIS, &fsIdx);CHKERRQ(ierr);
417fe209ef9SBlaise Bourdin     for (fs=0; fs<num_fs; ++fs) {
418fe209ef9SBlaise Bourdin       ierr = DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS);CHKERRQ(ierr);;
419fe209ef9SBlaise Bourdin       ierr = ISGetSize(stratumIS, &fsSize);CHKERRQ(ierr);
420fe209ef9SBlaise Bourdin       elem_list_size += fsSize;
421fe209ef9SBlaise Bourdin       ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
422fe209ef9SBlaise Bourdin     }
423fe209ef9SBlaise Bourdin     ierr = PetscMalloc3(num_fs, &elem_ind,
424fe209ef9SBlaise Bourdin                         elem_list_size, &elem_list,
425fe209ef9SBlaise Bourdin                         elem_list_size, &side_list);CHKERRQ(ierr);
426fe209ef9SBlaise Bourdin     elem_ind[0] = 0;
427fe209ef9SBlaise Bourdin     for (fs=0; fs<num_fs; ++fs) {
428fe209ef9SBlaise Bourdin       ierr = DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS);CHKERRQ(ierr);
429fe209ef9SBlaise Bourdin       ierr = ISGetIndices(stratumIS, &faces);CHKERRQ(ierr);
430fe209ef9SBlaise Bourdin       ierr = ISGetSize(stratumIS, &fsSize);CHKERRQ(ierr);
431fe209ef9SBlaise Bourdin       /* Set Parameters */
432fe209ef9SBlaise Bourdin       PetscStackCallStandard(ex_put_set_param,(exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0));
433fe209ef9SBlaise Bourdin       /* Indices */
434fe209ef9SBlaise Bourdin       if (fs<num_fs-1) {
435fe209ef9SBlaise Bourdin         elem_ind[fs+1] = elem_ind[fs] + fsSize;
436fe209ef9SBlaise Bourdin       }
437fe209ef9SBlaise Bourdin 
438fe209ef9SBlaise Bourdin       for (i=0; i<fsSize; ++i) {
439fe209ef9SBlaise Bourdin         /* Element List */
440fe209ef9SBlaise Bourdin         points = NULL;
441fe209ef9SBlaise Bourdin         ierr = DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points);CHKERRQ(ierr);
442fe209ef9SBlaise Bourdin         elem_list[elem_ind[fs] + i] = points[2] +1;
443fe209ef9SBlaise Bourdin         ierr = DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points);CHKERRQ(ierr);
444fe209ef9SBlaise Bourdin 
445fe209ef9SBlaise Bourdin         /* Side List */
446fe209ef9SBlaise Bourdin         points = NULL;
447fe209ef9SBlaise Bourdin         ierr = DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
448fe209ef9SBlaise Bourdin         for (j=1; j<numPoints; ++j) {
449fe209ef9SBlaise Bourdin           if (points[j*2]==faces[i]) {break;}
450fe209ef9SBlaise Bourdin         }
451fe209ef9SBlaise Bourdin         /* Convert HEX sides */
452fe209ef9SBlaise Bourdin         if (numPoints == 27) {
453fe209ef9SBlaise Bourdin           if      (j == 1) {j = 5;}
454fe209ef9SBlaise Bourdin           else if (j == 2) {j = 6;}
455fe209ef9SBlaise Bourdin           else if (j == 3) {j = 1;}
456fe209ef9SBlaise Bourdin           else if (j == 4) {j = 3;}
457fe209ef9SBlaise Bourdin           else if (j == 5) {j = 2;}
458fe209ef9SBlaise Bourdin           else if (j == 6) {j = 4;}
459fe209ef9SBlaise Bourdin         }
460fe209ef9SBlaise Bourdin         /* Convert TET sides */
461fe209ef9SBlaise Bourdin         if (numPoints == 15) {
462fe209ef9SBlaise Bourdin           --j;
463fe209ef9SBlaise Bourdin           if (j == 0) {j = 4;}
464fe209ef9SBlaise Bourdin         }
465fe209ef9SBlaise Bourdin         side_list[elem_ind[fs] + i] = j;
466fe209ef9SBlaise Bourdin         ierr = DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
467fe209ef9SBlaise Bourdin 
468fe209ef9SBlaise Bourdin       }
469fe209ef9SBlaise Bourdin       ierr = ISRestoreIndices(stratumIS, &faces);CHKERRQ(ierr);
470fe209ef9SBlaise Bourdin       ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
471fe209ef9SBlaise Bourdin     }
472fe209ef9SBlaise Bourdin     ierr = ISRestoreIndices(fsIS, &fsIdx);CHKERRQ(ierr);
473fe209ef9SBlaise Bourdin     ierr = ISDestroy(&fsIS);CHKERRQ(ierr);
474fe209ef9SBlaise Bourdin 
475fe209ef9SBlaise Bourdin     /* Put side sets */
476fe209ef9SBlaise Bourdin     for (fs=0; fs<num_fs; ++fs) {
477fe209ef9SBlaise Bourdin       PetscStackCallStandard(ex_put_set,(exoid, EX_SIDE_SET, fsIdx[fs], &elem_list[elem_ind[fs]], &side_list[elem_ind[fs]]));
478fe209ef9SBlaise Bourdin     }
479fe209ef9SBlaise Bourdin     ierr = PetscFree3(elem_ind,elem_list,side_list);CHKERRQ(ierr);
480fe209ef9SBlaise Bourdin   }
48178569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
48278569bb4SMatthew G. Knepley }
48378569bb4SMatthew G. Knepley 
48478569bb4SMatthew G. Knepley /*
48578569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file
48678569bb4SMatthew G. Knepley 
48778569bb4SMatthew G. Knepley   Collective on v
48878569bb4SMatthew G. Knepley 
48978569bb4SMatthew G. Knepley   Input Parameters:
49078569bb4SMatthew G. Knepley + v  - The vector to be written
49178569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
49278569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1)
49378569bb4SMatthew G. Knepley 
49478569bb4SMatthew G. Knepley   Notes:
49578569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
49678569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
49778569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
49878569bb4SMatthew G. Knepley   amongst all nodal variables.
49978569bb4SMatthew G. Knepley 
50078569bb4SMatthew G. Knepley   Level: beginner
50178569bb4SMatthew G. Knepley 
50278569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecLoadNodal_PlexEXO(),VecViewZonal_PlexEXO(),VecLoadZonal_PlexEXO()
50378569bb4SMatthew G. Knepley @*/
50478569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step)
50578569bb4SMatthew G. Knepley {
50678569bb4SMatthew G. Knepley   MPI_Comm         comm;
50778569bb4SMatthew G. Knepley   PetscMPIInt      size;
50878569bb4SMatthew G. Knepley   DM               dm;
50978569bb4SMatthew G. Knepley   Vec              vNatural, vComp;
51022a7544eSVaclav Hapla   const PetscScalar *varray;
51178569bb4SMatthew G. Knepley   const char      *vecname;
51278569bb4SMatthew G. Knepley   PetscInt         xs, xe, bs;
51378569bb4SMatthew G. Knepley   PetscBool        useNatural;
51478569bb4SMatthew G. Knepley   int              offset;
51578569bb4SMatthew G. Knepley   PetscErrorCode   ierr;
51678569bb4SMatthew G. Knepley 
51778569bb4SMatthew G. Knepley   PetscFunctionBegin;
51878569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
51978569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
52078569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
52178569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
52278569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
52378569bb4SMatthew G. Knepley   if (useNatural) {
52478569bb4SMatthew G. Knepley     ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr);
52578569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr);
52678569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr);
52778569bb4SMatthew G. Knepley   } else {
52878569bb4SMatthew G. Knepley     vNatural = v;
52978569bb4SMatthew G. Knepley   }
53078569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
53178569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
53278569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr);
53378569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file\n", vecname);
53478569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
53578569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
53678569bb4SMatthew G. Knepley      We assume that they are stored sequentially */
53778569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
53878569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
53978569bb4SMatthew G. Knepley   if (bs == 1) {
54078569bb4SMatthew G. Knepley     ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr);
5416de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray));
54278569bb4SMatthew G. Knepley     ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr);
54378569bb4SMatthew G. Knepley   } else {
54478569bb4SMatthew G. Knepley     IS       compIS;
54578569bb4SMatthew G. Knepley     PetscInt c;
54678569bb4SMatthew G. Knepley 
54778569bb4SMatthew G. Knepley     ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);
54878569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
54978569bb4SMatthew G. Knepley       ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
55078569bb4SMatthew G. Knepley       ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
55178569bb4SMatthew G. Knepley       ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr);
5526de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray));
55378569bb4SMatthew G. Knepley       ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr);
55478569bb4SMatthew G. Knepley       ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
55578569bb4SMatthew G. Knepley     }
55678569bb4SMatthew G. Knepley     ierr = ISDestroy(&compIS);CHKERRQ(ierr);
55778569bb4SMatthew G. Knepley   }
55878569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);}
55978569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
56078569bb4SMatthew G. Knepley }
56178569bb4SMatthew G. Knepley 
56278569bb4SMatthew G. Knepley /*
56378569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file
56478569bb4SMatthew G. Knepley 
56578569bb4SMatthew G. Knepley   Collective on v
56678569bb4SMatthew G. Knepley 
56778569bb4SMatthew G. Knepley   Input Parameters:
56878569bb4SMatthew G. Knepley + v  - The vector to be written
56978569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
57078569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1)
57178569bb4SMatthew G. Knepley 
57278569bb4SMatthew G. Knepley   Notes:
57378569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
57478569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
57578569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
57678569bb4SMatthew G. Knepley   amongst all nodal variables.
57778569bb4SMatthew G. Knepley 
57878569bb4SMatthew G. Knepley   Level: beginner
57978569bb4SMatthew G. Knepley 
58078569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO()
58178569bb4SMatthew G. Knepley */
58278569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step)
58378569bb4SMatthew G. Knepley {
58478569bb4SMatthew G. Knepley   MPI_Comm       comm;
58578569bb4SMatthew G. Knepley   PetscMPIInt    size;
58678569bb4SMatthew G. Knepley   DM             dm;
58778569bb4SMatthew G. Knepley   Vec            vNatural, vComp;
58822a7544eSVaclav Hapla   PetscScalar   *varray;
58978569bb4SMatthew G. Knepley   const char    *vecname;
59078569bb4SMatthew G. Knepley   PetscInt       xs, xe, bs;
59178569bb4SMatthew G. Knepley   PetscBool      useNatural;
59278569bb4SMatthew G. Knepley   int            offset;
59378569bb4SMatthew G. Knepley   PetscErrorCode ierr;
59478569bb4SMatthew G. Knepley 
59578569bb4SMatthew G. Knepley   PetscFunctionBegin;
59678569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
59778569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
59878569bb4SMatthew G. Knepley   ierr = VecGetDM(v,&dm);CHKERRQ(ierr);
59978569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
60078569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
60178569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
60278569bb4SMatthew G. Knepley   else            {vNatural = v;}
60378569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
60478569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
60578569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr);
60678569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file", vecname);
60778569bb4SMatthew G. Knepley   /* Read local chunk from the file */
60878569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
60978569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
61078569bb4SMatthew G. Knepley   if (bs == 1) {
61178569bb4SMatthew G. Knepley     ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr);
6126de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray));
61378569bb4SMatthew G. Knepley     ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr);
61478569bb4SMatthew G. Knepley   } else {
61578569bb4SMatthew G. Knepley     IS       compIS;
61678569bb4SMatthew G. Knepley     PetscInt c;
61778569bb4SMatthew G. Knepley 
61878569bb4SMatthew G. Knepley     ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);
61978569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
62078569bb4SMatthew G. Knepley       ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
62178569bb4SMatthew G. Knepley       ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
62278569bb4SMatthew G. Knepley       ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr);
6236de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray));
62478569bb4SMatthew G. Knepley       ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr);
62578569bb4SMatthew G. Knepley       ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
62678569bb4SMatthew G. Knepley     }
62778569bb4SMatthew G. Knepley     ierr = ISDestroy(&compIS);CHKERRQ(ierr);
62878569bb4SMatthew G. Knepley   }
62978569bb4SMatthew G. Knepley   if (useNatural) {
63078569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr);
63178569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr);
63278569bb4SMatthew G. Knepley     ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);
63378569bb4SMatthew G. Knepley   }
63478569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
63578569bb4SMatthew G. Knepley }
63678569bb4SMatthew G. Knepley 
63778569bb4SMatthew G. Knepley /*
63878569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file
63978569bb4SMatthew G. Knepley 
64078569bb4SMatthew G. Knepley   Collective on v
64178569bb4SMatthew G. Knepley 
64278569bb4SMatthew G. Knepley   Input Parameters:
64378569bb4SMatthew G. Knepley + v  - The vector to be written
64478569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
64578569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1)
64678569bb4SMatthew G. Knepley 
64778569bb4SMatthew G. Knepley   Notes:
64878569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
64978569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
65078569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
65178569bb4SMatthew G. Knepley   amongst all zonal variables.
65278569bb4SMatthew G. Knepley 
65378569bb4SMatthew G. Knepley   Level: beginner
65478569bb4SMatthew G. Knepley 
65578569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Zonal_Internal()
65678569bb4SMatthew G. Knepley */
65778569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step)
65878569bb4SMatthew G. Knepley {
65978569bb4SMatthew G. Knepley   MPI_Comm          comm;
66078569bb4SMatthew G. Knepley   PetscMPIInt       size;
66178569bb4SMatthew G. Knepley   DM                dm;
66278569bb4SMatthew G. Knepley   Vec               vNatural, vComp;
66322a7544eSVaclav Hapla   const PetscScalar *varray;
66478569bb4SMatthew G. Knepley   const char       *vecname;
66578569bb4SMatthew G. Knepley   PetscInt          xs, xe, bs;
66678569bb4SMatthew G. Knepley   PetscBool         useNatural;
66778569bb4SMatthew G. Knepley   int               offset;
66878569bb4SMatthew G. Knepley   IS                compIS;
66978569bb4SMatthew G. Knepley   PetscInt         *csSize, *csID;
67078569bb4SMatthew G. Knepley   PetscInt          numCS, set, csxs = 0;
67178569bb4SMatthew G. Knepley   PetscErrorCode    ierr;
67278569bb4SMatthew G. Knepley 
67378569bb4SMatthew G. Knepley   PetscFunctionBegin;
67478569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)v, &comm);CHKERRQ(ierr);
67578569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
67678569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
67778569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
67878569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
67978569bb4SMatthew G. Knepley   if (useNatural) {
68078569bb4SMatthew G. Knepley     ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr);
68178569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr);
68278569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr);
68378569bb4SMatthew G. Knepley   } else {
68478569bb4SMatthew G. Knepley     vNatural = v;
68578569bb4SMatthew G. Knepley   }
68678569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
68778569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
68878569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr);
68978569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname);
69078569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
69178569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
69278569bb4SMatthew G. Knepley      We assume that they are stored sequentially
69378569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
69478569bb4SMatthew G. Knepley      but once the vector has been reordered to natural size, we cannot use the label informations
69578569bb4SMatthew G. Knepley      to figure out what to save where. */
69678569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
69778569bb4SMatthew G. Knepley   ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr);
6986de834b4SMatthew G. Knepley   PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID));
69978569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
70078569bb4SMatthew G. Knepley     ex_block block;
70178569bb4SMatthew G. Knepley 
70278569bb4SMatthew G. Knepley     block.id   = csID[set];
70378569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
7046de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_block_param,(exoid, &block));
70578569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
70678569bb4SMatthew G. Knepley   }
70778569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
70878569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
70978569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);}
71078569bb4SMatthew G. Knepley   for (set = 0; set < numCS; set++) {
71178569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
71278569bb4SMatthew G. Knepley 
71378569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
71478569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
71578569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
71678569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
71778569bb4SMatthew G. Knepley     if (bs == 1) {
71878569bb4SMatthew G. Knepley       ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr);
7196de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]));
72078569bb4SMatthew G. Knepley       ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr);
72178569bb4SMatthew G. Knepley     } else {
72278569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
72378569bb4SMatthew G. Knepley         ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
72478569bb4SMatthew G. Knepley         ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
72578569bb4SMatthew G. Knepley         ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr);
7266de834b4SMatthew G. Knepley         PetscStackCallStandard(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)]));
72778569bb4SMatthew G. Knepley         ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr);
72878569bb4SMatthew G. Knepley         ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
72978569bb4SMatthew G. Knepley       }
73078569bb4SMatthew G. Knepley     }
73178569bb4SMatthew G. Knepley     csxs += csSize[set];
73278569bb4SMatthew G. Knepley   }
73378569bb4SMatthew G. Knepley   ierr = PetscFree2(csID, csSize);CHKERRQ(ierr);
73478569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);}
73578569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMRestoreGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
73678569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
73778569bb4SMatthew G. Knepley }
73878569bb4SMatthew G. Knepley 
73978569bb4SMatthew G. Knepley /*
74078569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file
74178569bb4SMatthew G. Knepley 
74278569bb4SMatthew G. Knepley   Collective on v
74378569bb4SMatthew G. Knepley 
74478569bb4SMatthew G. Knepley   Input Parameters:
74578569bb4SMatthew G. Knepley + v  - The vector to be written
74678569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
74778569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1)
74878569bb4SMatthew G. Knepley 
74978569bb4SMatthew G. Knepley   Notes:
75078569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
75178569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
75278569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
75378569bb4SMatthew G. Knepley   amongst all zonal variables.
75478569bb4SMatthew G. Knepley 
75578569bb4SMatthew G. Knepley   Level: beginner
75678569bb4SMatthew G. Knepley 
75778569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Zonal_Internal()
75878569bb4SMatthew G. Knepley */
75978569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step)
76078569bb4SMatthew G. Knepley {
76178569bb4SMatthew G. Knepley   MPI_Comm          comm;
76278569bb4SMatthew G. Knepley   PetscMPIInt       size;
76378569bb4SMatthew G. Knepley   DM                dm;
76478569bb4SMatthew G. Knepley   Vec               vNatural, vComp;
76522a7544eSVaclav Hapla   PetscScalar      *varray;
76678569bb4SMatthew G. Knepley   const char       *vecname;
76778569bb4SMatthew G. Knepley   PetscInt          xs, xe, bs;
76878569bb4SMatthew G. Knepley   PetscBool         useNatural;
76978569bb4SMatthew G. Knepley   int               offset;
77078569bb4SMatthew G. Knepley   IS                compIS;
77178569bb4SMatthew G. Knepley   PetscInt         *csSize, *csID;
77278569bb4SMatthew G. Knepley   PetscInt          numCS, set, csxs = 0;
77378569bb4SMatthew G. Knepley   PetscErrorCode    ierr;
77478569bb4SMatthew G. Knepley 
77578569bb4SMatthew G. Knepley   PetscFunctionBegin;
77678569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
77778569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
77878569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
77978569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
78078569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
78178569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
78278569bb4SMatthew G. Knepley   else            {vNatural = v;}
78378569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
78478569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
78578569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr);
78678569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname);
78778569bb4SMatthew G. Knepley   /* Read local chunk of the result in the exodus file
78878569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
78978569bb4SMatthew G. Knepley      We assume that they are stored sequentially
79078569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
79178569bb4SMatthew G. Knepley      but once the vector has been reordered to natural size, we cannot use the label informations
79278569bb4SMatthew G. Knepley      to figure out what to save where. */
79378569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
79478569bb4SMatthew G. Knepley   ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr);
7956de834b4SMatthew G. Knepley   PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID));
79678569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
79778569bb4SMatthew G. Knepley     ex_block block;
79878569bb4SMatthew G. Knepley 
79978569bb4SMatthew G. Knepley     block.id   = csID[set];
80078569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
8016de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_block_param,(exoid, &block));
80278569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
80378569bb4SMatthew G. Knepley   }
80478569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
80578569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
80678569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);}
80778569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
80878569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
80978569bb4SMatthew G. Knepley 
81078569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
81178569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
81278569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
81378569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
81478569bb4SMatthew G. Knepley     if (bs == 1) {
81578569bb4SMatthew G. Knepley       ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr);
8166de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]));
81778569bb4SMatthew G. Knepley       ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr);
81878569bb4SMatthew G. Knepley     } else {
81978569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
82078569bb4SMatthew G. Knepley         ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
82178569bb4SMatthew G. Knepley         ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
82278569bb4SMatthew G. Knepley         ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr);
8236de834b4SMatthew G. Knepley         PetscStackCallStandard(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)]));
82478569bb4SMatthew G. Knepley         ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr);
82578569bb4SMatthew G. Knepley         ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
82678569bb4SMatthew G. Knepley       }
82778569bb4SMatthew G. Knepley     }
82878569bb4SMatthew G. Knepley     csxs += csSize[set];
82978569bb4SMatthew G. Knepley   }
83078569bb4SMatthew G. Knepley   ierr = PetscFree2(csID, csSize);CHKERRQ(ierr);
83178569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);}
83278569bb4SMatthew G. Knepley   if (useNatural) {
83378569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr);
83478569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr);
83578569bb4SMatthew G. Knepley     ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);
83678569bb4SMatthew G. Knepley   }
83778569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
83878569bb4SMatthew G. Knepley }
839b53c8456SSatish Balay #endif
84078569bb4SMatthew G. Knepley 
84133751fbdSMichael Lange /*@C
842eaf898f9SPatrick Sanan   DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file.
84333751fbdSMichael Lange 
844d083f849SBarry Smith   Collective
84533751fbdSMichael Lange 
84633751fbdSMichael Lange   Input Parameters:
84733751fbdSMichael Lange + comm  - The MPI communicator
84833751fbdSMichael Lange . filename - The name of the ExodusII file
84933751fbdSMichael Lange - interpolate - Create faces and edges in the mesh
85033751fbdSMichael Lange 
85133751fbdSMichael Lange   Output Parameter:
85233751fbdSMichael Lange . dm  - The DM object representing the mesh
85333751fbdSMichael Lange 
85433751fbdSMichael Lange   Level: beginner
85533751fbdSMichael Lange 
85633751fbdSMichael Lange .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus()
85733751fbdSMichael Lange @*/
85833751fbdSMichael Lange PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
85933751fbdSMichael Lange {
86033751fbdSMichael Lange   PetscMPIInt    rank;
86133751fbdSMichael Lange   PetscErrorCode ierr;
86233751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
863e84f7031SBlaise Bourdin   int   CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1;
86433751fbdSMichael Lange   float version;
86533751fbdSMichael Lange #endif
86633751fbdSMichael Lange 
86733751fbdSMichael Lange   PetscFunctionBegin;
86833751fbdSMichael Lange   PetscValidCharPointer(filename, 2);
86933751fbdSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
87033751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
87133751fbdSMichael Lange   if (!rank) {
87233751fbdSMichael Lange     exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version);
87333751fbdSMichael Lange     if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename);
87433751fbdSMichael Lange   }
87533751fbdSMichael Lange   ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr);
8766de834b4SMatthew G. Knepley   if (!rank) {PetscStackCallStandard(ex_close,(exoid));}
87733751fbdSMichael Lange #else
87833751fbdSMichael Lange   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
87933751fbdSMichael Lange #endif
88033751fbdSMichael Lange   PetscFunctionReturn(0);
88133751fbdSMichael Lange }
88233751fbdSMichael Lange 
883552f7358SJed Brown /*@
88433751fbdSMichael Lange   DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID.
885552f7358SJed Brown 
886d083f849SBarry Smith   Collective
887552f7358SJed Brown 
888552f7358SJed Brown   Input Parameters:
889552f7358SJed Brown + comm  - The MPI communicator
890552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
891552f7358SJed Brown - interpolate - Create faces and edges in the mesh
892552f7358SJed Brown 
893552f7358SJed Brown   Output Parameter:
894552f7358SJed Brown . dm  - The DM object representing the mesh
895552f7358SJed Brown 
896552f7358SJed Brown   Level: beginner
897552f7358SJed Brown 
898e84e7769SJed Brown .seealso: DMPLEX, DMCreate()
899552f7358SJed Brown @*/
900552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm)
901552f7358SJed Brown {
902552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
903552f7358SJed Brown   PetscMPIInt    num_proc, rank;
904552f7358SJed Brown   PetscSection   coordSection;
905552f7358SJed Brown   Vec            coordinates;
906552f7358SJed Brown   PetscScalar    *coords;
907552f7358SJed Brown   PetscInt       coordSize, v;
908552f7358SJed Brown   PetscErrorCode ierr;
909552f7358SJed Brown   /* Read from ex_get_init() */
910552f7358SJed Brown   char title[PETSC_MAX_PATH_LEN+1];
911*e8f6893fSMatthew G. Knepley   int  dim    = 0, numVertices = 0, numCells = 0, numHybridCells = 0;
912552f7358SJed Brown   int  num_cs = 0, num_vs = 0, num_fs = 0;
913552f7358SJed Brown #endif
914552f7358SJed Brown 
915552f7358SJed Brown   PetscFunctionBegin;
916552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
917552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
918552f7358SJed Brown   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
919552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
920552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
921552f7358SJed Brown   /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */
922552f7358SJed Brown   if (!rank) {
923580bdb30SBarry Smith     ierr = PetscMemzero(title,PETSC_MAX_PATH_LEN+1);CHKERRQ(ierr);
9246de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_init,(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs));
925552f7358SJed Brown     if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n");
926552f7358SJed Brown   }
927552f7358SJed Brown   ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
928552f7358SJed Brown   ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
929552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr);
930c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
931552f7358SJed Brown   ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
932552f7358SJed Brown 
933552f7358SJed Brown   /* Read cell sets information */
934552f7358SJed Brown   if (!rank) {
935552f7358SJed Brown     PetscInt *cone;
936*e8f6893fSMatthew G. Knepley     int      c, cs, ncs, c_loc, v, v_loc;
937552f7358SJed Brown     /* Read from ex_get_elem_blk_ids() */
938*e8f6893fSMatthew G. Knepley     int *cs_id, *cs_order;
939552f7358SJed Brown     /* Read from ex_get_elem_block() */
940552f7358SJed Brown     char buffer[PETSC_MAX_PATH_LEN+1];
941*e8f6893fSMatthew G. Knepley     int  num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr;
942552f7358SJed Brown     /* Read from ex_get_elem_conn() */
943552f7358SJed Brown     int *cs_connect;
944552f7358SJed Brown 
945552f7358SJed Brown     /* Get cell sets IDs */
946*e8f6893fSMatthew G. Knepley     ierr = PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order);CHKERRQ(ierr);
9476de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, cs_id));
948552f7358SJed Brown     /* Read the cell set connectivity table and build mesh topology
949552f7358SJed Brown        EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
950*e8f6893fSMatthew G. Knepley     /* Check for a hybrid mesh */
951*e8f6893fSMatthew G. Knepley     for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) {
952*e8f6893fSMatthew G. Knepley       PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr));
953*e8f6893fSMatthew G. Knepley       switch (dim) {
954*e8f6893fSMatthew G. Knepley         case 3:
955*e8f6893fSMatthew G. Knepley         switch (num_vertex_per_cell) {
956*e8f6893fSMatthew G. Knepley           case 6:
957*e8f6893fSMatthew G. Knepley             cs_order[cs] = cs;
958*e8f6893fSMatthew G. Knepley             numHybridCells += num_cell_in_set;
959*e8f6893fSMatthew G. Knepley             ++num_hybrid;
960*e8f6893fSMatthew G. Knepley           break;
961*e8f6893fSMatthew G. Knepley           default:
962*e8f6893fSMatthew G. Knepley             for (c = cs; c > cs-num_hybrid; --c) cs_order[c] = cs_order[c-1];
963*e8f6893fSMatthew G. Knepley             cs_order[cs-num_hybrid] = cs;
964*e8f6893fSMatthew G. Knepley         }
965*e8f6893fSMatthew G. Knepley         break;
966*e8f6893fSMatthew G. Knepley       default:
967*e8f6893fSMatthew G. Knepley         for (c = cs; c > cs-num_hybrid; --c) cs_order[c] = cs_order[c-1];
968*e8f6893fSMatthew G. Knepley         cs_order[cs-num_hybrid] = cs;
969*e8f6893fSMatthew G. Knepley       }
970*e8f6893fSMatthew G. Knepley     }
971*e8f6893fSMatthew G. Knepley     if (num_hybrid) {ierr = DMPlexSetHybridBounds(*dm, numCells-numHybridCells, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);}
972552f7358SJed Brown     /* First set sizes */
973*e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
974*e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
9756de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr));
976552f7358SJed Brown       for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
977552f7358SJed Brown         ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr);
978552f7358SJed Brown       }
979552f7358SJed Brown     }
980552f7358SJed Brown     ierr = DMSetUp(*dm);CHKERRQ(ierr);
981*e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
982*e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
9836de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr));
984dcca6d9dSJed Brown       ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr);
9856de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_conn,(exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL));
986eaf898f9SPatrick Sanan       /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
987552f7358SJed Brown       for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
988552f7358SJed Brown         for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) {
989552f7358SJed Brown           cone[v_loc] = cs_connect[v]+numCells-1;
990552f7358SJed Brown         }
991731dcdf4SMatthew G. Knepley         if (dim == 3) {
9922e1b13c2SMatthew G. Knepley           /* Tetrahedra are inverted */
9932e1b13c2SMatthew G. Knepley           if (num_vertex_per_cell == 4) {
9942e1b13c2SMatthew G. Knepley             PetscInt tmp = cone[0];
9952e1b13c2SMatthew G. Knepley             cone[0] = cone[1];
9962e1b13c2SMatthew G. Knepley             cone[1] = tmp;
9972e1b13c2SMatthew G. Knepley           }
9982e1b13c2SMatthew G. Knepley           /* Hexahedra are inverted */
9992e1b13c2SMatthew G. Knepley           if (num_vertex_per_cell == 8) {
10002e1b13c2SMatthew G. Knepley             PetscInt tmp = cone[1];
10012e1b13c2SMatthew G. Knepley             cone[1] = cone[3];
10022e1b13c2SMatthew G. Knepley             cone[3] = tmp;
10032e1b13c2SMatthew G. Knepley           }
1004731dcdf4SMatthew G. Knepley         }
1005552f7358SJed Brown         ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
1006c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr);
1007552f7358SJed Brown       }
1008552f7358SJed Brown       ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr);
1009552f7358SJed Brown     }
1010*e8f6893fSMatthew G. Knepley     ierr = PetscFree2(cs_id, cs_order);CHKERRQ(ierr);
1011552f7358SJed Brown   }
1012552f7358SJed Brown   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
1013552f7358SJed Brown   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
1014552f7358SJed Brown   if (interpolate) {
10155fd9971aSMatthew G. Knepley     DM idm;
1016552f7358SJed Brown 
10179f074e33SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
1018552f7358SJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
1019552f7358SJed Brown     *dm  = idm;
1020552f7358SJed Brown   }
1021552f7358SJed Brown 
1022552f7358SJed Brown   /* Create vertex set label */
1023552f7358SJed Brown   if (!rank && (num_vs > 0)) {
1024552f7358SJed Brown     int vs, v;
1025552f7358SJed Brown     /* Read from ex_get_node_set_ids() */
1026552f7358SJed Brown     int *vs_id;
1027552f7358SJed Brown     /* Read from ex_get_node_set_param() */
1028f958083aSBlaise Bourdin     int num_vertex_in_set;
1029552f7358SJed Brown     /* Read from ex_get_node_set() */
1030552f7358SJed Brown     int *vs_vertex_list;
1031552f7358SJed Brown 
1032552f7358SJed Brown     /* Get vertex set ids */
1033785e854fSJed Brown     ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr);
10346de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_ids,(exoid, EX_NODE_SET, vs_id));
1035552f7358SJed Brown     for (vs = 0; vs < num_vs; ++vs) {
10366de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_set_param,(exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL));
1037785e854fSJed Brown       ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr);
10386de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_set,(exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL));
1039552f7358SJed Brown       for (v = 0; v < num_vertex_in_set; ++v) {
1040c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr);
1041552f7358SJed Brown       }
1042552f7358SJed Brown       ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr);
1043552f7358SJed Brown     }
1044552f7358SJed Brown     ierr = PetscFree(vs_id);CHKERRQ(ierr);
1045552f7358SJed Brown   }
1046552f7358SJed Brown   /* Read coordinates */
104769d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
1048552f7358SJed Brown   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
1049552f7358SJed Brown   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
1050552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
1051552f7358SJed Brown   for (v = numCells; v < numCells+numVertices; ++v) {
1052552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
1053552f7358SJed Brown     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
1054552f7358SJed Brown   }
1055552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1056552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
10578b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
1058552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1059552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
10604e90ef8eSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
10612eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
1062552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
10630aba08f6SKarl Rupp   if (!rank) {
1064e84f7031SBlaise Bourdin     PetscReal *x, *y, *z;
1065552f7358SJed Brown 
1066dcca6d9dSJed Brown     ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr);
10676de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_coord,(exoid, x, y, z));
10680d644c17SKarl Rupp     if (dim > 0) {
10690d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v];
10700d644c17SKarl Rupp     }
10710d644c17SKarl Rupp     if (dim > 1) {
10720d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v];
10730d644c17SKarl Rupp     }
10740d644c17SKarl Rupp     if (dim > 2) {
10750d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v];
10760d644c17SKarl Rupp     }
1077552f7358SJed Brown     ierr = PetscFree3(x,y,z);CHKERRQ(ierr);
1078552f7358SJed Brown   }
1079552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1080552f7358SJed Brown   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
1081552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1082552f7358SJed Brown 
1083552f7358SJed Brown   /* Create side set label */
1084552f7358SJed Brown   if (!rank && interpolate && (num_fs > 0)) {
1085552f7358SJed Brown     int fs, f, voff;
1086552f7358SJed Brown     /* Read from ex_get_side_set_ids() */
1087552f7358SJed Brown     int *fs_id;
1088552f7358SJed Brown     /* Read from ex_get_side_set_param() */
1089f958083aSBlaise Bourdin     int num_side_in_set;
1090552f7358SJed Brown     /* Read from ex_get_side_set_node_list() */
1091552f7358SJed Brown     int *fs_vertex_count_list, *fs_vertex_list;
1092ef073283Smichael_afanasiev     /* Read side set labels */
1093ef073283Smichael_afanasiev     char fs_name[MAX_STR_LENGTH+1];
1094552f7358SJed Brown 
1095552f7358SJed Brown     /* Get side set ids */
1096785e854fSJed Brown     ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr);
10976de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_ids,(exoid, EX_SIDE_SET, fs_id));
1098552f7358SJed Brown     for (fs = 0; fs < num_fs; ++fs) {
10996de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_set_param,(exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL));
1100dcca6d9dSJed Brown       ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr);
11016de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_side_set_node_list,(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list));
1102ef073283Smichael_afanasiev       /* Get the specific name associated with this side set ID. */
1103ef073283Smichael_afanasiev       int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name);
1104552f7358SJed Brown       for (f = 0, voff = 0; f < num_side_in_set; ++f) {
11050298fd71SBarry Smith         const PetscInt *faces   = NULL;
1106552f7358SJed Brown         PetscInt       faceSize = fs_vertex_count_list[f], numFaces;
1107552f7358SJed Brown         PetscInt       faceVertices[4], v;
1108552f7358SJed Brown 
1109552f7358SJed Brown         if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize);
1110552f7358SJed Brown         for (v = 0; v < faceSize; ++v, ++voff) {
1111552f7358SJed Brown           faceVertices[v] = fs_vertex_list[voff]+numCells-1;
1112552f7358SJed Brown         }
1113552f7358SJed Brown         ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
1114552f7358SJed Brown         if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces);
1115c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr);
1116ef073283Smichael_afanasiev         /* Only add the label if one has been detected for this side set. */
1117ef073283Smichael_afanasiev         if (!fs_name_err) {
1118ef073283Smichael_afanasiev           ierr = DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]);CHKERRQ(ierr);
1119ef073283Smichael_afanasiev         }
1120552f7358SJed Brown         ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
1121552f7358SJed Brown       }
1122552f7358SJed Brown       ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr);
1123552f7358SJed Brown     }
1124552f7358SJed Brown     ierr = PetscFree(fs_id);CHKERRQ(ierr);
1125552f7358SJed Brown   }
1126552f7358SJed Brown #else
1127552f7358SJed Brown   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
1128552f7358SJed Brown #endif
1129552f7358SJed Brown   PetscFunctionReturn(0);
1130552f7358SJed Brown }
1131