xref: /petsc/src/dm/impls/plex/plexexodusii.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
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 #include <petsc/private/viewerimpl.h>
10 #include <petsc/private/viewerexodusiiimpl.h>
11 #if defined(PETSC_HAVE_EXODUSII)
12 /*
13   PETSC_VIEWER_EXODUSII_ - Creates an ExodusII PetscViewer shared by all processors in a communicator.
14 
15   Collective
16 
17   Input Parameter:
18 . comm - the MPI communicator to share the ExodusII PetscViewer
19 
20   Level: intermediate
21 
22   Notes:
23     misses Fortran bindings
24 
25   Notes:
26   Unlike almost all other PETSc routines, PETSC_VIEWER_EXODUSII_ does not return
27   an error code.  The GLVIS PetscViewer is usually used in the form
28 $       XXXView(XXX object, PETSC_VIEWER_EXODUSII_(comm));
29 
30 .seealso: `PetscViewerExodusIIOpen()`, `PetscViewerType`, `PetscViewerCreate()`, `PetscViewerDestroy()`
31 */
32 PetscViewer PETSC_VIEWER_EXODUSII_(MPI_Comm comm)
33 {
34   PetscViewer    viewer;
35   PetscErrorCode ierr;
36 
37   PetscFunctionBegin;
38   ierr = PetscViewerExodusIIOpen(comm, "mesh.exo", FILE_MODE_WRITE, &viewer);
39   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_EXODUSII_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
40   ierr = PetscObjectRegisterDestroy((PetscObject) viewer);
41   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_EXODUSII_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
42   PetscFunctionReturn(viewer);
43 }
44 
45 static PetscErrorCode PetscViewerView_ExodusII(PetscViewer v, PetscViewer viewer)
46 {
47   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) v->data;
48 
49   PetscFunctionBegin;
50   if (exo->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename:    %s\n", exo->filename));
51   if (exo->exoid)    PetscCall(PetscViewerASCIIPrintf(viewer, "exoid:       %d\n", exo->exoid));
52   if (exo->btype)    PetscCall(PetscViewerASCIIPrintf(viewer, "IO Mode:     %d\n", exo->btype));
53   if (exo->order)    PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh order:  %d\n", exo->order));
54   PetscFunctionReturn(0);
55 }
56 
57 static PetscErrorCode PetscViewerSetFromOptions_ExodusII(PetscOptionItems *PetscOptionsObject, PetscViewer v)
58 {
59   PetscFunctionBegin;
60   PetscOptionsHeadBegin(PetscOptionsObject, "ExodusII PetscViewer Options");
61   PetscOptionsHeadEnd();
62   PetscFunctionReturn(0);
63 }
64 
65 static PetscErrorCode PetscViewerSetUp_ExodusII(PetscViewer viewer)
66 {
67   PetscFunctionBegin;
68   PetscFunctionReturn(0);
69 }
70 
71 static PetscErrorCode PetscViewerDestroy_ExodusII(PetscViewer viewer)
72 {
73   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
74 
75   PetscFunctionBegin;
76   if (exo->exoid >= 0) {PetscStackCallStandard(ex_close,exo->exoid);}
77   PetscCall(PetscFree(exo->filename));
78   PetscCall(PetscFree(exo));
79   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL));
80   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL));
81   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL));
82   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetMode_C",NULL));
83   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGetId_C",NULL));
84   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGetOrder_C",NULL));
85   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerSetOrder_C",NULL));
86   PetscFunctionReturn(0);
87 }
88 
89 static PetscErrorCode PetscViewerFileSetName_ExodusII(PetscViewer viewer, const char name[])
90 {
91   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
92   PetscMPIInt           rank;
93   int                   CPU_word_size, IO_word_size, EXO_mode;
94   MPI_Info              mpi_info = MPI_INFO_NULL;
95   float                 EXO_version;
96 
97   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) viewer), &rank));
98   CPU_word_size = sizeof(PetscReal);
99   IO_word_size  = sizeof(PetscReal);
100 
101   PetscFunctionBegin;
102   if (exo->exoid >= 0) {
103     PetscStackCallStandard(ex_close,exo->exoid);
104     exo->exoid = -1;
105   }
106   if (exo->filename) PetscCall(PetscFree(exo->filename));
107   PetscCall(PetscStrallocpy(name, &exo->filename));
108   switch (exo->btype) {
109   case FILE_MODE_READ:
110     EXO_mode = EX_READ;
111     break;
112   case FILE_MODE_APPEND:
113   case FILE_MODE_UPDATE:
114   case FILE_MODE_APPEND_UPDATE:
115     /* Will fail if the file does not already exist */
116     EXO_mode = EX_WRITE;
117     break;
118   case FILE_MODE_WRITE:
119     /*
120       exodus only allows writing geometry upon file creation, so we will let DMView create the file.
121     */
122     PetscFunctionReturn(0);
123     break;
124   default:
125     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
126   }
127   #if defined(PETSC_USE_64BIT_INDICES)
128   EXO_mode += EX_ALL_INT64_API;
129   #endif
130   exo->exoid = ex_open_par(name,EXO_mode,&CPU_word_size,&IO_word_size,&EXO_version,PETSC_COMM_WORLD,mpi_info);
131   PetscCheck(exo->exoid >= 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", name);
132   PetscFunctionReturn(0);
133 }
134 
135 static PetscErrorCode PetscViewerFileGetName_ExodusII(PetscViewer viewer, const char **name)
136 {
137   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
138 
139   PetscFunctionBegin;
140   *name = exo->filename;
141   PetscFunctionReturn(0);
142 }
143 
144 static PetscErrorCode PetscViewerFileSetMode_ExodusII(PetscViewer viewer, PetscFileMode type)
145 {
146   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
147 
148   PetscFunctionBegin;
149   exo->btype = type;
150   PetscFunctionReturn(0);
151 }
152 
153 static PetscErrorCode PetscViewerFileGetMode_ExodusII(PetscViewer viewer, PetscFileMode *type)
154 {
155   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
156 
157   PetscFunctionBegin;
158   *type = exo->btype;
159   PetscFunctionReturn(0);
160 }
161 
162 static PetscErrorCode PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer, int *exoid)
163 {
164   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
165 
166   PetscFunctionBegin;
167   *exoid = exo->exoid;
168   PetscFunctionReturn(0);
169 }
170 
171 static PetscErrorCode PetscViewerExodusIIGetOrder_ExodusII(PetscViewer viewer, PetscInt *order)
172 {
173   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
174 
175   PetscFunctionBegin;
176   *order = exo->order;
177   PetscFunctionReturn(0);
178 }
179 
180 static PetscErrorCode PetscViewerExodusIISetOrder_ExodusII(PetscViewer viewer, PetscInt order)
181 {
182   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
183 
184   PetscFunctionBegin;
185   exo->order = order;
186   PetscFunctionReturn(0);
187 }
188 
189 /*MC
190    PETSCVIEWEREXODUSII - A viewer that writes to an Exodus II file
191 
192 .seealso: `PetscViewerExodusIIOpen()`, `PetscViewerCreate()`, `PETSCVIEWERBINARY`, `PETSCVIEWERHDF5`, `DMView()`,
193           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
194 
195   Level: beginner
196 M*/
197 
198 PETSC_EXTERN PetscErrorCode PetscViewerCreate_ExodusII(PetscViewer v)
199 {
200   PetscViewer_ExodusII *exo;
201 
202   PetscFunctionBegin;
203   PetscCall(PetscNewLog(v,&exo));
204 
205   v->data                = (void*) exo;
206   v->ops->destroy        = PetscViewerDestroy_ExodusII;
207   v->ops->setfromoptions = PetscViewerSetFromOptions_ExodusII;
208   v->ops->setup          = PetscViewerSetUp_ExodusII;
209   v->ops->view           = PetscViewerView_ExodusII;
210   v->ops->flush          = 0;
211   exo->btype             = (PetscFileMode) -1;
212   exo->filename          = 0;
213   exo->exoid             = -1;
214 
215   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_ExodusII));
216   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_ExodusII));
217   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_ExodusII));
218   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_ExodusII));
219   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerGetId_C",PetscViewerExodusIIGetId_ExodusII));
220   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerSetOrder_C",PetscViewerExodusIISetOrder_ExodusII));
221   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerGetOrder_C",PetscViewerExodusIIGetOrder_ExodusII));
222   PetscFunctionReturn(0);
223 }
224 
225 /*
226   EXOGetVarIndex - Locate a result in an exodus file based on its name
227 
228   Collective
229 
230   Input Parameters:
231 + exoid    - the exodus id of a file (obtained from ex_open or ex_create for instance)
232 . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK
233 - name     - the name of the result
234 
235   Output Parameters:
236 . varIndex - the location in the exodus file of the result
237 
238   Notes:
239   The exodus variable index is obtained by comparing name and the
240   names of zonal variables declared in the exodus file. For instance if name is "V"
241   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
242   amongst all variables of type obj_type.
243 
244   Level: beginner
245 
246 .seealso: `EXOGetVarIndex()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadNodal_PlexEXO()`, `VecLoadZonal_PlexEXO()`
247 */
248 PetscErrorCode EXOGetVarIndex_Internal(int exoid, ex_entity_type obj_type, const char name[], int *varIndex)
249 {
250   int            num_vars, i, j;
251   char           ext_name[MAX_STR_LENGTH+1], var_name[MAX_STR_LENGTH+1];
252   const int      num_suffix = 5;
253   char          *suffix[5];
254   PetscBool      flg;
255 
256   PetscFunctionBegin;
257   suffix[0] = (char *) "";
258   suffix[1] = (char *) "_X";
259   suffix[2] = (char *) "_XX";
260   suffix[3] = (char *) "_1";
261   suffix[4] = (char *) "_11";
262 
263   *varIndex = -1;
264   PetscStackCallStandard(ex_get_variable_param,exoid, obj_type, &num_vars);
265   for (i = 0; i < num_vars; ++i) {
266     PetscStackCallStandard(ex_get_variable_name,exoid, obj_type, i+1, var_name);
267     for (j = 0; j < num_suffix; ++j) {
268       PetscCall(PetscStrncpy(ext_name, name, MAX_STR_LENGTH));
269       PetscCall(PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH));
270       PetscCall(PetscStrcasecmp(ext_name, var_name, &flg));
271       if (flg) {
272         *varIndex = i+1;
273       }
274     }
275   }
276   PetscFunctionReturn(0);
277 }
278 
279 /*
280   DMView_PlexExodusII - Write a DM to disk in exodus format
281 
282   Collective on dm
283 
284   Input Parameters:
285 + dm  - The dm to be written
286 . viewer - an exodusII viewer
287 
288   Notes:
289   Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels)
290   consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted.
291 
292   If the dm has been distributed, only the part of the DM on MPI rank 0 (including "ghost" cells and vertices)
293   will be written.
294 
295   DMPlex only represents geometry while most post-processing software expect that a mesh also provides information
296   on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2.
297   The order of the mesh shall be set using PetscViewerExodusIISetOrder
298   It should be extended to use PetscFE objects.
299 
300   This function will only handle TRI, TET, QUAD, and HEX cells.
301   Level: beginner
302 
303 .seealso:
304 */
305 PetscErrorCode DMView_PlexExodusII(DM dm, PetscViewer viewer)
306 {
307   enum ElemType {TRI, QUAD, TET, HEX};
308   MPI_Comm        comm;
309   PetscInt        degree; /* the order of the mesh */
310   /* Connectivity Variables */
311   PetscInt        cellsNotInConnectivity;
312   /* Cell Sets */
313   DMLabel         csLabel;
314   IS              csIS;
315   const PetscInt *csIdx;
316   PetscInt        num_cs, cs;
317   enum ElemType  *type;
318   PetscBool       hasLabel;
319   /* Coordinate Variables */
320   DM              cdm;
321   PetscSection    coordSection;
322   Vec             coord;
323   PetscInt      **nodes;
324   PetscInt        depth, d, dim, skipCells = 0;
325   PetscInt        pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes;
326   PetscInt        num_vs, num_fs;
327   PetscMPIInt     rank, size;
328   const char     *dmName;
329   PetscInt        nodesTriP1[4]  = {3,0,0,0};
330   PetscInt        nodesTriP2[4]  = {3,3,0,0};
331   PetscInt        nodesQuadP1[4] = {4,0,0,0};
332   PetscInt        nodesQuadP2[4] = {4,4,0,1};
333   PetscInt        nodesTetP1[4]  = {4,0,0,0};
334   PetscInt        nodesTetP2[4]  = {4,6,0,0};
335   PetscInt        nodesHexP1[4]  = {8,0,0,0};
336   PetscInt        nodesHexP2[4]  = {8,12,6,1};
337   int             CPU_word_size, IO_word_size, EXO_mode;
338   float           EXO_version;
339 
340   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
341 
342   PetscFunctionBegin;
343   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
344   PetscCallMPI(MPI_Comm_rank(comm, &rank));
345   PetscCallMPI(MPI_Comm_size(comm, &size));
346 
347   /*
348     Creating coordSection is a collective operation so we do it somewhat out of sequence
349   */
350   PetscCall(PetscSectionCreate(comm, &coordSection));
351   PetscCall(DMGetCoordinatesLocalSetUp(dm));
352   if (rank == 0) {
353     switch (exo->btype) {
354     case FILE_MODE_READ:
355     case FILE_MODE_APPEND:
356     case FILE_MODE_UPDATE:
357     case FILE_MODE_APPEND_UPDATE:
358       /* exodusII does not allow writing geometry to an existing file */
359       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "cannot add geometry to existing file %s", exo->filename);
360     case FILE_MODE_WRITE:
361       /* Create an empty file if one already exists*/
362       EXO_mode = EX_CLOBBER;
363 #if defined(PETSC_USE_64BIT_INDICES)
364       EXO_mode += EX_ALL_INT64_API;
365 #endif
366         CPU_word_size = sizeof(PetscReal);
367         IO_word_size  = sizeof(PetscReal);
368         exo->exoid = ex_create(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size);
369         PetscCheck(exo->exoid >= 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", exo->filename);
370 
371       break;
372     default:
373       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
374     }
375 
376     /* --- Get DM info --- */
377     PetscCall(PetscObjectGetName((PetscObject) dm, &dmName));
378     PetscCall(DMPlexGetDepth(dm, &depth));
379     PetscCall(DMGetDimension(dm, &dim));
380     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
381     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
382     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
383     PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
384     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
385     numCells    = cEnd - cStart;
386     numEdges    = eEnd - eStart;
387     numVertices = vEnd - vStart;
388     if (depth == 3) {numFaces = fEnd - fStart;}
389     else            {numFaces = 0;}
390     PetscCall(DMGetLabelSize(dm, "Cell Sets", &num_cs));
391     PetscCall(DMGetLabelSize(dm, "Vertex Sets", &num_vs));
392     PetscCall(DMGetLabelSize(dm, "Face Sets", &num_fs));
393     PetscCall(DMGetCoordinatesLocal(dm, &coord));
394     PetscCall(DMGetCoordinateDM(dm, &cdm));
395     if (num_cs > 0) {
396       PetscCall(DMGetLabel(dm, "Cell Sets", &csLabel));
397       PetscCall(DMLabelGetValueIS(csLabel, &csIS));
398       PetscCall(ISGetIndices(csIS, &csIdx));
399     }
400     PetscCall(PetscMalloc1(num_cs, &nodes));
401     /* Set element type for each block and compute total number of nodes */
402     PetscCall(PetscMalloc1(num_cs, &type));
403     numNodes = numVertices;
404 
405     PetscCall(PetscViewerExodusIIGetOrder(viewer, &degree));
406     if (degree == 2) {numNodes += numEdges;}
407     cellsNotInConnectivity = numCells;
408     for (cs = 0; cs < num_cs; ++cs) {
409       IS              stratumIS;
410       const PetscInt *cells;
411       PetscScalar    *xyz = NULL;
412       PetscInt        csSize, closureSize;
413 
414       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
415       PetscCall(ISGetIndices(stratumIS, &cells));
416       PetscCall(ISGetSize(stratumIS, &csSize));
417       PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz));
418       switch (dim) {
419         case 2:
420           if      (closureSize == 3*dim) {type[cs] = TRI;}
421           else if (closureSize == 4*dim) {type[cs] = QUAD;}
422           else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %" PetscInt_FMT " in dimension %" PetscInt_FMT " has no ExodusII type", closureSize/dim, dim);
423           break;
424         case 3:
425           if      (closureSize == 4*dim) {type[cs] = TET;}
426           else if (closureSize == 8*dim) {type[cs] = HEX;}
427           else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %" PetscInt_FMT " in dimension %" PetscInt_FMT " has no ExodusII type", closureSize/dim, dim);
428           break;
429         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not handled by ExodusII viewer", dim);
430       }
431       if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;}
432       if ((degree == 2) && (type[cs] == HEX))  {numNodes += csSize; numNodes += numFaces;}
433       PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz));
434       /* Set nodes and Element type */
435       if (type[cs] == TRI) {
436         if      (degree == 1) nodes[cs] = nodesTriP1;
437         else if (degree == 2) nodes[cs] = nodesTriP2;
438       } else if (type[cs] == QUAD) {
439         if      (degree == 1) nodes[cs] = nodesQuadP1;
440         else if (degree == 2) nodes[cs] = nodesQuadP2;
441       } else if (type[cs] == TET) {
442         if      (degree == 1) nodes[cs] = nodesTetP1;
443         else if (degree == 2) nodes[cs] = nodesTetP2;
444       } else if (type[cs] == HEX) {
445         if      (degree == 1) nodes[cs] = nodesHexP1;
446         else if (degree == 2) nodes[cs] = nodesHexP2;
447       }
448       /* Compute the number of cells not in the connectivity table */
449       cellsNotInConnectivity -= nodes[cs][3]*csSize;
450 
451       PetscCall(ISRestoreIndices(stratumIS, &cells));
452       PetscCall(ISDestroy(&stratumIS));
453     }
454     if (num_cs > 0) {PetscStackCallStandard(ex_put_init,exo->exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs);}
455     /* --- Connectivity --- */
456     for (cs = 0; cs < num_cs; ++cs) {
457       IS              stratumIS;
458       const PetscInt *cells;
459       PetscInt       *connect, off = 0;
460       PetscInt        edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0;
461       PetscInt        csSize, c, connectSize, closureSize;
462       char           *elem_type = NULL;
463       char            elem_type_tri3[]  = "TRI3",  elem_type_quad4[] = "QUAD4";
464       char            elem_type_tri6[]  = "TRI6",  elem_type_quad9[] = "QUAD9";
465       char            elem_type_tet4[]  = "TET4",  elem_type_hex8[]  = "HEX8";
466       char            elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27";
467 
468       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
469       PetscCall(ISGetIndices(stratumIS, &cells));
470       PetscCall(ISGetSize(stratumIS, &csSize));
471       /* Set Element type */
472       if (type[cs] == TRI) {
473         if      (degree == 1) elem_type = elem_type_tri3;
474         else if (degree == 2) elem_type = elem_type_tri6;
475       } else if (type[cs] == QUAD) {
476         if      (degree == 1) elem_type = elem_type_quad4;
477         else if (degree == 2) elem_type = elem_type_quad9;
478       } else if (type[cs] == TET) {
479         if      (degree == 1) elem_type = elem_type_tet4;
480         else if (degree == 2) elem_type = elem_type_tet10;
481       } else if (type[cs] == HEX) {
482         if      (degree == 1) elem_type = elem_type_hex8;
483         else if (degree == 2) elem_type = elem_type_hex27;
484       }
485       connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3];
486       PetscCall(PetscMalloc1(PetscMax(27,connectSize)*csSize, &connect));
487       PetscStackCallStandard(ex_put_block,exo->exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1);
488       /* Find number of vertices, edges, and faces in the closure */
489       verticesInClosure = nodes[cs][0];
490       if (depth > 1) {
491         if (dim == 2) {
492           PetscCall(DMPlexGetConeSize(dm, cells[0], &edgesInClosure));
493         } else if (dim == 3) {
494           PetscInt *closure = NULL;
495 
496           PetscCall(DMPlexGetConeSize(dm, cells[0], &facesInClosure));
497           PetscCall(DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure));
498           edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure;
499           PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure));
500         }
501       }
502       /* Get connectivity for each cell */
503       for (c = 0; c < csSize; ++c) {
504         PetscInt *closure = NULL;
505         PetscInt  temp, i;
506 
507         PetscCall(DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure));
508         for (i = 0; i < connectSize; ++i) {
509           if (i < nodes[cs][0]) {/* Vertices */
510             connect[i+off] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1;
511             connect[i+off] -= cellsNotInConnectivity;
512           } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */
513             connect[i+off] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1;
514             if (nodes[cs][2] == 0) connect[i+off] -= numFaces;
515             connect[i+off] -= cellsNotInConnectivity;
516           } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]) { /* Cells */
517             connect[i+off] = closure[0] + 1;
518             connect[i+off] -= skipCells;
519           } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]) { /* Faces */
520             connect[i+off] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1;
521             connect[i+off] -= cellsNotInConnectivity;
522           } else {
523             connect[i+off] = -1;
524           }
525         }
526         /* Tetrahedra are inverted */
527         if (type[cs] == TET) {
528           temp = connect[0+off]; connect[0+off] = connect[1+off]; connect[1+off] = temp;
529           if (degree == 2) {
530             temp = connect[5+off]; connect[5+off] = connect[6+off]; connect[6+off] = temp;
531             temp = connect[7+off]; connect[7+off] = connect[8+off]; connect[8+off] = temp;
532           }
533         }
534         /* Hexahedra are inverted */
535         if (type[cs] == HEX) {
536           temp = connect[1+off]; connect[1+off] = connect[3+off]; connect[3+off] = temp;
537           if (degree == 2) {
538             temp = connect[8+off];  connect[8+off]  = connect[11+off]; connect[11+off] = temp;
539             temp = connect[9+off];  connect[9+off]  = connect[10+off]; connect[10+off] = temp;
540             temp = connect[16+off]; connect[16+off] = connect[17+off]; connect[17+off] = temp;
541             temp = connect[18+off]; connect[18+off] = connect[19+off]; connect[19+off] = temp;
542 
543             temp = connect[12+off]; connect[12+off] = connect[16+off]; connect[16+off] = temp;
544             temp = connect[13+off]; connect[13+off] = connect[17+off]; connect[17+off] = temp;
545             temp = connect[14+off]; connect[14+off] = connect[18+off]; connect[18+off] = temp;
546             temp = connect[15+off]; connect[15+off] = connect[19+off]; connect[19+off] = temp;
547 
548             temp = connect[23+off]; connect[23+off] = connect[26+off]; connect[26+off] = temp;
549             temp = connect[24+off]; connect[24+off] = connect[25+off]; connect[25+off] = temp;
550             temp = connect[25+off]; connect[25+off] = connect[26+off]; connect[26+off] = temp;
551           }
552         }
553         off += connectSize;
554         PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure));
555       }
556       PetscStackCallStandard(ex_put_conn,exo->exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0);
557       skipCells += (nodes[cs][3] == 0)*csSize;
558       PetscCall(PetscFree(connect));
559       PetscCall(ISRestoreIndices(stratumIS, &cells));
560       PetscCall(ISDestroy(&stratumIS));
561     }
562     PetscCall(PetscFree(type));
563     /* --- Coordinates --- */
564     PetscCall(PetscSectionSetChart(coordSection, pStart, pEnd));
565     if (num_cs) {
566       for (d = 0; d < depth; ++d) {
567         PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
568         for (p = pStart; p < pEnd; ++p) {
569           PetscCall(PetscSectionSetDof(coordSection, p, nodes[0][d] > 0));
570         }
571       }
572     }
573     for (cs = 0; cs < num_cs; ++cs) {
574       IS              stratumIS;
575       const PetscInt *cells;
576       PetscInt        csSize, c;
577 
578       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
579       PetscCall(ISGetIndices(stratumIS, &cells));
580       PetscCall(ISGetSize(stratumIS, &csSize));
581       for (c = 0; c < csSize; ++c) {
582         PetscCall(PetscSectionSetDof(coordSection, cells[c], nodes[cs][3] > 0));
583       }
584       PetscCall(ISRestoreIndices(stratumIS, &cells));
585       PetscCall(ISDestroy(&stratumIS));
586     }
587     if (num_cs > 0) {
588       PetscCall(ISRestoreIndices(csIS, &csIdx));
589       PetscCall(ISDestroy(&csIS));
590     }
591     PetscCall(PetscFree(nodes));
592     PetscCall(PetscSectionSetUp(coordSection));
593     if (numNodes > 0) {
594       const char  *coordNames[3] = {"x", "y", "z"};
595       PetscScalar *closure, *cval;
596       PetscReal   *coords;
597       PetscInt     hasDof, n = 0;
598 
599       /* There can't be more than 24 values in the closure of a point for the coord coordSection */
600       PetscCall(PetscCalloc3(numNodes*3, &coords, dim, &cval, 24, &closure));
601       PetscCall(DMGetCoordinatesLocalNoncollective(dm, &coord));
602       PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
603       for (p = pStart; p < pEnd; ++p) {
604         PetscCall(PetscSectionGetDof(coordSection, p, &hasDof));
605         if (hasDof) {
606           PetscInt closureSize = 24, j;
607 
608           PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure));
609           for (d = 0; d < dim; ++d) {
610             cval[d] = 0.0;
611             for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d];
612             coords[d*numNodes+n] = PetscRealPart(cval[d]) * dim / closureSize;
613           }
614           ++n;
615         }
616       }
617       PetscStackCallStandard(ex_put_coord,exo->exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes]);
618       PetscCall(PetscFree3(coords, cval, closure));
619       PetscStackCallStandard(ex_put_coord_names,exo->exoid, (char **) coordNames);
620     }
621 
622     /* --- Node Sets/Vertex Sets --- */
623     DMHasLabel(dm, "Vertex Sets", &hasLabel);
624     if (hasLabel) {
625       PetscInt        i, vs, vsSize;
626       const PetscInt *vsIdx, *vertices;
627       PetscInt       *nodeList;
628       IS              vsIS, stratumIS;
629       DMLabel         vsLabel;
630       PetscCall(DMGetLabel(dm, "Vertex Sets", &vsLabel));
631       PetscCall(DMLabelGetValueIS(vsLabel, &vsIS));
632       PetscCall(ISGetIndices(vsIS, &vsIdx));
633       for (vs=0; vs<num_vs; ++vs) {
634         PetscCall(DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS));
635         PetscCall(ISGetIndices(stratumIS, &vertices));
636         PetscCall(ISGetSize(stratumIS, &vsSize));
637         PetscCall(PetscMalloc1(vsSize, &nodeList));
638         for (i=0; i<vsSize; ++i) {
639           nodeList[i] = vertices[i] - skipCells + 1;
640         }
641         PetscStackCallStandard(ex_put_set_param,exo->exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0);
642         PetscStackCallStandard(ex_put_set,exo->exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL);
643         PetscCall(ISRestoreIndices(stratumIS, &vertices));
644         PetscCall(ISDestroy(&stratumIS));
645         PetscCall(PetscFree(nodeList));
646       }
647       PetscCall(ISRestoreIndices(vsIS, &vsIdx));
648       PetscCall(ISDestroy(&vsIS));
649     }
650     /* --- Side Sets/Face Sets --- */
651     PetscCall(DMHasLabel(dm, "Face Sets", &hasLabel));
652     if (hasLabel) {
653       PetscInt        i, j, fs, fsSize;
654       const PetscInt *fsIdx, *faces;
655       IS              fsIS, stratumIS;
656       DMLabel         fsLabel;
657       PetscInt        numPoints, *points;
658       PetscInt        elem_list_size = 0;
659       PetscInt       *elem_list, *elem_ind, *side_list;
660 
661       PetscCall(DMGetLabel(dm, "Face Sets", &fsLabel));
662       /* Compute size of Node List and Element List */
663       PetscCall(DMLabelGetValueIS(fsLabel, &fsIS));
664       PetscCall(ISGetIndices(fsIS, &fsIdx));
665       for (fs=0; fs<num_fs; ++fs) {
666         PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS));
667         PetscCall(ISGetSize(stratumIS, &fsSize));
668         elem_list_size += fsSize;
669         PetscCall(ISDestroy(&stratumIS));
670       }
671       PetscCall(PetscMalloc3(num_fs, &elem_ind,elem_list_size, &elem_list,elem_list_size, &side_list));
672       elem_ind[0] = 0;
673       for (fs=0; fs<num_fs; ++fs) {
674         PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS));
675         PetscCall(ISGetIndices(stratumIS, &faces));
676         PetscCall(ISGetSize(stratumIS, &fsSize));
677         /* Set Parameters */
678         PetscStackCallStandard(ex_put_set_param,exo->exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0);
679         /* Indices */
680         if (fs<num_fs-1) {
681           elem_ind[fs+1] = elem_ind[fs] + fsSize;
682         }
683 
684         for (i=0; i<fsSize; ++i) {
685           /* Element List */
686           points = NULL;
687           PetscCall(DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points));
688           elem_list[elem_ind[fs] + i] = points[2] +1;
689           PetscCall(DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points));
690 
691           /* Side List */
692           points = NULL;
693           PetscCall(DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points));
694           for (j=1; j<numPoints; ++j) {
695             if (points[j*2]==faces[i]) {break;}
696           }
697           /* Convert HEX sides */
698           if (numPoints == 27) {
699             if      (j == 1) {j = 5;}
700             else if (j == 2) {j = 6;}
701             else if (j == 3) {j = 1;}
702             else if (j == 4) {j = 3;}
703             else if (j == 5) {j = 2;}
704             else if (j == 6) {j = 4;}
705           }
706           /* Convert TET sides */
707           if (numPoints == 15) {
708             --j;
709             if (j == 0) {j = 4;}
710           }
711           side_list[elem_ind[fs] + i] = j;
712           PetscCall(DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points));
713 
714         }
715         PetscCall(ISRestoreIndices(stratumIS, &faces));
716         PetscCall(ISDestroy(&stratumIS));
717       }
718       PetscCall(ISRestoreIndices(fsIS, &fsIdx));
719       PetscCall(ISDestroy(&fsIS));
720 
721       /* Put side sets */
722       for (fs=0; fs<num_fs; ++fs) {
723         PetscStackCallStandard(ex_put_set,exo->exoid, EX_SIDE_SET, fsIdx[fs], &elem_list[elem_ind[fs]], &side_list[elem_ind[fs]]);
724       }
725       PetscCall(PetscFree3(elem_ind,elem_list,side_list));
726     }
727     /*
728       close the exodus file
729     */
730     ex_close(exo->exoid);
731     exo->exoid = -1;
732   }
733   PetscCall(PetscSectionDestroy(&coordSection));
734 
735   /*
736     reopen the file in parallel
737   */
738   EXO_mode = EX_WRITE;
739 #if defined(PETSC_USE_64BIT_INDICES)
740   EXO_mode += EX_ALL_INT64_API;
741 #endif
742   CPU_word_size = sizeof(PetscReal);
743   IO_word_size  = sizeof(PetscReal);
744   exo->exoid = ex_open_par(exo->filename,EXO_mode,&CPU_word_size,&IO_word_size,&EXO_version,comm,MPI_INFO_NULL);
745   PetscCheck(exo->exoid >= 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", exo->filename);
746   PetscFunctionReturn(0);
747 }
748 
749 /*
750   VecView_PlexExodusII_Internal - Write a Vec corresponding in an exodus file
751 
752   Collective on v
753 
754   Input Parameters:
755 + v  - The vector to be written
756 - viewer - The PetscViewerExodusII viewer associate to an exodus file
757 
758   Notes:
759   The exodus result variable index is obtained by comparing the Vec name and the
760   names of variables declared in the exodus file. For instance for a Vec named "V"
761   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
762   amongst all variables.
763   In the event where a nodal and zonal variable both match, the function will return an error instead of
764   possibly corrupting the file
765 
766   Level: beginner
767 
768 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()`
769 @*/
770 PetscErrorCode VecView_PlexExodusII_Internal(Vec v, PetscViewer viewer)
771 {
772   DM                 dm;
773   MPI_Comm           comm;
774   PetscMPIInt        rank;
775 
776   int                exoid,offsetN = 0, offsetZ = 0;
777   const char        *vecname;
778   PetscInt           step;
779 
780   PetscFunctionBegin;
781   PetscCall(PetscObjectGetComm((PetscObject) v, &comm));
782   PetscCallMPI(MPI_Comm_rank(comm, &rank));
783   PetscCall(PetscViewerExodusIIGetId(viewer,&exoid));
784   PetscCall(VecGetDM(v, &dm));
785   PetscCall(PetscObjectGetName((PetscObject) v, &vecname));
786 
787   PetscCall(DMGetOutputSequenceNumber(dm,&step,NULL));
788   PetscCall(EXOGetVarIndex_Internal(exoid,EX_NODAL,vecname,&offsetN));
789   PetscCall(EXOGetVarIndex_Internal(exoid,EX_ELEM_BLOCK,vecname,&offsetZ));
790   PetscCheck(offsetN > 0 || offsetZ > 0,comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname);
791   if (offsetN > 0) {
792     PetscCall(VecViewPlex_ExodusII_Nodal_Internal(v,exoid,(int) step+1,offsetN));
793   } else if (offsetZ > 0) {
794     PetscCall(VecViewPlex_ExodusII_Zonal_Internal(v,exoid,(int) step+1,offsetZ));
795   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname);
796   PetscFunctionReturn(0);
797 }
798 
799 /*
800   VecLoad_PlexExodusII_Internal - Write a Vec corresponding in an exodus file
801 
802   Collective on v
803 
804   Input Parameters:
805 + v  - The vector to be written
806 - viewer - The PetscViewerExodusII viewer associate to an exodus file
807 
808   Notes:
809   The exodus result variable index is obtained by comparing the Vec name and the
810   names of variables declared in the exodus file. For instance for a Vec named "V"
811   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
812   amongst all variables.
813   In the event where a nodal and zonal variable both match, the function will return an error instead of
814   possibly corrupting the file
815 
816   Level: beginner
817 
818 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()`
819 @*/
820 PetscErrorCode VecLoad_PlexExodusII_Internal(Vec v, PetscViewer viewer)
821 {
822   DM                 dm;
823   MPI_Comm           comm;
824   PetscMPIInt        rank;
825 
826   int                exoid,offsetN = 0, offsetZ = 0;
827   const char        *vecname;
828   PetscInt           step;
829 
830   PetscFunctionBegin;
831   PetscCall(PetscObjectGetComm((PetscObject) v, &comm));
832   PetscCallMPI(MPI_Comm_rank(comm, &rank));
833   PetscCall(PetscViewerExodusIIGetId(viewer,&exoid));
834   PetscCall(VecGetDM(v, &dm));
835   PetscCall(PetscObjectGetName((PetscObject) v, &vecname));
836 
837   PetscCall(DMGetOutputSequenceNumber(dm,&step,NULL));
838   PetscCall(EXOGetVarIndex_Internal(exoid,EX_NODAL,vecname,&offsetN));
839   PetscCall(EXOGetVarIndex_Internal(exoid,EX_ELEM_BLOCK,vecname,&offsetZ));
840   PetscCheck(offsetN > 0 || offsetZ > 0,comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname);
841   if (offsetN > 0) PetscCall(VecLoadPlex_ExodusII_Nodal_Internal(v,exoid,(int) step+1,offsetN));
842   else if (offsetZ > 0) PetscCall(VecLoadPlex_ExodusII_Zonal_Internal(v,exoid,(int) step+1,offsetZ));
843   else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname);
844   PetscFunctionReturn(0);
845 }
846 
847 /*
848   VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file
849 
850   Collective on v
851 
852   Input Parameters:
853 + v  - The vector to be written
854 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
855 . step - the time step to write at (exodus steps are numbered starting from 1)
856 - offset - the location of the variable in the file
857 
858   Notes:
859   The exodus result nodal variable index is obtained by comparing the Vec name and the
860   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
861   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
862   amongst all nodal variables.
863 
864   Level: beginner
865 
866 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecLoadNodal_PlexEXO()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()`
867 @*/
868 PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset)
869 {
870   MPI_Comm           comm;
871   PetscMPIInt        size;
872   DM                 dm;
873   Vec                vNatural, vComp;
874   const PetscScalar *varray;
875   PetscInt           xs, xe, bs;
876   PetscBool          useNatural;
877 
878   PetscFunctionBegin;
879   PetscCall(PetscObjectGetComm((PetscObject) v, &comm));
880   PetscCallMPI(MPI_Comm_size(comm, &size));
881   PetscCall(VecGetDM(v, &dm));
882   PetscCall(DMGetUseNatural(dm, &useNatural));
883   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
884   if (useNatural) {
885     PetscCall(DMGetGlobalVector(dm, &vNatural));
886     PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural));
887     PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural));
888   } else {
889     vNatural = v;
890   }
891 
892   /* Write local chunk of the result in the exodus file
893      exodus stores each component of a vector-valued field as a separate variable.
894      We assume that they are stored sequentially */
895   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
896   PetscCall(VecGetBlockSize(vNatural, &bs));
897   if (bs == 1) {
898     PetscCall(VecGetArrayRead(vNatural, &varray));
899     PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);
900     PetscCall(VecRestoreArrayRead(vNatural, &varray));
901   } else {
902     IS       compIS;
903     PetscInt c;
904 
905     PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS));
906     for (c = 0; c < bs; ++c) {
907       PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs));
908       PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
909       PetscCall(VecGetArrayRead(vComp, &varray));
910       PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);
911       PetscCall(VecRestoreArrayRead(vComp, &varray));
912       PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
913     }
914     PetscCall(ISDestroy(&compIS));
915   }
916   if (useNatural) PetscCall(DMRestoreGlobalVector(dm, &vNatural));
917   PetscFunctionReturn(0);
918 }
919 
920 /*
921   VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file
922 
923   Collective on v
924 
925   Input Parameters:
926 + v  - The vector to be written
927 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
928 . step - the time step to read at (exodus steps are numbered starting from 1)
929 - offset - the location of the variable in the file
930 
931   Notes:
932   The exodus result nodal variable index is obtained by comparing the Vec name and the
933   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
934   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
935   amongst all nodal variables.
936 
937   Level: beginner
938 
939 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()`
940 */
941 PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset)
942 {
943   MPI_Comm       comm;
944   PetscMPIInt    size;
945   DM             dm;
946   Vec            vNatural, vComp;
947   PetscScalar   *varray;
948   PetscInt       xs, xe, bs;
949   PetscBool      useNatural;
950 
951   PetscFunctionBegin;
952   PetscCall(PetscObjectGetComm((PetscObject) v, &comm));
953   PetscCallMPI(MPI_Comm_size(comm, &size));
954   PetscCall(VecGetDM(v,&dm));
955   PetscCall(DMGetUseNatural(dm, &useNatural));
956   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
957   if (useNatural) PetscCall(DMGetGlobalVector(dm,&vNatural));
958   else            {vNatural = v;}
959 
960   /* Read local chunk from the file */
961   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
962   PetscCall(VecGetBlockSize(vNatural, &bs));
963   if (bs == 1) {
964     PetscCall(VecGetArray(vNatural, &varray));
965     PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);
966     PetscCall(VecRestoreArray(vNatural, &varray));
967   } else {
968     IS       compIS;
969     PetscInt c;
970 
971     PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS));
972     for (c = 0; c < bs; ++c) {
973       PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs));
974       PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
975       PetscCall(VecGetArray(vComp, &varray));
976       PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);
977       PetscCall(VecRestoreArray(vComp, &varray));
978       PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
979     }
980     PetscCall(ISDestroy(&compIS));
981   }
982   if (useNatural) {
983     PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v));
984     PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v));
985     PetscCall(DMRestoreGlobalVector(dm, &vNatural));
986   }
987   PetscFunctionReturn(0);
988 }
989 
990 /*
991   VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file
992 
993   Collective on v
994 
995   Input Parameters:
996 + v  - The vector to be written
997 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
998 . step - the time step to write at (exodus steps are numbered starting from 1)
999 - offset - the location of the variable in the file
1000 
1001   Notes:
1002   The exodus result zonal variable index is obtained by comparing the Vec name and the
1003   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
1004   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
1005   amongst all zonal variables.
1006 
1007   Level: beginner
1008 
1009 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()`
1010 */
1011 PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset)
1012 {
1013   MPI_Comm          comm;
1014   PetscMPIInt       size;
1015   DM                dm;
1016   Vec               vNatural, vComp;
1017   const PetscScalar *varray;
1018   PetscInt          xs, xe, bs;
1019   PetscBool         useNatural;
1020   IS                compIS;
1021   PetscInt         *csSize, *csID;
1022   PetscInt          numCS, set, csxs = 0;
1023 
1024   PetscFunctionBegin;
1025   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
1026   PetscCallMPI(MPI_Comm_size(comm, &size));
1027   PetscCall(VecGetDM(v, &dm));
1028   PetscCall(DMGetUseNatural(dm, &useNatural));
1029   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
1030   if (useNatural) {
1031     PetscCall(DMGetGlobalVector(dm, &vNatural));
1032     PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural));
1033     PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural));
1034   } else {
1035     vNatural = v;
1036   }
1037 
1038   /* Write local chunk of the result in the exodus file
1039      exodus stores each component of a vector-valued field as a separate variable.
1040      We assume that they are stored sequentially
1041      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
1042      but once the vector has been reordered to natural size, we cannot use the label information
1043      to figure out what to save where. */
1044   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
1045   PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize));
1046   PetscStackCallStandard(ex_get_ids,exoid, EX_ELEM_BLOCK, csID);
1047   for (set = 0; set < numCS; ++set) {
1048     ex_block block;
1049 
1050     block.id   = csID[set];
1051     block.type = EX_ELEM_BLOCK;
1052     PetscStackCallStandard(ex_get_block_param,exoid, &block);
1053     csSize[set] = block.num_entry;
1054   }
1055   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
1056   PetscCall(VecGetBlockSize(vNatural, &bs));
1057   if (bs > 1) PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS));
1058   for (set = 0; set < numCS; set++) {
1059     PetscInt csLocalSize, c;
1060 
1061     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
1062        local slice of zonal values:         xs/bs,xm/bs-1
1063        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
1064     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
1065     if (bs == 1) {
1066       PetscCall(VecGetArrayRead(vNatural, &varray));
1067       PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]);
1068       PetscCall(VecRestoreArrayRead(vNatural, &varray));
1069     } else {
1070       for (c = 0; c < bs; ++c) {
1071         PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs));
1072         PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
1073         PetscCall(VecGetArrayRead(vComp, &varray));
1074         PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]);
1075         PetscCall(VecRestoreArrayRead(vComp, &varray));
1076         PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
1077       }
1078     }
1079     csxs += csSize[set];
1080   }
1081   PetscCall(PetscFree2(csID, csSize));
1082   if (bs > 1) PetscCall(ISDestroy(&compIS));
1083   if (useNatural) PetscCall(DMRestoreGlobalVector(dm,&vNatural));
1084   PetscFunctionReturn(0);
1085 }
1086 
1087 /*
1088   VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file
1089 
1090   Collective on v
1091 
1092   Input Parameters:
1093 + v  - The vector to be written
1094 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
1095 . step - the time step to read at (exodus steps are numbered starting from 1)
1096 - offset - the location of the variable in the file
1097 
1098   Notes:
1099   The exodus result zonal variable index is obtained by comparing the Vec name and the
1100   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
1101   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
1102   amongst all zonal variables.
1103 
1104   Level: beginner
1105 
1106 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()`
1107 */
1108 PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset)
1109 {
1110   MPI_Comm          comm;
1111   PetscMPIInt       size;
1112   DM                dm;
1113   Vec               vNatural, vComp;
1114   PetscScalar      *varray;
1115   PetscInt          xs, xe, bs;
1116   PetscBool         useNatural;
1117   IS                compIS;
1118   PetscInt         *csSize, *csID;
1119   PetscInt          numCS, set, csxs = 0;
1120 
1121   PetscFunctionBegin;
1122   PetscCall(PetscObjectGetComm((PetscObject)v,&comm));
1123   PetscCallMPI(MPI_Comm_size(comm, &size));
1124   PetscCall(VecGetDM(v, &dm));
1125   PetscCall(DMGetUseNatural(dm, &useNatural));
1126   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
1127   if (useNatural) PetscCall(DMGetGlobalVector(dm,&vNatural));
1128   else            {vNatural = v;}
1129 
1130   /* Read local chunk of the result in the exodus file
1131      exodus stores each component of a vector-valued field as a separate variable.
1132      We assume that they are stored sequentially
1133      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
1134      but once the vector has been reordered to natural size, we cannot use the label information
1135      to figure out what to save where. */
1136   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
1137   PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize));
1138   PetscStackCallStandard(ex_get_ids,exoid, EX_ELEM_BLOCK, csID);
1139   for (set = 0; set < numCS; ++set) {
1140     ex_block block;
1141 
1142     block.id   = csID[set];
1143     block.type = EX_ELEM_BLOCK;
1144     PetscStackCallStandard(ex_get_block_param,exoid, &block);
1145     csSize[set] = block.num_entry;
1146   }
1147   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
1148   PetscCall(VecGetBlockSize(vNatural, &bs));
1149   if (bs > 1) PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS));
1150   for (set = 0; set < numCS; ++set) {
1151     PetscInt csLocalSize, c;
1152 
1153     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
1154        local slice of zonal values:         xs/bs,xm/bs-1
1155        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
1156     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
1157     if (bs == 1) {
1158       PetscCall(VecGetArray(vNatural, &varray));
1159       PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]);
1160       PetscCall(VecRestoreArray(vNatural, &varray));
1161     } else {
1162       for (c = 0; c < bs; ++c) {
1163         PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs));
1164         PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
1165         PetscCall(VecGetArray(vComp, &varray));
1166         PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]);
1167         PetscCall(VecRestoreArray(vComp, &varray));
1168         PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
1169       }
1170     }
1171     csxs += csSize[set];
1172   }
1173   PetscCall(PetscFree2(csID, csSize));
1174   if (bs > 1) PetscCall(ISDestroy(&compIS));
1175   if (useNatural) {
1176     PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v));
1177     PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v));
1178     PetscCall(DMRestoreGlobalVector(dm, &vNatural));
1179   }
1180   PetscFunctionReturn(0);
1181 }
1182 #endif
1183 
1184 /*@
1185   PetscViewerExodusIIGetId - Get the file id of the ExodusII file
1186 
1187   Logically Collective on PetscViewer
1188 
1189   Input Parameter:
1190 .  viewer - the PetscViewer
1191 
1192   Output Parameter:
1193 .  exoid - The ExodusII file id
1194 
1195   Level: intermediate
1196 
1197 .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
1198 @*/
1199 PetscErrorCode PetscViewerExodusIIGetId(PetscViewer viewer, int *exoid)
1200 {
1201   PetscFunctionBegin;
1202   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1203   PetscTryMethod(viewer, "PetscViewerGetId_C",(PetscViewer,int*),(viewer,exoid));
1204   PetscFunctionReturn(0);
1205 }
1206 
1207 /*@
1208    PetscViewerExodusIISetOrder - Set the elements order in the exodusII file.
1209 
1210    Collective
1211 
1212    Input Parameters:
1213 +  viewer - the viewer
1214 -  order - elements order
1215 
1216    Output Parameter:
1217 
1218    Level: beginner
1219 
1220    Note:
1221 
1222 .seealso: `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()`
1223 @*/
1224 PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer viewer, PetscInt order)
1225 {
1226   PetscFunctionBegin;
1227   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1228   PetscTryMethod(viewer, "PetscViewerSetOrder_C",(PetscViewer,PetscInt),(viewer,order));
1229   PetscFunctionReturn(0);
1230 }
1231 
1232 /*@
1233    PetscViewerExodusIIGetOrder - Get the elements order in the exodusII file.
1234 
1235    Collective
1236 
1237    Input Parameters:
1238 +  viewer - the viewer
1239 -  order - elements order
1240 
1241    Output Parameter:
1242 
1243    Level: beginner
1244 
1245    Note:
1246 
1247 .seealso: `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()`
1248 @*/
1249 PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer viewer, PetscInt *order)
1250 {
1251   PetscFunctionBegin;
1252   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1253   PetscTryMethod(viewer, "PetscViewerGetOrder_C",(PetscViewer,PetscInt*),(viewer,order));
1254   PetscFunctionReturn(0);
1255 }
1256 
1257 /*@C
1258    PetscViewerExodusIIOpen - Opens a file for ExodusII input/output.
1259 
1260    Collective
1261 
1262    Input Parameters:
1263 +  comm - MPI communicator
1264 .  name - name of file
1265 -  type - type of file
1266 $    FILE_MODE_WRITE - create new file for binary output
1267 $    FILE_MODE_READ - open existing file for binary input
1268 $    FILE_MODE_APPEND - open existing file for binary output
1269 
1270    Output Parameter:
1271 .  exo - PetscViewer for Exodus II input/output to use with the specified file
1272 
1273    Level: beginner
1274 
1275    Note:
1276    This PetscViewer should be destroyed with PetscViewerDestroy().
1277 
1278 .seealso: `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1279           `DMLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
1280 @*/
1281 PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo)
1282 {
1283   PetscFunctionBegin;
1284   PetscCall(PetscViewerCreate(comm, exo));
1285   PetscCall(PetscViewerSetType(*exo, PETSCVIEWEREXODUSII));
1286   PetscCall(PetscViewerFileSetMode(*exo, type));
1287   PetscCall(PetscViewerFileSetName(*exo, name));
1288   PetscCall(PetscViewerSetFromOptions(*exo));
1289   PetscFunctionReturn(0);
1290 }
1291 
1292 /*@C
1293   DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file.
1294 
1295   Collective
1296 
1297   Input Parameters:
1298 + comm  - The MPI communicator
1299 . filename - The name of the ExodusII file
1300 - interpolate - Create faces and edges in the mesh
1301 
1302   Output Parameter:
1303 . dm  - The DM object representing the mesh
1304 
1305   Level: beginner
1306 
1307 .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateExodus()`
1308 @*/
1309 PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
1310 {
1311   PetscMPIInt    rank;
1312 #if defined(PETSC_HAVE_EXODUSII)
1313   int   CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1;
1314   float version;
1315 #endif
1316 
1317   PetscFunctionBegin;
1318   PetscValidCharPointer(filename, 2);
1319   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1320 #if defined(PETSC_HAVE_EXODUSII)
1321   if (rank == 0) {
1322     exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version);
1323     PetscCheck(exoid > 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename);
1324   }
1325   PetscCall(DMPlexCreateExodus(comm, exoid, interpolate, dm));
1326   if (rank == 0) {PetscStackCallStandard(ex_close,exoid);}
1327   PetscFunctionReturn(0);
1328 #else
1329   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
1330 #endif
1331 }
1332 
1333 #if defined(PETSC_HAVE_EXODUSII)
1334 static PetscErrorCode ExodusGetCellType_Internal(const char *elem_type, DMPolytopeType *ct)
1335 {
1336   PetscBool      flg;
1337 
1338   PetscFunctionBegin;
1339   *ct = DM_POLYTOPE_UNKNOWN;
1340   PetscCall(PetscStrcmp(elem_type, "TRI", &flg));
1341   if (flg) {*ct = DM_POLYTOPE_TRIANGLE; goto done;}
1342   PetscCall(PetscStrcmp(elem_type, "TRI3", &flg));
1343   if (flg) {*ct = DM_POLYTOPE_TRIANGLE; goto done;}
1344   PetscCall(PetscStrcmp(elem_type, "QUAD", &flg));
1345   if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;}
1346   PetscCall(PetscStrcmp(elem_type, "QUAD4", &flg));
1347   if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;}
1348   PetscCall(PetscStrcmp(elem_type, "SHELL4", &flg));
1349   if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;}
1350   PetscCall(PetscStrcmp(elem_type, "TETRA", &flg));
1351   if (flg) {*ct = DM_POLYTOPE_TETRAHEDRON; goto done;}
1352   PetscCall(PetscStrcmp(elem_type, "TET4", &flg));
1353   if (flg) {*ct = DM_POLYTOPE_TETRAHEDRON; goto done;}
1354   PetscCall(PetscStrcmp(elem_type, "WEDGE", &flg));
1355   if (flg) {*ct = DM_POLYTOPE_TRI_PRISM; goto done;}
1356   PetscCall(PetscStrcmp(elem_type, "HEX", &flg));
1357   if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;}
1358   PetscCall(PetscStrcmp(elem_type, "HEX8", &flg));
1359   if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;}
1360   PetscCall(PetscStrcmp(elem_type, "HEXAHEDRON", &flg));
1361   if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;}
1362   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized element type %s", elem_type);
1363   done:
1364   PetscFunctionReturn(0);
1365 }
1366 #endif
1367 
1368 /*@
1369   DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID.
1370 
1371   Collective
1372 
1373   Input Parameters:
1374 + comm  - The MPI communicator
1375 . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
1376 - interpolate - Create faces and edges in the mesh
1377 
1378   Output Parameter:
1379 . dm  - The DM object representing the mesh
1380 
1381   Level: beginner
1382 
1383 .seealso: `DMPLEX`, `DMCreate()`
1384 @*/
1385 PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm)
1386 {
1387 #if defined(PETSC_HAVE_EXODUSII)
1388   PetscMPIInt    num_proc, rank;
1389   DMLabel        cellSets = NULL, faceSets = NULL, vertSets = NULL;
1390   PetscSection   coordSection;
1391   Vec            coordinates;
1392   PetscScalar    *coords;
1393   PetscInt       coordSize, v;
1394   /* Read from ex_get_init() */
1395   char title[PETSC_MAX_PATH_LEN+1];
1396   int  dim    = 0, dimEmbed = 0, numVertices = 0, numCells = 0;
1397   int  num_cs = 0, num_vs = 0, num_fs = 0;
1398 #endif
1399 
1400   PetscFunctionBegin;
1401 #if defined(PETSC_HAVE_EXODUSII)
1402   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1403   PetscCallMPI(MPI_Comm_size(comm, &num_proc));
1404   PetscCall(DMCreate(comm, dm));
1405   PetscCall(DMSetType(*dm, DMPLEX));
1406   /* Open EXODUS II file and read basic information on rank 0, then broadcast to all processors */
1407   if (rank == 0) {
1408     PetscCall(PetscMemzero(title,PETSC_MAX_PATH_LEN+1));
1409     PetscStackCallStandard(ex_get_init,exoid, title, &dimEmbed, &numVertices, &numCells, &num_cs, &num_vs, &num_fs);
1410     PetscCheck(num_cs,PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set");
1411   }
1412   PetscCallMPI(MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm));
1413   PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm));
1414   PetscCall(PetscObjectSetName((PetscObject) *dm, title));
1415   PetscCall(DMPlexSetChart(*dm, 0, numCells+numVertices));
1416   /*   We do not want this label automatically computed, instead we compute it here */
1417   PetscCall(DMCreateLabel(*dm, "celltype"));
1418 
1419   /* Read cell sets information */
1420   if (rank == 0) {
1421     PetscInt *cone;
1422     int      c, cs, ncs, c_loc, v, v_loc;
1423     /* Read from ex_get_elem_blk_ids() */
1424     int *cs_id, *cs_order;
1425     /* Read from ex_get_elem_block() */
1426     char buffer[PETSC_MAX_PATH_LEN+1];
1427     int  num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr;
1428     /* Read from ex_get_elem_conn() */
1429     int *cs_connect;
1430 
1431     /* Get cell sets IDs */
1432     PetscCall(PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order));
1433     PetscStackCallStandard(ex_get_ids,exoid, EX_ELEM_BLOCK, cs_id);
1434     /* Read the cell set connectivity table and build mesh topology
1435        EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
1436     /* Check for a hybrid mesh */
1437     for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) {
1438       DMPolytopeType ct;
1439       char           elem_type[PETSC_MAX_PATH_LEN];
1440 
1441       PetscCall(PetscArrayzero(elem_type, sizeof(elem_type)));
1442       PetscStackCallStandard(ex_get_elem_type,exoid, cs_id[cs], elem_type);
1443       PetscCall(ExodusGetCellType_Internal(elem_type, &ct));
1444       dim  = PetscMax(dim, DMPolytopeTypeGetDim(ct));
1445       PetscStackCallStandard(ex_get_block,exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr);
1446       switch (ct) {
1447         case DM_POLYTOPE_TRI_PRISM:
1448           cs_order[cs] = cs;
1449           ++num_hybrid;
1450           break;
1451         default:
1452           for (c = cs; c > cs-num_hybrid; --c) cs_order[c] = cs_order[c-1];
1453           cs_order[cs-num_hybrid] = cs;
1454       }
1455     }
1456     /* First set sizes */
1457     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
1458       DMPolytopeType ct;
1459       char           elem_type[PETSC_MAX_PATH_LEN];
1460       const PetscInt cs = cs_order[ncs];
1461 
1462       PetscCall(PetscArrayzero(elem_type, sizeof(elem_type)));
1463       PetscStackCallStandard(ex_get_elem_type,exoid, cs_id[cs], elem_type);
1464       PetscCall(ExodusGetCellType_Internal(elem_type, &ct));
1465       PetscStackCallStandard(ex_get_block,exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr);
1466       for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
1467         PetscCall(DMPlexSetConeSize(*dm, c, num_vertex_per_cell));
1468         PetscCall(DMPlexSetCellType(*dm, c, ct));
1469       }
1470     }
1471     for (v = numCells; v < numCells+numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
1472     PetscCall(DMSetUp(*dm));
1473     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
1474       const PetscInt cs = cs_order[ncs];
1475       PetscStackCallStandard(ex_get_block,exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr);
1476       PetscCall(PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone));
1477       PetscStackCallStandard(ex_get_conn,exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL);
1478       /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
1479       for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
1480         DMPolytopeType ct;
1481 
1482         for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) {
1483           cone[v_loc] = cs_connect[v]+numCells-1;
1484         }
1485         PetscCall(DMPlexGetCellType(*dm, c, &ct));
1486         PetscCall(DMPlexInvertCell(ct, cone));
1487         PetscCall(DMPlexSetCone(*dm, c, cone));
1488         PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", c, cs_id[cs]));
1489       }
1490       PetscCall(PetscFree2(cs_connect,cone));
1491     }
1492     PetscCall(PetscFree2(cs_id, cs_order));
1493   }
1494   {
1495     PetscInt ints[] = {dim, dimEmbed};
1496 
1497     PetscCallMPI(MPI_Bcast(ints, 2, MPIU_INT, 0, comm));
1498     PetscCall(DMSetDimension(*dm, ints[0]));
1499     PetscCall(DMSetCoordinateDim(*dm, ints[1]));
1500     dim      = ints[0];
1501     dimEmbed = ints[1];
1502   }
1503   PetscCall(DMPlexSymmetrize(*dm));
1504   PetscCall(DMPlexStratify(*dm));
1505   if (interpolate) {
1506     DM idm;
1507 
1508     PetscCall(DMPlexInterpolate(*dm, &idm));
1509     PetscCall(DMDestroy(dm));
1510     *dm  = idm;
1511   }
1512 
1513   /* Create vertex set label */
1514   if (rank == 0 && (num_vs > 0)) {
1515     int vs, v;
1516     /* Read from ex_get_node_set_ids() */
1517     int *vs_id;
1518     /* Read from ex_get_node_set_param() */
1519     int num_vertex_in_set;
1520     /* Read from ex_get_node_set() */
1521     int *vs_vertex_list;
1522 
1523     /* Get vertex set ids */
1524     PetscCall(PetscMalloc1(num_vs, &vs_id));
1525     PetscStackCallStandard(ex_get_ids,exoid, EX_NODE_SET, vs_id);
1526     for (vs = 0; vs < num_vs; ++vs) {
1527       PetscStackCallStandard(ex_get_set_param,exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL);
1528       PetscCall(PetscMalloc1(num_vertex_in_set, &vs_vertex_list));
1529       PetscStackCallStandard(ex_get_set,exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL);
1530       for (v = 0; v < num_vertex_in_set; ++v) {
1531         PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]));
1532       }
1533       PetscCall(PetscFree(vs_vertex_list));
1534     }
1535     PetscCall(PetscFree(vs_id));
1536   }
1537   /* Read coordinates */
1538   PetscCall(DMGetCoordinateSection(*dm, &coordSection));
1539   PetscCall(PetscSectionSetNumFields(coordSection, 1));
1540   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dimEmbed));
1541   PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
1542   for (v = numCells; v < numCells+numVertices; ++v) {
1543     PetscCall(PetscSectionSetDof(coordSection, v, dimEmbed));
1544     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed));
1545   }
1546   PetscCall(PetscSectionSetUp(coordSection));
1547   PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
1548   PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
1549   PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates"));
1550   PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
1551   PetscCall(VecSetBlockSize(coordinates, dimEmbed));
1552   PetscCall(VecSetType(coordinates,VECSTANDARD));
1553   PetscCall(VecGetArray(coordinates, &coords));
1554   if (rank == 0) {
1555     PetscReal *x, *y, *z;
1556 
1557     PetscCall(PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z));
1558     PetscStackCallStandard(ex_get_coord,exoid, x, y, z);
1559     if (dimEmbed > 0) {
1560       for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+0] = x[v];
1561     }
1562     if (dimEmbed > 1) {
1563       for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+1] = y[v];
1564     }
1565     if (dimEmbed > 2) {
1566       for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+2] = z[v];
1567     }
1568     PetscCall(PetscFree3(x,y,z));
1569   }
1570   PetscCall(VecRestoreArray(coordinates, &coords));
1571   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
1572   PetscCall(VecDestroy(&coordinates));
1573 
1574   /* Create side set label */
1575   if (rank == 0 && interpolate && (num_fs > 0)) {
1576     int fs, f, voff;
1577     /* Read from ex_get_side_set_ids() */
1578     int *fs_id;
1579     /* Read from ex_get_side_set_param() */
1580     int num_side_in_set;
1581     /* Read from ex_get_side_set_node_list() */
1582     int *fs_vertex_count_list, *fs_vertex_list;
1583     /* Read side set labels */
1584     char fs_name[MAX_STR_LENGTH+1];
1585 
1586     /* Get side set ids */
1587     PetscCall(PetscMalloc1(num_fs, &fs_id));
1588     PetscStackCallStandard(ex_get_ids,exoid, EX_SIDE_SET, fs_id);
1589     for (fs = 0; fs < num_fs; ++fs) {
1590       PetscStackCallStandard(ex_get_set_param,exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL);
1591       PetscCall(PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list));
1592       PetscStackCallStandard(ex_get_side_set_node_list,exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);
1593       /* Get the specific name associated with this side set ID. */
1594       int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name);
1595       for (f = 0, voff = 0; f < num_side_in_set; ++f) {
1596         const PetscInt *faces   = NULL;
1597         PetscInt       faceSize = fs_vertex_count_list[f], numFaces;
1598         PetscInt       faceVertices[4], v;
1599 
1600         PetscCheck(faceSize <= 4,comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize);
1601         for (v = 0; v < faceSize; ++v, ++voff) {
1602           faceVertices[v] = fs_vertex_list[voff]+numCells-1;
1603         }
1604         PetscCall(DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces));
1605         PetscCheck(numFaces == 1,comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces);
1606         PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", faces[0], fs_id[fs]));
1607         /* Only add the label if one has been detected for this side set. */
1608         if (!fs_name_err) {
1609           PetscCall(DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]));
1610         }
1611         PetscCall(DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces));
1612       }
1613       PetscCall(PetscFree2(fs_vertex_count_list,fs_vertex_list));
1614     }
1615     PetscCall(PetscFree(fs_id));
1616   }
1617 
1618   { /* Create Cell/Face/Vertex Sets labels at all processes */
1619     enum {n = 3};
1620     PetscBool flag[n];
1621 
1622     flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE;
1623     flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE;
1624     flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE;
1625     PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm));
1626     if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets"));
1627     if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets"));
1628     if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets"));
1629   }
1630   PetscFunctionReturn(0);
1631 #else
1632   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
1633 #endif
1634 }
1635