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