xref: /petsc/src/dm/impls/plex/plexexodusii.c (revision e02669251c099ae2faea4c9cee5397e8f99c2eb2)
1 #define PETSCDM_DLL
2 #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3 
4 #if defined(PETSC_HAVE_EXODUSII)
5 #include <netcdf.h>
6 #include <exodusII.h>
7 #endif
8 
9 #if defined(PETSC_HAVE_EXODUSII)
10 /*
11   EXOGetVarIndex - Locate a result in an exodus file based on its name
12 
13   Not collective
14 
15   Input Parameters:
16 + exoid    - the exodus id of a file (obtained from ex_open or ex_create for instance)
17 . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK
18 - name     - the name of the result
19 
20   Output Parameters:
21 . varIndex - the location in the exodus file of the result
22 
23   Notes:
24   The exodus variable index is obtained by comparing name and the
25   names of zonal variables declared in the exodus file. For instance if name is "V"
26   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
27   amongst all variables of type obj_type.
28 
29   Level: beginner
30 
31 .keywords: mesh,ExodusII
32 .seealso: EXOGetVarIndex(),DMView_Exodus(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadNodal_PlexEXO(),VecLoadZonal_PlexEXO()
33 */
34 static PetscErrorCode EXOGetVarIndex_Private(int exoid, ex_entity_type obj_type, const char name[], int *varIndex)
35 {
36   int            exoerr, num_vars, i, j;
37   char           ext_name[MAX_STR_LENGTH+1], var_name[MAX_STR_LENGTH+1];
38   const int      num_suffix = 5;
39   char          *suffix[5];
40   PetscBool      flg;
41   PetscErrorCode ierr;
42 
43   PetscFunctionBegin;
44   suffix[0] = (char *) "";
45   suffix[1] = (char *) "_X";
46   suffix[2] = (char *) "_XX";
47   suffix[3] = (char *) "_1";
48   suffix[4] = (char *) "_11";
49 
50   *varIndex = 0;
51   exoerr = ex_get_variable_param(exoid, obj_type, &num_vars);
52   for (i = 0; i < num_vars; ++i) {
53     exoerr = ex_get_variable_name(exoid, obj_type, i+1, var_name);
54     for (j = 0; j < num_suffix; ++j){
55       ierr = PetscStrncpy(ext_name, name, MAX_STR_LENGTH);CHKERRQ(ierr);
56       ierr = PetscStrncat(ext_name, suffix[j], MAX_STR_LENGTH);CHKERRQ(ierr);
57       ierr = PetscStrcasecmp(ext_name, var_name, &flg);
58       if (flg) {
59         *varIndex = i+1;
60         PetscFunctionReturn(0);
61       }
62     }
63   }
64   SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to locate variable %s in ExodusII file.", name);
65  PetscFunctionReturn(-1);
66 }
67 #endif
68 
69 /*
70   DMPlexView_ExodusII_Internal - Write a DM to disk in exodus format
71 
72   Collective on dm
73 
74   Input Parameters:
75 + dm  - The dm to be written
76 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
77 - degree - the degree of the interpolation space
78 
79   Notes:
80   Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels)
81   consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted.
82 
83   If the dm has been distributed and exoid points to different files on each MPI rank, only the "local" part of the DM
84   (including "ghost" cells and vertices) will be written. If all exoid points to the same file, the resulting file will
85   probably be corrupted.
86 
87   DMPlex only represents geometry while most post-processing software expect that a mesh also provides information
88   on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2.
89   It should be extended to use PetscFE objects.
90 
91   This function will only handle TRI, TET, QUAD and HEX cells.
92   Level: beginner
93 
94 .keywords: mesh, ExodusII
95 .seealso: EXOGetVarIndex_Private(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadNodal_PlexEXO(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO()
96 */
97 PetscErrorCode DMPlexView_ExodusII_Internal(DM dm, int exoid, PetscInt degree)
98 {
99   enum ElemType {TRI, QUAD, TET, HEX};
100   MPI_Comm        comm;
101   /* Connectivity Variables */
102   PetscInt        cellsNotInConnectivity;
103   /* Cell Sets */
104   DMLabel         csLabel;
105   IS              csIS;
106   const PetscInt *csIdx;
107   PetscInt        num_cs, cs;
108   enum ElemType  *type;
109   /* Coordinate Variables */
110   DM              cdm;
111   PetscSection    section;
112   Vec             coord;
113   PetscInt      **nodes;
114   PetscInt        depth, d, dim, skipCells = 0;
115   PetscInt        pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes;
116   PetscInt        num_vs, num_fs;
117   PetscMPIInt     rank, size;
118   const char     *dmName;
119   PetscErrorCode  ierr;
120 
121   PetscFunctionBegin;
122   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
123   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
124   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
125   /* --- Get DM info --- */
126   ierr = PetscObjectGetName((PetscObject) dm, &dmName);CHKERRQ(ierr);
127   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
128   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
129   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
130   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
131   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
132   ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);
133   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
134   numCells    = cEnd - cStart;
135   numEdges    = eEnd - eStart;
136   numVertices = vEnd - vStart;
137   if (depth == 3) {numFaces = fEnd - fStart;}
138   else            {numFaces = 0;}
139   ierr = DMGetLabelSize(dm, "Cell Sets", &num_cs);CHKERRQ(ierr);
140   ierr = DMGetLabelSize(dm, "Vertex Sets", &num_vs);CHKERRQ(ierr);
141   ierr = DMGetLabelSize(dm, "Face Sets", &num_fs);CHKERRQ(ierr);
142   ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr);
143   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
144   if (num_cs > 0) {
145     ierr = DMGetLabel(dm, "Cell Sets", &csLabel);CHKERRQ(ierr);
146     ierr = DMLabelGetValueIS(csLabel, &csIS);CHKERRQ(ierr);
147     ierr = ISGetIndices(csIS, &csIdx);CHKERRQ(ierr);
148   }
149   ierr = PetscMalloc1(num_cs, &nodes);CHKERRQ(ierr);
150   /* Set element type for each block and compute total number of nodes */
151   ierr = PetscMalloc1(num_cs, &type);CHKERRQ(ierr);
152   numNodes = numVertices;
153   if (degree == 2) {numNodes += numEdges;}
154   cellsNotInConnectivity = numCells;
155   for (cs = 0; cs < num_cs; ++cs) {
156     IS              stratumIS;
157     const PetscInt *cells;
158     PetscScalar    *xyz = NULL;
159     PetscInt        csSize, closureSize;
160     PetscInt        nodesTriP1[4]  = {3,0,0,0};
161     PetscInt        nodesTriP2[4]  = {3,3,0,0};
162     PetscInt        nodesQuadP1[4] = {4,0,0,0};
163     PetscInt        nodesQuadP2[4] = {4,4,0,1};
164     PetscInt        nodesTetP1[4]  = {4,0,0,0};
165     PetscInt        nodesTetP2[4]  = {4,6,0,0};
166     PetscInt        nodesHexP1[4]  = {8,0,0,0};
167     PetscInt        nodesHexP2[4]  = {8,12,6,1};
168 
169     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
170     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
171     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
172     ierr = DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr);
173     switch (dim) {
174     case 2:
175       if      (closureSize == 3*dim) {type[cs] = TRI;}
176       else if (closureSize == 4*dim) {type[cs] = QUAD;}
177       else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim);
178       break;
179     case 3:
180       if      (closureSize == 4*dim) {type[cs] = TET;}
181       else if (closureSize == 8*dim) {type[cs] = HEX;}
182       else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim);
183       break;
184     default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not handled by ExodusII viewer", dim);
185     }
186     if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;}
187     if ((degree == 2) && (type[cs] == HEX))  {numNodes += csSize; numNodes += numFaces;}
188     ierr = DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr);
189     /* Set nodes and Element type */
190     if (type[cs] == TRI) {
191       if      (degree == 1) nodes[cs] = nodesTriP1;
192       else if (degree == 2) nodes[cs] = nodesTriP2;
193     } else if (type[cs] == QUAD) {
194       if      (degree == 1) nodes[cs] = nodesQuadP1;
195       else if (degree == 2) nodes[cs] = nodesQuadP2;
196     } else if (type[cs] == TET) {
197       if      (degree == 1) nodes[cs] = nodesTetP1;
198       else if (degree == 2) nodes[cs] = nodesTetP2;
199     } else if (type[cs] == HEX) {
200       if      (degree == 1) nodes[cs] = nodesHexP1;
201       else if (degree == 2) nodes[cs] = nodesHexP2;
202     }
203     /* Compute the number of cells not in the connectivity table */
204     cellsNotInConnectivity -= nodes[cs][3]*csSize;
205 
206     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
207     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
208   }
209   if (num_cs > 0) {ierr = ex_put_init(exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs);CHKERRQ(ierr);}
210   /* --- Connectivity --- */
211   for (cs = 0; cs < num_cs; ++cs) {
212     IS              stratumIS;
213     const PetscInt *cells;
214     PetscInt       *connect;
215     PetscInt        edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0;
216     PetscInt        csSize, c, connectSize, closureSize;
217     char           *elem_type = NULL;
218     char            elem_type_tri3[]  = "TRI3",  elem_type_quad4[] = "QUAD4";
219     char            elem_type_tri6[]  = "TRI6",  elem_type_quad9[] = "QUAD9";
220     char            elem_type_tet4[]  = "TET4",  elem_type_hex8[]  = "HEX8";
221     char            elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27";
222 
223     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
224     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
225     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
226     /* Set Element type */
227     if (type[cs] == TRI) {
228       if      (degree == 1) elem_type = elem_type_tri3;
229       else if (degree == 2) elem_type = elem_type_tri6;
230     } else if (type[cs] == QUAD) {
231       if      (degree == 1) elem_type = elem_type_quad4;
232       else if (degree == 2) elem_type = elem_type_quad9;
233     } else if (type[cs] == TET) {
234       if      (degree == 1) elem_type = elem_type_tet4;
235       else if (degree == 2) elem_type = elem_type_tet10;
236     } else if (type[cs] == HEX) {
237       if      (degree == 1) elem_type = elem_type_hex8;
238       else if (degree == 2) elem_type = elem_type_hex27;
239     }
240     connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3];
241     ierr = PetscMalloc1(connectSize, &connect);CHKERRQ(ierr);
242     ierr = ex_put_block(exoid, EX_ELEM_BLOCK, cs, elem_type, csSize, connectSize, 0, 0, 1);
243     /* Find number of vertices, edges, and faces in the closure */
244     verticesInClosure = nodes[cs][0];
245     if (depth > 1) {
246       if (dim == 2) {
247         ierr = DMPlexGetConeSize(dm, cells[0], &edgesInClosure);CHKERRQ(ierr);
248       } else if (dim == 3) {
249         PetscInt *closure = NULL;
250 
251         ierr = DMPlexGetConeSize(dm, cells[0], &facesInClosure);CHKERRQ(ierr);
252         ierr = DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
253         edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure;
254         ierr = DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
255       }
256     }
257     /* Get connectivity for each cell */
258     for (c = 0; c < csSize; ++c) {
259       PetscInt *closure = NULL;
260       PetscInt  temp, i;
261 
262       ierr = DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
263       for (i = 0; i < connectSize; ++i) {
264         if (i < nodes[cs][0]) {/* Vertices */
265           connect[i] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1;
266           connect[i] -= cellsNotInConnectivity;
267         } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */
268           connect[i] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1;
269           if (nodes[cs][2] == 0) connect[i] -= numFaces;
270           connect[i] -= cellsNotInConnectivity;
271         } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]){ /* Cells */
272           connect[i] = closure[0] + 1;
273           connect[i] -= skipCells;
274         } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]){ /* Faces */
275           connect[i] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1;
276           connect[i] -= cellsNotInConnectivity;
277         } else {
278           connect[i] = -1;
279         }
280       }
281       /* Tetrahedra are inverted */
282       if (type[cs] == TET) {
283         temp = connect[0]; connect[0] = connect[1]; connect[1] = temp;
284         if (degree == 2) {
285           temp = connect[5]; connect[5] = connect[6]; connect[6] = temp;
286           temp = connect[7]; connect[7] = connect[8]; connect[8] = temp;
287         }
288       }
289       /* Hexahedra are inverted */
290       if (type[cs] == HEX) {
291         temp = connect[1]; connect[1] = connect[3]; connect[3] = temp;
292         if (degree == 2) {
293           temp = connect[8];  connect[8]  = connect[11]; connect[11] = temp;
294           temp = connect[9];  connect[9]  = connect[10]; connect[10] = temp;
295           temp = connect[16]; connect[16] = connect[17]; connect[17] = temp;
296           temp = connect[18]; connect[18] = connect[19]; connect[19] = temp;
297 
298           temp = connect[12]; connect[12] = connect[16]; connect[16] = temp;
299           temp = connect[13]; connect[13] = connect[17]; connect[17] = temp;
300           temp = connect[14]; connect[14] = connect[18]; connect[18] = temp;
301           temp = connect[15]; connect[15] = connect[19]; connect[19] = temp;
302 
303           temp = connect[23]; connect[23] = connect[26]; connect[26] = temp;
304           temp = connect[24]; connect[24] = connect[25]; connect[25] = temp;
305           temp = connect[25]; connect[25] = connect[26]; connect[26] = temp;
306         }
307       }
308       ierr = ex_put_partial_conn(exoid, EX_ELEM_BLOCK, cs, c+1, 1, connect, NULL, NULL);
309       ierr = DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
310     }
311     skipCells += (nodes[cs][3] == 0)*csSize;
312     ierr = PetscFree(connect);CHKERRQ(ierr);
313     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
314     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
315   }
316   ierr = PetscFree(type);CHKERRQ(ierr);
317   /* --- Coordinates --- */
318   ierr = PetscSectionCreate(comm, &section);CHKERRQ(ierr);
319   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
320   for (d = 0; d < depth; ++d) {
321     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
322     for (p = pStart; p < pEnd; ++p) {
323       ierr = PetscSectionSetDof(section, p, nodes[0][d] > 0);CHKERRQ(ierr);
324     }
325   }
326   for (cs = 0; cs < num_cs; ++cs) {
327     IS              stratumIS;
328     const PetscInt *cells;
329     PetscInt        csSize, c;
330 
331     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
332     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
333     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
334     for (c = 0; c < csSize; ++c) {
335       ierr = PetscSectionSetDof(section, cells[c], nodes[cs][3] > 0);CHKERRQ(ierr);
336     }
337     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
338     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
339   }
340   if (num_cs > 0) {
341     ierr = ISRestoreIndices(csIS, &csIdx);CHKERRQ(ierr);
342     ierr = ISDestroy(&csIS);CHKERRQ(ierr);
343   }
344   ierr = PetscFree(nodes);CHKERRQ(ierr);
345   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
346   if (numNodes > 0) {
347     const char  *coordNames[3] = {"x", "y", "z"};
348     PetscScalar *coords, *closure;
349     PetscReal   *cval;
350     PetscInt     hasDof, n = 0;
351 
352     /* There can't be more than 24 values in the closure of a point for the coord section */
353     ierr = PetscMalloc3(numNodes*3, &coords, dim, &cval, 24, &closure);CHKERRQ(ierr);
354     ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr);
355     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
356     for (p = pStart; p < pEnd; ++p) {
357       ierr = PetscSectionGetDof(section, p, &hasDof);
358       if (hasDof) {
359         PetscInt closureSize = 24, j;
360 
361         ierr = DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure);CHKERRQ(ierr);
362         for (d = 0; d < dim; ++d) {
363           cval[d] = 0.0;
364           for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d];
365           coords[d*numNodes+n] = cval[d] * dim / closureSize;
366         }
367         ++n;
368       }
369     }
370     ierr = ex_put_coord(exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes]);CHKERRQ(ierr);
371     ierr = PetscFree3(coords, cval, closure);CHKERRQ(ierr);
372     ierr = ex_put_coord_names(exoid, (char **) coordNames);CHKERRQ(ierr);
373   }
374   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
375   PetscFunctionReturn(0);
376 }
377 
378 /*
379   VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file
380 
381   Collective on v
382 
383   Input Parameters:
384 + v  - The vector to be written
385 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
386 - step - the time step to write at (exodus steps are numbered starting from 1)
387 
388   Notes:
389   The exodus result nodal variable index is obtained by comparing the Vec name and the
390   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
391   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
392   amongst all nodal variables.
393 
394   Level: beginner
395 
396 .keywords: mesh,ExodusII
397 .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecLoadNodal_PlexEXO(),VecViewZonal_PlexEXO(),VecLoadZonal_PlexEXO()
398 @*/
399 PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step)
400 {
401   MPI_Comm         comm;
402   PetscMPIInt      size;
403 #if defined(PETSC_HAVE_EXODUSII)
404   DM               dm;
405   Vec              vNatural, vComp;
406   const PetscReal *varray;
407   const char      *vecname;
408   PetscInt         xs, xe, bs;
409   PetscBool        useNatural;
410   int              offset;
411   PetscErrorCode   ierr;
412 #endif
413 
414   PetscFunctionBegin;
415   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
416   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
417 #if defined(PETSC_HAVE_EXODUSII)
418   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
419   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
420   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
421   if (useNatural) {
422     ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr);
423     ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr);
424     ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr);
425   } else {
426     vNatural = v;
427   }
428   /* Get the location of the variable in the exodus file */
429   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
430   ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr);
431   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file\n", vecname);
432   /* Write local chunk of the result in the exodus file
433      exodus stores each component of a vector-valued field as a separate variable.
434      We assume that they are stored sequentially */
435   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
436   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
437   if (bs == 1) {
438     ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr);
439     ierr = ex_put_partial_var(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);CHKERRQ(ierr);
440     ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr);
441   } else {
442     IS       compIS;
443     PetscInt c;
444 
445     ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);
446     for (c = 0; c < bs; ++c) {
447       ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
448       ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
449       ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr);
450       ierr = ex_put_partial_var(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);CHKERRQ(ierr);
451       ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr);
452       ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
453     }
454     ierr = ISDestroy(&compIS);CHKERRQ(ierr);
455   }
456   if (useNatural) {ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);}
457 #else
458   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
459 #endif
460   PetscFunctionReturn(0);
461 }
462 
463 /*
464   VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file
465 
466   Collective on v
467 
468   Input Parameters:
469 + v  - The vector to be written
470 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
471 - step - the time step to read at (exodus steps are numbered starting from 1)
472 
473   Notes:
474   The exodus result nodal variable index is obtained by comparing the Vec name and the
475   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
476   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
477   amongst all nodal variables.
478 
479   Level: beginner
480 
481 .keywords: mesh, ExodusII
482 .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO()
483 */
484 PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step)
485 {
486   MPI_Comm       comm;
487   PetscMPIInt    size;
488 #if defined(PETSC_HAVE_EXODUSII)
489   DM             dm;
490   Vec            vNatural, vComp;
491   PetscReal     *varray;
492   const char    *vecname;
493   PetscInt       xs, xe, bs;
494   PetscBool      useNatural;
495   int            offset;
496   PetscErrorCode ierr;
497 #endif
498 
499   PetscFunctionBegin;
500   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
501   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
502 #if defined(PETSC_HAVE_EXODUSII)
503   ierr = VecGetDM(v,&dm);CHKERRQ(ierr);
504   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
505   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
506   if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
507   else            {vNatural = v;}
508   /* Get the location of the variable in the exodus file */
509   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
510   ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr);
511   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file", vecname);
512   /* Read local chunk from the file */
513   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
514   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
515   if (bs == 1) {
516     ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr);
517     ierr = ex_get_partial_var(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);CHKERRQ(ierr);
518     ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr);
519   } else {
520     IS       compIS;
521     PetscInt c;
522 
523     ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);
524     for (c = 0; c < bs; ++c) {
525       ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
526       ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
527       ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr);
528       ierr = ex_get_partial_var(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);CHKERRQ(ierr);
529       ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr);
530       ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
531     }
532     ierr = ISDestroy(&compIS);CHKERRQ(ierr);
533   }
534   if (useNatural) {
535     ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr);
536     ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr);
537     ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);
538   }
539 #else
540   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
541 #endif
542   PetscFunctionReturn(0);
543 }
544 
545 /*
546   VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file
547 
548   Collective on v
549 
550   Input Parameters:
551 + v  - The vector to be written
552 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
553 - step - the time step to write at (exodus steps are numbered starting from 1)
554 
555   Notes:
556   The exodus result zonal variable index is obtained by comparing the Vec name and the
557   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
558   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
559   amongst all zonal variables.
560 
561   Level: beginner
562 
563 .keywords: mesh,ExodusII
564 .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Zonal_Internal()
565 */
566 PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step)
567 {
568   MPI_Comm          comm;
569   PetscMPIInt       size;
570 #if defined(PETSC_HAVE_EXODUSII)
571   DM                dm;
572   Vec               vNatural, vComp;
573   const PetscReal  *varray;
574   const char       *vecname;
575   PetscInt          xs, xe, bs;
576   PetscBool         useNatural;
577   int               offset;
578   IS                compIS;
579   PetscInt         *csSize, *csID;
580   PetscInt          numCS, set, csxs = 0;
581   PetscErrorCode    ierr;
582 #endif
583 
584   PetscFunctionBegin;
585   ierr = PetscObjectGetComm((PetscObject)v, &comm);CHKERRQ(ierr);
586   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
587 #if defined(PETSC_HAVE_EXODUSII)
588   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
589   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
590   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
591   if (useNatural) {
592     ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr);
593     ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr);
594     ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr);
595   } else {
596     vNatural = v;
597   }
598   /* Get the location of the variable in the exodus file */
599   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
600   ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr);
601   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname);
602   /* Write local chunk of the result in the exodus file
603      exodus stores each component of a vector-valued field as a separate variable.
604      We assume that they are stored sequentially
605      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
606      but once the vector has been reordered to natural size, we cannot use the label informations
607      to figure out what to save where. */
608   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
609   ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr);
610   ierr = ex_get_ids(exoid, EX_ELEM_BLOCK, csID);CHKERRQ(ierr);
611   for (set = 0; set < numCS; ++set) {
612     ex_block block;
613 
614     block.id   = csID[set];
615     block.type = EX_ELEM_BLOCK;
616     ierr = ex_get_block_param(exoid, &block);CHKERRQ(ierr);
617     csSize[set] = block.num_entry;
618   }
619   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
620   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
621   if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);}
622   for (set = 0; set < numCS; set++) {
623     PetscInt csLocalSize, c;
624 
625     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
626        local slice of zonal values:         xs/bs,xm/bs-1
627        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
628     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
629     if (bs == 1) {
630       ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr);
631       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);
632       ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr);
633     } else {
634       for (c = 0; c < bs; ++c) {
635         ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
636         ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
637         ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr);
638         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)]);
639         ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr);
640         ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
641       }
642     }
643     csxs += csSize[set];
644   }
645   ierr = PetscFree2(csID, csSize);CHKERRQ(ierr);
646   if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);}
647   if (useNatural) {ierr = DMRestoreGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
648 #else
649   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
650 #endif
651   PetscFunctionReturn(0);
652 }
653 
654 /*
655   VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file
656 
657   Collective on v
658 
659   Input Parameters:
660 + v  - The vector to be written
661 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
662 - step - the time step to read at (exodus steps are numbered starting from 1)
663 
664   Notes:
665   The exodus result zonal variable index is obtained by comparing the Vec name and the
666   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
667   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
668   amongst all zonal variables.
669 
670   Level: beginner
671 
672 .keywords: mesh,ExodusII
673 .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Zonal_Internal()
674 */
675 PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step)
676 {
677   MPI_Comm          comm;
678   PetscMPIInt       size;
679 #if defined(PETSC_HAVE_EXODUSII)
680   DM                dm;
681   Vec               vNatural, vComp;
682   PetscReal        *varray;
683   const char       *vecname;
684   PetscInt          xs, xe, bs;
685   PetscBool         useNatural;
686   int               offset;
687   IS                compIS;
688   PetscInt         *csSize, *csID;
689   PetscInt          numCS, set, csxs = 0;
690   PetscErrorCode    ierr;
691 #endif
692 
693   PetscFunctionBegin;
694   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
695   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
696 #if defined(PETSC_HAVE_EXODUSII)
697   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
698   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
699   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
700   if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
701   else            {vNatural = v;}
702   /* Get the location of the variable in the exodus file */
703   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
704   ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr);
705   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname);
706   /* Read local chunk of the result in the exodus file
707      exodus stores each component of a vector-valued field as a separate variable.
708      We assume that they are stored sequentially
709      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
710      but once the vector has been reordered to natural size, we cannot use the label informations
711      to figure out what to save where. */
712   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
713   ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr);
714   ierr = ex_get_ids(exoid, EX_ELEM_BLOCK, csID);CHKERRQ(ierr);
715   for (set = 0; set < numCS; ++set) {
716     ex_block block;
717 
718     block.id   = csID[set];
719     block.type = EX_ELEM_BLOCK;
720     ierr = ex_get_block_param(exoid, &block);CHKERRQ(ierr);
721     csSize[set] = block.num_entry;
722   }
723   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
724   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
725   if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);}
726   for (set = 0; set < numCS; ++set) {
727     PetscInt csLocalSize, c;
728 
729     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
730        local slice of zonal values:         xs/bs,xm/bs-1
731        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
732     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
733     if (bs == 1) {
734       ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr);
735       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);
736       ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr);
737     } else {
738       for (c = 0; c < bs; ++c) {
739         ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
740         ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
741         ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr);
742         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)]);
743         ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr);
744         ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
745       }
746     }
747     csxs += csSize[set];
748   }
749   ierr = PetscFree2(csID, csSize);CHKERRQ(ierr);
750   if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);}
751   if (useNatural) {
752     ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr);
753     ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr);
754     ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);
755   }
756 #else
757   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
758 #endif
759   PetscFunctionReturn(0);
760 }
761 
762 /*@C
763   DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file.
764 
765   Collective on comm
766 
767   Input Parameters:
768 + comm  - The MPI communicator
769 . filename - The name of the ExodusII file
770 - interpolate - Create faces and edges in the mesh
771 
772   Output Parameter:
773 . dm  - The DM object representing the mesh
774 
775   Level: beginner
776 
777 .keywords: mesh,ExodusII
778 .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus()
779 @*/
780 PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
781 {
782   PetscMPIInt    rank;
783   PetscErrorCode ierr;
784 #if defined(PETSC_HAVE_EXODUSII)
785   int   CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1;
786   float version;
787 #endif
788 
789   PetscFunctionBegin;
790   PetscValidCharPointer(filename, 2);
791   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
792 #if defined(PETSC_HAVE_EXODUSII)
793   if (!rank) {
794     exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version);
795     if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename);
796   }
797   ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr);
798   if (!rank) {ierr = ex_close(exoid);CHKERRQ(ierr);}
799 #else
800   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
801 #endif
802   PetscFunctionReturn(0);
803 }
804 
805 /*@
806   DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID.
807 
808   Collective on comm
809 
810   Input Parameters:
811 + comm  - The MPI communicator
812 . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
813 - interpolate - Create faces and edges in the mesh
814 
815   Output Parameter:
816 . dm  - The DM object representing the mesh
817 
818   Level: beginner
819 
820 .keywords: mesh,ExodusII
821 .seealso: DMPLEX, DMCreate()
822 @*/
823 PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm)
824 {
825 #if defined(PETSC_HAVE_EXODUSII)
826   PetscMPIInt    num_proc, rank;
827   PetscSection   coordSection;
828   Vec            coordinates;
829   PetscScalar    *coords;
830   PetscInt       coordSize, v;
831   PetscErrorCode ierr;
832   /* Read from ex_get_init() */
833   char title[PETSC_MAX_PATH_LEN+1];
834   int  dim    = 0, numVertices = 0, numCells = 0;
835   int  num_cs = 0, num_vs = 0, num_fs = 0;
836 #endif
837 
838   PetscFunctionBegin;
839 #if defined(PETSC_HAVE_EXODUSII)
840   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
841   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
842   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
843   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
844   /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */
845   if (!rank) {
846     ierr = PetscMemzero(title,(PETSC_MAX_PATH_LEN+1)*sizeof(char));CHKERRQ(ierr);
847     ierr = ex_get_init(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs);
848     if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ExodusII ex_get_init() failed with error code %D\n",ierr);
849     if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n");
850   }
851   ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
852   ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
853   ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr);
854   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
855   ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
856 
857   /* Read cell sets information */
858   if (!rank) {
859     PetscInt *cone;
860     int      c, cs, c_loc, v, v_loc;
861     /* Read from ex_get_elem_blk_ids() */
862     int *cs_id;
863     /* Read from ex_get_elem_block() */
864     char buffer[PETSC_MAX_PATH_LEN+1];
865     int  num_cell_in_set, num_vertex_per_cell, num_attr;
866     /* Read from ex_get_elem_conn() */
867     int *cs_connect;
868 
869     /* Get cell sets IDs */
870     ierr = PetscMalloc1(num_cs, &cs_id);CHKERRQ(ierr);
871     ierr = ex_get_ids (exoid, EX_ELEM_BLOCK, cs_id);CHKERRQ(ierr);
872     /* Read the cell set connectivity table and build mesh topology
873        EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
874     /* First set sizes */
875     for (cs = 0, c = 0; cs < num_cs; cs++) {
876       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);
877       for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
878         ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr);
879       }
880     }
881     ierr = DMSetUp(*dm);CHKERRQ(ierr);
882     for (cs = 0, c = 0; cs < num_cs; cs++) {
883       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);
884       ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr);
885       ierr = ex_get_conn (exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL);CHKERRQ(ierr);
886       /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
887       for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
888         for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) {
889           cone[v_loc] = cs_connect[v]+numCells-1;
890         }
891         if (dim == 3) {
892           /* Tetrahedra are inverted */
893           if (num_vertex_per_cell == 4) {
894             PetscInt tmp = cone[0];
895             cone[0] = cone[1];
896             cone[1] = tmp;
897           }
898           /* Hexahedra are inverted */
899           if (num_vertex_per_cell == 8) {
900             PetscInt tmp = cone[1];
901             cone[1] = cone[3];
902             cone[3] = tmp;
903           }
904         }
905         ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
906         ierr = DMSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr);
907       }
908       ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr);
909     }
910     ierr = PetscFree(cs_id);CHKERRQ(ierr);
911   }
912   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
913   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
914   if (interpolate) {
915     DM idm;
916 
917     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
918     ierr = DMDestroy(dm);CHKERRQ(ierr);
919     *dm  = idm;
920   }
921 
922   /* Create vertex set label */
923   if (!rank && (num_vs > 0)) {
924     int vs, v;
925     /* Read from ex_get_node_set_ids() */
926     int *vs_id;
927     /* Read from ex_get_node_set_param() */
928     int num_vertex_in_set;
929     /* Read from ex_get_node_set() */
930     int *vs_vertex_list;
931 
932     /* Get vertex set ids */
933     ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr);
934     ierr = ex_get_ids(exoid, EX_NODE_SET, vs_id);CHKERRQ(ierr);
935     for (vs = 0; vs < num_vs; ++vs) {
936       ierr = ex_get_set_param(exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL);
937       ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr);
938       ierr = ex_get_set(exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL);CHKERRQ(ierr);
939       for (v = 0; v < num_vertex_in_set; ++v) {
940         ierr = DMSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr);
941       }
942       ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr);
943     }
944     ierr = PetscFree(vs_id);CHKERRQ(ierr);
945   }
946   /* Read coordinates */
947   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
948   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
949   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
950   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
951   for (v = numCells; v < numCells+numVertices; ++v) {
952     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
953     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
954   }
955   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
956   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
957   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
958   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
959   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
960   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
961   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
962   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
963   if (!rank) {
964     PetscReal *x, *y, *z;
965 
966     ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr);
967     ierr = ex_get_coord(exoid, x, y, z);CHKERRQ(ierr);
968     if (dim > 0) {
969       for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v];
970     }
971     if (dim > 1) {
972       for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v];
973     }
974     if (dim > 2) {
975       for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v];
976     }
977     ierr = PetscFree3(x,y,z);CHKERRQ(ierr);
978   }
979   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
980   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
981   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
982 
983   /* Create side set label */
984   if (!rank && interpolate && (num_fs > 0)) {
985     int fs, f, voff;
986     /* Read from ex_get_side_set_ids() */
987     int *fs_id;
988     /* Read from ex_get_side_set_param() */
989     int num_side_in_set;
990     /* Read from ex_get_side_set_node_list() */
991     int *fs_vertex_count_list, *fs_vertex_list;
992     /* Read side set labels */
993     char fs_name[MAX_STR_LENGTH+1];
994 
995     /* Get side set ids */
996     ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr);
997     ierr = ex_get_ids(exoid, EX_SIDE_SET, fs_id);CHKERRQ(ierr);
998     for (fs = 0; fs < num_fs; ++fs) {
999       ierr = ex_get_set_param(exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL);CHKERRQ(ierr);
1000       ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr);
1001       ierr = ex_get_side_set_node_list(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);CHKERRQ(ierr);
1002       /* Get the specific name associated with this side set ID. */
1003       int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name);
1004       for (f = 0, voff = 0; f < num_side_in_set; ++f) {
1005         const PetscInt *faces   = NULL;
1006         PetscInt       faceSize = fs_vertex_count_list[f], numFaces;
1007         PetscInt       faceVertices[4], v;
1008 
1009         if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize);
1010         for (v = 0; v < faceSize; ++v, ++voff) {
1011           faceVertices[v] = fs_vertex_list[voff]+numCells-1;
1012         }
1013         ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
1014         if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces);
1015         ierr = DMSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr);
1016         /* Only add the label if one has been detected for this side set. */
1017         if (!fs_name_err) {
1018           ierr = DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]);CHKERRQ(ierr);
1019         }
1020         ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
1021       }
1022       ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr);
1023     }
1024     ierr = PetscFree(fs_id);CHKERRQ(ierr);
1025   }
1026 #else
1027   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
1028 #endif
1029   PetscFunctionReturn(0);
1030 }
1031