xref: /petsc/src/dm/impls/plex/plexexodusii.c (revision 1e50132f2fca1f8b659af561d8e08fe7445712a1)
1552f7358SJed Brown #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown 
4519f805aSKarl Rupp #if defined(PETSC_HAVE_EXODUSII)
5552f7358SJed Brown #include <netcdf.h>
6552f7358SJed Brown #include <exodusII.h>
7552f7358SJed Brown #endif
8552f7358SJed Brown 
9*1e50132fSMatthew G. Knepley #include <petsc/private/viewerimpl.h>
10*1e50132fSMatthew G. Knepley #include <petsc/private/viewerexodusiiimpl.h>
11e0266925SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
1278569bb4SMatthew G. Knepley /*
13*1e50132fSMatthew G. Knepley   PETSC_VIEWER_EXODUSII_ - Creates an ExodusII PetscViewer shared by all processors in a communicator.
14*1e50132fSMatthew G. Knepley 
15*1e50132fSMatthew G. Knepley   Collective
16*1e50132fSMatthew G. Knepley 
17*1e50132fSMatthew G. Knepley   Input Parameter:
18*1e50132fSMatthew G. Knepley . comm - the MPI communicator to share the ExodusII PetscViewer
19*1e50132fSMatthew G. Knepley 
20*1e50132fSMatthew G. Knepley   Level: intermediate
21*1e50132fSMatthew G. Knepley 
22*1e50132fSMatthew G. Knepley   Notes:
23*1e50132fSMatthew G. Knepley     misses Fortran bindings
24*1e50132fSMatthew G. Knepley 
25*1e50132fSMatthew G. Knepley   Notes:
26*1e50132fSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_EXODUSII_ does not return
27*1e50132fSMatthew G. Knepley   an error code.  The GLVIS PetscViewer is usually used in the form
28*1e50132fSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_EXODUSII_(comm));
29*1e50132fSMatthew G. Knepley 
30*1e50132fSMatthew G. Knepley .seealso: PetscViewerGLVISOpen(), PetscViewerGLVisType, PetscViewerCreate(), PetscViewerDestroy()
31*1e50132fSMatthew G. Knepley */
32*1e50132fSMatthew G. Knepley PetscViewer PETSC_VIEWER_EXODUSII_(MPI_Comm comm)
33*1e50132fSMatthew G. Knepley {
34*1e50132fSMatthew G. Knepley   PetscViewer    viewer;
35*1e50132fSMatthew G. Knepley   PetscErrorCode ierr;
36*1e50132fSMatthew G. Knepley 
37*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
38*1e50132fSMatthew G. Knepley   ierr = PetscViewerExodusIIOpen(comm, "mesh.exo", FILE_MODE_WRITE, &viewer);
39*1e50132fSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
40*1e50132fSMatthew G. Knepley   ierr = PetscObjectRegisterDestroy((PetscObject) viewer);
41*1e50132fSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
42*1e50132fSMatthew G. Knepley   PetscFunctionReturn(viewer);
43*1e50132fSMatthew G. Knepley }
44*1e50132fSMatthew G. Knepley 
45*1e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerView_ExodusII(PetscViewer v, PetscViewer viewer)
46*1e50132fSMatthew G. Knepley {
47*1e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
48*1e50132fSMatthew G. Knepley   PetscErrorCode        ierr;
49*1e50132fSMatthew G. Knepley 
50*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
51*1e50132fSMatthew G. Knepley   if (exo->filename) {ierr = PetscViewerASCIIPrintf(viewer, "Filename: %s\n", exo->filename);CHKERRQ(ierr);}
52*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
53*1e50132fSMatthew G. Knepley }
54*1e50132fSMatthew G. Knepley 
55*1e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerSetFromOptions_ExodusII(PetscOptionItems *PetscOptionsObject, PetscViewer v)
56*1e50132fSMatthew G. Knepley {
57*1e50132fSMatthew G. Knepley   PetscErrorCode ierr;
58*1e50132fSMatthew G. Knepley 
59*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
60*1e50132fSMatthew G. Knepley   ierr = PetscOptionsHead(PetscOptionsObject, "ExodusII PetscViewer Options");CHKERRQ(ierr);
61*1e50132fSMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
62*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
63*1e50132fSMatthew G. Knepley }
64*1e50132fSMatthew G. Knepley 
65*1e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerSetUp_ExodusII(PetscViewer viewer)
66*1e50132fSMatthew G. Knepley {
67*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
68*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
69*1e50132fSMatthew G. Knepley }
70*1e50132fSMatthew G. Knepley 
71*1e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerDestroy_ExodusII(PetscViewer viewer)
72*1e50132fSMatthew G. Knepley {
73*1e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
74*1e50132fSMatthew G. Knepley   PetscErrorCode        ierr;
75*1e50132fSMatthew G. Knepley 
76*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
77*1e50132fSMatthew G. Knepley   if (exo->exoid >= 0) {ierr = ex_close(exo->exoid);CHKERRQ(ierr);}
78*1e50132fSMatthew G. Knepley   ierr = PetscFree(exo);CHKERRQ(ierr);
79*1e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
80*1e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
81*1e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
82*1e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerExodusIIGetId",NULL);CHKERRQ(ierr);
83*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
84*1e50132fSMatthew G. Knepley }
85*1e50132fSMatthew G. Knepley 
86*1e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileSetName_ExodusII(PetscViewer viewer, const char name[])
87*1e50132fSMatthew G. Knepley {
88*1e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
89*1e50132fSMatthew G. Knepley   PetscMPIInt           rank;
90*1e50132fSMatthew G. Knepley   int                   CPU_word_size, IO_word_size, EXO_mode;
91*1e50132fSMatthew G. Knepley   PetscErrorCode        ierr;
92*1e50132fSMatthew G. Knepley 
93*1e50132fSMatthew G. Knepley 
94*1e50132fSMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) viewer), &rank);CHKERRQ(ierr);
95*1e50132fSMatthew G. Knepley   CPU_word_size = sizeof(PetscReal);
96*1e50132fSMatthew G. Knepley   IO_word_size  = sizeof(PetscReal);
97*1e50132fSMatthew G. Knepley 
98*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
99*1e50132fSMatthew G. Knepley   if (exo->exoid >= 0) {ex_close(exo->exoid); exo->exoid = -1;}
100*1e50132fSMatthew G. Knepley   if (exo->filename) {ierr = PetscFree(exo->filename);CHKERRQ(ierr);}
101*1e50132fSMatthew G. Knepley   ierr = PetscStrallocpy(name, &exo->filename);CHKERRQ(ierr);
102*1e50132fSMatthew G. Knepley   /* Create or open the file collectively */
103*1e50132fSMatthew G. Knepley   switch (exo->btype) {
104*1e50132fSMatthew G. Knepley   case FILE_MODE_READ:
105*1e50132fSMatthew G. Knepley     EXO_mode = EX_CLOBBER;
106*1e50132fSMatthew G. Knepley     break;
107*1e50132fSMatthew G. Knepley   case FILE_MODE_APPEND:
108*1e50132fSMatthew G. Knepley     EXO_mode = EX_CLOBBER;
109*1e50132fSMatthew G. Knepley     break;
110*1e50132fSMatthew G. Knepley   case FILE_MODE_WRITE:
111*1e50132fSMatthew G. Knepley     EXO_mode = EX_CLOBBER;
112*1e50132fSMatthew G. Knepley     break;
113*1e50132fSMatthew G. Knepley   default:
114*1e50132fSMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
115*1e50132fSMatthew G. Knepley   }
116*1e50132fSMatthew G. Knepley   #if defined(PETSC_USE_64BIT_INDICES)
117*1e50132fSMatthew G. Knepley   EXO_mode += EX_ALL_INT64_API;
118*1e50132fSMatthew G. Knepley   #endif
119*1e50132fSMatthew G. Knepley   exo->exoid = ex_create(name, EXO_mode, &CPU_word_size, &IO_word_size);
120*1e50132fSMatthew G. Knepley   if (exo->exoid < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", name);
121*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
122*1e50132fSMatthew G. Knepley }
123*1e50132fSMatthew G. Knepley 
124*1e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileGetName_ExodusII(PetscViewer viewer, const char **name)
125*1e50132fSMatthew G. Knepley {
126*1e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
127*1e50132fSMatthew G. Knepley 
128*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
129*1e50132fSMatthew G. Knepley   *name = exo->filename;
130*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
131*1e50132fSMatthew G. Knepley }
132*1e50132fSMatthew G. Knepley 
133*1e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileSetMode_ExodusII(PetscViewer viewer, PetscFileMode type)
134*1e50132fSMatthew G. Knepley {
135*1e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
136*1e50132fSMatthew G. Knepley 
137*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
138*1e50132fSMatthew G. Knepley   exo->btype = type;
139*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
140*1e50132fSMatthew G. Knepley }
141*1e50132fSMatthew G. Knepley 
142*1e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileGetMode_ExodusII(PetscViewer viewer, PetscFileMode *type)
143*1e50132fSMatthew G. Knepley {
144*1e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
145*1e50132fSMatthew G. Knepley 
146*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
147*1e50132fSMatthew G. Knepley   *type = exo->btype;
148*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
149*1e50132fSMatthew G. Knepley }
150*1e50132fSMatthew G. Knepley 
151*1e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer, int *exoid)
152*1e50132fSMatthew G. Knepley {
153*1e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
154*1e50132fSMatthew G. Knepley 
155*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
156*1e50132fSMatthew G. Knepley   *exoid = exo->exoid;
157*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
158*1e50132fSMatthew G. Knepley }
159*1e50132fSMatthew G. Knepley 
160*1e50132fSMatthew G. Knepley /*@C
161*1e50132fSMatthew G. Knepley    PetscViewerExodusIIOpen - Opens a file for ExodusII input/output.
162*1e50132fSMatthew G. Knepley 
163*1e50132fSMatthew G. Knepley    Collective
164*1e50132fSMatthew G. Knepley 
165*1e50132fSMatthew G. Knepley    Input Parameters:
166*1e50132fSMatthew G. Knepley +  comm - MPI communicator
167*1e50132fSMatthew G. Knepley .  name - name of file
168*1e50132fSMatthew G. Knepley -  type - type of file
169*1e50132fSMatthew G. Knepley $    FILE_MODE_WRITE - create new file for binary output
170*1e50132fSMatthew G. Knepley $    FILE_MODE_READ - open existing file for binary input
171*1e50132fSMatthew G. Knepley $    FILE_MODE_APPEND - open existing file for binary output
172*1e50132fSMatthew G. Knepley 
173*1e50132fSMatthew G. Knepley    Output Parameter:
174*1e50132fSMatthew G. Knepley .  exo - PetscViewer for HDF5 input/output to use with the specified file
175*1e50132fSMatthew G. Knepley 
176*1e50132fSMatthew G. Knepley    Level: beginner
177*1e50132fSMatthew G. Knepley 
178*1e50132fSMatthew G. Knepley    Note:
179*1e50132fSMatthew G. Knepley    This PetscViewer should be destroyed with PetscViewerDestroy().
180*1e50132fSMatthew G. Knepley 
181*1e50132fSMatthew G. Knepley 
182*1e50132fSMatthew G. Knepley .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
183*1e50132fSMatthew G. Knepley           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
184*1e50132fSMatthew G. Knepley           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
185*1e50132fSMatthew G. Knepley @*/
186*1e50132fSMatthew G. Knepley PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo)
187*1e50132fSMatthew G. Knepley {
188*1e50132fSMatthew G. Knepley   PetscErrorCode ierr;
189*1e50132fSMatthew G. Knepley 
190*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
191*1e50132fSMatthew G. Knepley   ierr = PetscViewerCreate(comm, exo);CHKERRQ(ierr);
192*1e50132fSMatthew G. Knepley   ierr = PetscViewerSetType(*exo, PETSCVIEWEREXODUSII);CHKERRQ(ierr);
193*1e50132fSMatthew G. Knepley   ierr = PetscViewerFileSetMode(*exo, type);CHKERRQ(ierr);
194*1e50132fSMatthew G. Knepley   ierr = PetscViewerFileSetName(*exo, name);CHKERRQ(ierr);
195*1e50132fSMatthew G. Knepley   ierr = PetscViewerSetFromOptions(*exo);CHKERRQ(ierr);
196*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
197*1e50132fSMatthew G. Knepley }
198*1e50132fSMatthew G. Knepley 
199*1e50132fSMatthew G. Knepley /*MC
200*1e50132fSMatthew G. Knepley    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
201*1e50132fSMatthew G. Knepley 
202*1e50132fSMatthew G. Knepley 
203*1e50132fSMatthew G. Knepley .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
204*1e50132fSMatthew G. Knepley            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
205*1e50132fSMatthew G. Knepley            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
206*1e50132fSMatthew G. Knepley            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
207*1e50132fSMatthew G. Knepley 
208*1e50132fSMatthew G. Knepley   Level: beginner
209*1e50132fSMatthew G. Knepley M*/
210*1e50132fSMatthew G. Knepley 
211*1e50132fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerCreate_ExodusII(PetscViewer v)
212*1e50132fSMatthew G. Knepley {
213*1e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo;
214*1e50132fSMatthew G. Knepley   PetscErrorCode        ierr;
215*1e50132fSMatthew G. Knepley 
216*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
217*1e50132fSMatthew G. Knepley   ierr = PetscNewLog(v,&exo);CHKERRQ(ierr);
218*1e50132fSMatthew G. Knepley 
219*1e50132fSMatthew G. Knepley   v->data                = (void*) exo;
220*1e50132fSMatthew G. Knepley   v->ops->destroy        = PetscViewerDestroy_ExodusII;
221*1e50132fSMatthew G. Knepley   v->ops->setfromoptions = PetscViewerSetFromOptions_ExodusII;
222*1e50132fSMatthew G. Knepley   v->ops->setup          = PetscViewerSetUp_ExodusII;
223*1e50132fSMatthew G. Knepley   v->ops->view           = PetscViewerView_ExodusII;
224*1e50132fSMatthew G. Knepley   v->ops->flush          = 0;
225*1e50132fSMatthew G. Knepley   exo->btype             = (PetscFileMode) -1;
226*1e50132fSMatthew G. Knepley   exo->filename          = 0;
227*1e50132fSMatthew G. Knepley   exo->exoid             = -1;
228*1e50132fSMatthew G. Knepley 
229*1e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_ExodusII);CHKERRQ(ierr);
230*1e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_ExodusII);CHKERRQ(ierr);
231*1e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_ExodusII);CHKERRQ(ierr);
232*1e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_ExodusII);CHKERRQ(ierr);
233*1e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerExodusIIGetId_C",PetscViewerExodusIIGetId_ExodusII);CHKERRQ(ierr);
234*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
235*1e50132fSMatthew G. Knepley }
236*1e50132fSMatthew G. Knepley 
237*1e50132fSMatthew G. Knepley /*
23878569bb4SMatthew G. Knepley   EXOGetVarIndex - Locate a result in an exodus file based on its name
23978569bb4SMatthew G. Knepley 
24078569bb4SMatthew G. Knepley   Not collective
24178569bb4SMatthew G. Knepley 
24278569bb4SMatthew G. Knepley   Input Parameters:
24378569bb4SMatthew G. Knepley + exoid    - the exodus id of a file (obtained from ex_open or ex_create for instance)
24478569bb4SMatthew G. Knepley . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK
24578569bb4SMatthew G. Knepley - name     - the name of the result
24678569bb4SMatthew G. Knepley 
24778569bb4SMatthew G. Knepley   Output Parameters:
24878569bb4SMatthew G. Knepley . varIndex - the location in the exodus file of the result
24978569bb4SMatthew G. Knepley 
25078569bb4SMatthew G. Knepley   Notes:
25178569bb4SMatthew G. Knepley   The exodus variable index is obtained by comparing name and the
25278569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance if name is "V"
25378569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
25478569bb4SMatthew G. Knepley   amongst all variables of type obj_type.
25578569bb4SMatthew G. Knepley 
25678569bb4SMatthew G. Knepley   Level: beginner
25778569bb4SMatthew G. Knepley 
258cd921712SStefano Zampini .seealso: EXOGetVarIndex(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadNodal_PlexEXO(),VecLoadZonal_PlexEXO()
25978569bb4SMatthew G. Knepley */
26078569bb4SMatthew G. Knepley static PetscErrorCode EXOGetVarIndex_Private(int exoid, ex_entity_type obj_type, const char name[], int *varIndex)
26178569bb4SMatthew G. Knepley {
2626de834b4SMatthew G. Knepley   int            num_vars, i, j;
26378569bb4SMatthew G. Knepley   char           ext_name[MAX_STR_LENGTH+1], var_name[MAX_STR_LENGTH+1];
26478569bb4SMatthew G. Knepley   const int      num_suffix = 5;
26578569bb4SMatthew G. Knepley   char          *suffix[5];
26678569bb4SMatthew G. Knepley   PetscBool      flg;
26778569bb4SMatthew G. Knepley   PetscErrorCode ierr;
26878569bb4SMatthew G. Knepley 
26978569bb4SMatthew G. Knepley   PetscFunctionBegin;
27078569bb4SMatthew G. Knepley   suffix[0] = (char *) "";
27178569bb4SMatthew G. Knepley   suffix[1] = (char *) "_X";
27278569bb4SMatthew G. Knepley   suffix[2] = (char *) "_XX";
27378569bb4SMatthew G. Knepley   suffix[3] = (char *) "_1";
27478569bb4SMatthew G. Knepley   suffix[4] = (char *) "_11";
27578569bb4SMatthew G. Knepley 
27678569bb4SMatthew G. Knepley   *varIndex = 0;
2776de834b4SMatthew G. Knepley   PetscStackCallStandard(ex_get_variable_param,(exoid, obj_type, &num_vars));
27878569bb4SMatthew G. Knepley   for (i = 0; i < num_vars; ++i) {
2796de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_variable_name,(exoid, obj_type, i+1, var_name));
28078569bb4SMatthew G. Knepley     for (j = 0; j < num_suffix; ++j){
28178569bb4SMatthew G. Knepley       ierr = PetscStrncpy(ext_name, name, MAX_STR_LENGTH);CHKERRQ(ierr);
282a126751eSBarry Smith       ierr = PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH);CHKERRQ(ierr);
28378569bb4SMatthew G. Knepley       ierr = PetscStrcasecmp(ext_name, var_name, &flg);
28478569bb4SMatthew G. Knepley       if (flg) {
28578569bb4SMatthew G. Knepley         *varIndex = i+1;
28678569bb4SMatthew G. Knepley         PetscFunctionReturn(0);
28778569bb4SMatthew G. Knepley       }
28878569bb4SMatthew G. Knepley     }
28978569bb4SMatthew G. Knepley   }
29078569bb4SMatthew G. Knepley   SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to locate variable %s in ExodusII file.", name);
29178569bb4SMatthew G. Knepley  PetscFunctionReturn(-1);
29278569bb4SMatthew G. Knepley }
29378569bb4SMatthew G. Knepley 
29478569bb4SMatthew G. Knepley /*
29578569bb4SMatthew G. Knepley   DMPlexView_ExodusII_Internal - Write a DM to disk in exodus format
29678569bb4SMatthew G. Knepley 
29778569bb4SMatthew G. Knepley   Collective on dm
29878569bb4SMatthew G. Knepley 
29978569bb4SMatthew G. Knepley   Input Parameters:
30078569bb4SMatthew G. Knepley + dm  - The dm to be written
30178569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
30278569bb4SMatthew G. Knepley - degree - the degree of the interpolation space
30378569bb4SMatthew G. Knepley 
30478569bb4SMatthew G. Knepley   Notes:
30578569bb4SMatthew G. Knepley   Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels)
30678569bb4SMatthew G. Knepley   consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted.
30778569bb4SMatthew G. Knepley 
30878569bb4SMatthew G. Knepley   If the dm has been distributed and exoid points to different files on each MPI rank, only the "local" part of the DM
30978569bb4SMatthew G. Knepley   (including "ghost" cells and vertices) will be written. If all exoid points to the same file, the resulting file will
31078569bb4SMatthew G. Knepley   probably be corrupted.
31178569bb4SMatthew G. Knepley 
31278569bb4SMatthew G. Knepley   DMPlex only represents geometry while most post-processing software expect that a mesh also provides information
31378569bb4SMatthew G. Knepley   on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2.
31478569bb4SMatthew G. Knepley   It should be extended to use PetscFE objects.
31578569bb4SMatthew G. Knepley 
31678569bb4SMatthew G. Knepley   This function will only handle TRI, TET, QUAD and HEX cells.
31778569bb4SMatthew G. Knepley   Level: beginner
31878569bb4SMatthew G. Knepley 
31978569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadNodal_PlexEXO(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO()
32078569bb4SMatthew G. Knepley */
32178569bb4SMatthew G. Knepley PetscErrorCode DMPlexView_ExodusII_Internal(DM dm, int exoid, PetscInt degree)
32278569bb4SMatthew G. Knepley {
32378569bb4SMatthew G. Knepley   enum ElemType {TRI, QUAD, TET, HEX};
32478569bb4SMatthew G. Knepley   MPI_Comm        comm;
32578569bb4SMatthew G. Knepley   /* Connectivity Variables */
32678569bb4SMatthew G. Knepley   PetscInt        cellsNotInConnectivity;
32778569bb4SMatthew G. Knepley   /* Cell Sets */
32878569bb4SMatthew G. Knepley   DMLabel         csLabel;
32978569bb4SMatthew G. Knepley   IS              csIS;
33078569bb4SMatthew G. Knepley   const PetscInt *csIdx;
33178569bb4SMatthew G. Knepley   PetscInt        num_cs, cs;
33278569bb4SMatthew G. Knepley   enum ElemType  *type;
333fe209ef9SBlaise Bourdin   PetscBool       hasLabel;
33478569bb4SMatthew G. Knepley   /* Coordinate Variables */
33578569bb4SMatthew G. Knepley   DM              cdm;
33678569bb4SMatthew G. Knepley   PetscSection    section;
33778569bb4SMatthew G. Knepley   Vec             coord;
33878569bb4SMatthew G. Knepley   PetscInt      **nodes;
339eae8a387SMatthew G. Knepley   PetscInt        depth, d, dim, skipCells = 0;
34078569bb4SMatthew G. Knepley   PetscInt        pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes;
34178569bb4SMatthew G. Knepley   PetscInt        num_vs, num_fs;
34278569bb4SMatthew G. Knepley   PetscMPIInt     rank, size;
34378569bb4SMatthew G. Knepley   const char     *dmName;
344fe209ef9SBlaise Bourdin   PetscInt        nodesTriP1[4]  = {3,0,0,0};
345fe209ef9SBlaise Bourdin   PetscInt        nodesTriP2[4]  = {3,3,0,0};
346fe209ef9SBlaise Bourdin   PetscInt        nodesQuadP1[4] = {4,0,0,0};
347fe209ef9SBlaise Bourdin   PetscInt        nodesQuadP2[4] = {4,4,0,1};
348fe209ef9SBlaise Bourdin   PetscInt        nodesTetP1[4]  = {4,0,0,0};
349fe209ef9SBlaise Bourdin   PetscInt        nodesTetP2[4]  = {4,6,0,0};
350fe209ef9SBlaise Bourdin   PetscInt        nodesHexP1[4]  = {8,0,0,0};
351fe209ef9SBlaise Bourdin   PetscInt        nodesHexP2[4]  = {8,12,6,1};
35278569bb4SMatthew G. Knepley   PetscErrorCode  ierr;
35378569bb4SMatthew G. Knepley 
35478569bb4SMatthew G. Knepley   PetscFunctionBegin;
35578569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
35678569bb4SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
35778569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
35878569bb4SMatthew G. Knepley   /* --- Get DM info --- */
35978569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) dm, &dmName);CHKERRQ(ierr);
36078569bb4SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
36178569bb4SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
36278569bb4SMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
36378569bb4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
36478569bb4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
36578569bb4SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);
36678569bb4SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
36778569bb4SMatthew G. Knepley   numCells    = cEnd - cStart;
36878569bb4SMatthew G. Knepley   numEdges    = eEnd - eStart;
36978569bb4SMatthew G. Knepley   numVertices = vEnd - vStart;
37078569bb4SMatthew G. Knepley   if (depth == 3) {numFaces = fEnd - fStart;}
37178569bb4SMatthew G. Knepley   else            {numFaces = 0;}
37278569bb4SMatthew G. Knepley   ierr = DMGetLabelSize(dm, "Cell Sets", &num_cs);CHKERRQ(ierr);
37378569bb4SMatthew G. Knepley   ierr = DMGetLabelSize(dm, "Vertex Sets", &num_vs);CHKERRQ(ierr);
37478569bb4SMatthew G. Knepley   ierr = DMGetLabelSize(dm, "Face Sets", &num_fs);CHKERRQ(ierr);
37578569bb4SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr);
37678569bb4SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
37778569bb4SMatthew G. Knepley   if (num_cs > 0) {
37878569bb4SMatthew G. Knepley     ierr = DMGetLabel(dm, "Cell Sets", &csLabel);CHKERRQ(ierr);
37978569bb4SMatthew G. Knepley     ierr = DMLabelGetValueIS(csLabel, &csIS);CHKERRQ(ierr);
38078569bb4SMatthew G. Knepley     ierr = ISGetIndices(csIS, &csIdx);CHKERRQ(ierr);
38178569bb4SMatthew G. Knepley   }
38278569bb4SMatthew G. Knepley   ierr = PetscMalloc1(num_cs, &nodes);CHKERRQ(ierr);
38378569bb4SMatthew G. Knepley   /* Set element type for each block and compute total number of nodes */
38478569bb4SMatthew G. Knepley   ierr = PetscMalloc1(num_cs, &type);CHKERRQ(ierr);
38578569bb4SMatthew G. Knepley   numNodes = numVertices;
38678569bb4SMatthew G. Knepley   if (degree == 2) {numNodes += numEdges;}
38778569bb4SMatthew G. Knepley   cellsNotInConnectivity = numCells;
38878569bb4SMatthew G. Knepley   for (cs = 0; cs < num_cs; ++cs) {
38978569bb4SMatthew G. Knepley     IS              stratumIS;
39078569bb4SMatthew G. Knepley     const PetscInt *cells;
39178569bb4SMatthew G. Knepley     PetscScalar    *xyz = NULL;
39278569bb4SMatthew G. Knepley     PetscInt        csSize, closureSize;
39378569bb4SMatthew G. Knepley 
39478569bb4SMatthew G. Knepley     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
39578569bb4SMatthew G. Knepley     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
39678569bb4SMatthew G. Knepley     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
39778569bb4SMatthew G. Knepley     ierr = DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr);
39878569bb4SMatthew G. Knepley     switch (dim) {
39978569bb4SMatthew G. Knepley       case 2:
40078569bb4SMatthew G. Knepley         if      (closureSize == 3*dim) {type[cs] = TRI;}
40178569bb4SMatthew G. Knepley         else if (closureSize == 4*dim) {type[cs] = QUAD;}
40278569bb4SMatthew G. Knepley         else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim);
40378569bb4SMatthew G. Knepley         break;
40478569bb4SMatthew G. Knepley       case 3:
40578569bb4SMatthew G. Knepley         if      (closureSize == 4*dim) {type[cs] = TET;}
40678569bb4SMatthew G. Knepley         else if (closureSize == 8*dim) {type[cs] = HEX;}
40778569bb4SMatthew G. Knepley         else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim);
40878569bb4SMatthew G. Knepley         break;
40978569bb4SMatthew G. Knepley       default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not handled by ExodusII viewer", dim);
41078569bb4SMatthew G. Knepley     }
41178569bb4SMatthew G. Knepley     if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;}
41278569bb4SMatthew G. Knepley     if ((degree == 2) && (type[cs] == HEX))  {numNodes += csSize; numNodes += numFaces;}
41378569bb4SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr);
41478569bb4SMatthew G. Knepley     /* Set nodes and Element type */
41578569bb4SMatthew G. Knepley     if (type[cs] == TRI) {
41678569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesTriP1;
41778569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesTriP2;
41878569bb4SMatthew G. Knepley     } else if (type[cs] == QUAD) {
41978569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesQuadP1;
42078569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesQuadP2;
42178569bb4SMatthew G. Knepley     } else if (type[cs] == TET) {
42278569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesTetP1;
42378569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesTetP2;
42478569bb4SMatthew G. Knepley     } else if (type[cs] == HEX) {
42578569bb4SMatthew G. Knepley       if      (degree == 1) nodes[cs] = nodesHexP1;
42678569bb4SMatthew G. Knepley       else if (degree == 2) nodes[cs] = nodesHexP2;
42778569bb4SMatthew G. Knepley     }
42878569bb4SMatthew G. Knepley     /* Compute the number of cells not in the connectivity table */
42978569bb4SMatthew G. Knepley     cellsNotInConnectivity -= nodes[cs][3]*csSize;
43078569bb4SMatthew G. Knepley 
43178569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
43278569bb4SMatthew G. Knepley     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
43378569bb4SMatthew G. Knepley   }
4346de834b4SMatthew G. Knepley   if (num_cs > 0) {PetscStackCallStandard(ex_put_init,(exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs));}
43578569bb4SMatthew G. Knepley   /* --- Connectivity --- */
43678569bb4SMatthew G. Knepley   for (cs = 0; cs < num_cs; ++cs) {
43778569bb4SMatthew G. Knepley     IS              stratumIS;
43878569bb4SMatthew G. Knepley     const PetscInt *cells;
439fe209ef9SBlaise Bourdin     PetscInt       *connect, off = 0;
440eae8a387SMatthew G. Knepley     PetscInt        edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0;
44178569bb4SMatthew G. Knepley     PetscInt        csSize, c, connectSize, closureSize;
44278569bb4SMatthew G. Knepley     char           *elem_type = NULL;
44378569bb4SMatthew G. Knepley     char            elem_type_tri3[]  = "TRI3",  elem_type_quad4[] = "QUAD4";
44478569bb4SMatthew G. Knepley     char            elem_type_tri6[]  = "TRI6",  elem_type_quad9[] = "QUAD9";
44578569bb4SMatthew G. Knepley     char            elem_type_tet4[]  = "TET4",  elem_type_hex8[]  = "HEX8";
44678569bb4SMatthew G. Knepley     char            elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27";
44778569bb4SMatthew G. Knepley 
44878569bb4SMatthew G. Knepley     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
44978569bb4SMatthew G. Knepley     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
45078569bb4SMatthew G. Knepley     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
45178569bb4SMatthew G. Knepley     /* Set Element type */
45278569bb4SMatthew G. Knepley     if (type[cs] == TRI) {
45378569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_tri3;
45478569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_tri6;
45578569bb4SMatthew G. Knepley     } else if (type[cs] == QUAD) {
45678569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_quad4;
45778569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_quad9;
45878569bb4SMatthew G. Knepley     } else if (type[cs] == TET) {
45978569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_tet4;
46078569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_tet10;
46178569bb4SMatthew G. Knepley     } else if (type[cs] == HEX) {
46278569bb4SMatthew G. Knepley       if      (degree == 1) elem_type = elem_type_hex8;
46378569bb4SMatthew G. Knepley       else if (degree == 2) elem_type = elem_type_hex27;
46478569bb4SMatthew G. Knepley     }
46578569bb4SMatthew G. Knepley     connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3];
466fe209ef9SBlaise Bourdin     ierr = PetscMalloc1(PetscMax(27,connectSize)*csSize, &connect);CHKERRQ(ierr);
467fe209ef9SBlaise Bourdin     PetscStackCallStandard(ex_put_block,(exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1));
46878569bb4SMatthew G. Knepley     /* Find number of vertices, edges, and faces in the closure */
46978569bb4SMatthew G. Knepley     verticesInClosure = nodes[cs][0];
47078569bb4SMatthew G. Knepley     if (depth > 1) {
47178569bb4SMatthew G. Knepley       if (dim == 2) {
47278569bb4SMatthew G. Knepley         ierr = DMPlexGetConeSize(dm, cells[0], &edgesInClosure);CHKERRQ(ierr);
47378569bb4SMatthew G. Knepley       } else if (dim == 3) {
47478569bb4SMatthew G. Knepley         PetscInt *closure = NULL;
47578569bb4SMatthew G. Knepley 
47678569bb4SMatthew G. Knepley         ierr = DMPlexGetConeSize(dm, cells[0], &facesInClosure);CHKERRQ(ierr);
47778569bb4SMatthew G. Knepley         ierr = DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
47878569bb4SMatthew G. Knepley         edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure;
47978569bb4SMatthew G. Knepley         ierr = DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
48078569bb4SMatthew G. Knepley       }
48178569bb4SMatthew G. Knepley     }
48278569bb4SMatthew G. Knepley     /* Get connectivity for each cell */
48378569bb4SMatthew G. Knepley     for (c = 0; c < csSize; ++c) {
48478569bb4SMatthew G. Knepley       PetscInt *closure = NULL;
48578569bb4SMatthew G. Knepley       PetscInt  temp, i;
48678569bb4SMatthew G. Knepley 
48778569bb4SMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
48878569bb4SMatthew G. Knepley       for (i = 0; i < connectSize; ++i) {
48978569bb4SMatthew G. Knepley         if (i < nodes[cs][0]) {/* Vertices */
490fe209ef9SBlaise Bourdin           connect[i+off] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1;
491fe209ef9SBlaise Bourdin           connect[i+off] -= cellsNotInConnectivity;
49278569bb4SMatthew G. Knepley         } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */
493fe209ef9SBlaise Bourdin           connect[i+off] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1;
494fe209ef9SBlaise Bourdin           if (nodes[cs][2] == 0) connect[i+off] -= numFaces;
495fe209ef9SBlaise Bourdin           connect[i+off] -= cellsNotInConnectivity;
49678569bb4SMatthew G. Knepley         } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]){ /* Cells */
497fe209ef9SBlaise Bourdin           connect[i+off] = closure[0] + 1;
498fe209ef9SBlaise Bourdin           connect[i+off] -= skipCells;
49978569bb4SMatthew G. Knepley         } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]){ /* Faces */
500fe209ef9SBlaise Bourdin           connect[i+off] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1;
501fe209ef9SBlaise Bourdin           connect[i+off] -= cellsNotInConnectivity;
50278569bb4SMatthew G. Knepley         } else {
503fe209ef9SBlaise Bourdin           connect[i+off] = -1;
50478569bb4SMatthew G. Knepley         }
50578569bb4SMatthew G. Knepley       }
50678569bb4SMatthew G. Knepley       /* Tetrahedra are inverted */
50778569bb4SMatthew G. Knepley       if (type[cs] == TET) {
508fe209ef9SBlaise Bourdin         temp = connect[0+off]; connect[0+off] = connect[1+off]; connect[1+off] = temp;
50978569bb4SMatthew G. Knepley         if (degree == 2) {
510e2c9586dSBlaise Bourdin           temp = connect[5+off]; connect[5+off] = connect[6+off]; connect[6+off] = temp;
511e2c9586dSBlaise Bourdin           temp = connect[7+off]; connect[7+off] = connect[8+off]; connect[8+off] = temp;
51278569bb4SMatthew G. Knepley         }
51378569bb4SMatthew G. Knepley       }
51478569bb4SMatthew G. Knepley       /* Hexahedra are inverted */
51578569bb4SMatthew G. Knepley       if (type[cs] == HEX) {
516fe209ef9SBlaise Bourdin         temp = connect[1+off]; connect[1+off] = connect[3+off]; connect[3+off] = temp;
51778569bb4SMatthew G. Knepley         if (degree == 2) {
518fe209ef9SBlaise Bourdin           temp = connect[8+off];  connect[8+off]  = connect[11+off]; connect[11+off] = temp;
519fe209ef9SBlaise Bourdin           temp = connect[9+off];  connect[9+off]  = connect[10+off]; connect[10+off] = temp;
520fe209ef9SBlaise Bourdin           temp = connect[16+off]; connect[16+off] = connect[17+off]; connect[17+off] = temp;
521fe209ef9SBlaise Bourdin           temp = connect[18+off]; connect[18+off] = connect[19+off]; connect[19+off] = temp;
52278569bb4SMatthew G. Knepley 
523fe209ef9SBlaise Bourdin           temp = connect[12+off]; connect[12+off] = connect[16+off]; connect[16+off] = temp;
524fe209ef9SBlaise Bourdin           temp = connect[13+off]; connect[13+off] = connect[17+off]; connect[17+off] = temp;
525fe209ef9SBlaise Bourdin           temp = connect[14+off]; connect[14+off] = connect[18+off]; connect[18+off] = temp;
526fe209ef9SBlaise Bourdin           temp = connect[15+off]; connect[15+off] = connect[19+off]; connect[19+off] = temp;
52778569bb4SMatthew G. Knepley 
528fe209ef9SBlaise Bourdin           temp = connect[23+off]; connect[23+off] = connect[26+off]; connect[26+off] = temp;
529fe209ef9SBlaise Bourdin           temp = connect[24+off]; connect[24+off] = connect[25+off]; connect[25+off] = temp;
530fe209ef9SBlaise Bourdin           temp = connect[25+off]; connect[25+off] = connect[26+off]; connect[26+off] = temp;
53178569bb4SMatthew G. Knepley         }
53278569bb4SMatthew G. Knepley       }
533fe209ef9SBlaise Bourdin       off += connectSize;
53478569bb4SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
53578569bb4SMatthew G. Knepley     }
536fe209ef9SBlaise Bourdin     PetscStackCallStandard(ex_put_conn,(exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0));
53778569bb4SMatthew G. Knepley     skipCells += (nodes[cs][3] == 0)*csSize;
53878569bb4SMatthew G. Knepley     ierr = PetscFree(connect);CHKERRQ(ierr);
53978569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
54078569bb4SMatthew G. Knepley     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
54178569bb4SMatthew G. Knepley   }
54278569bb4SMatthew G. Knepley   ierr = PetscFree(type);CHKERRQ(ierr);
54378569bb4SMatthew G. Knepley   /* --- Coordinates --- */
54478569bb4SMatthew G. Knepley   ierr = PetscSectionCreate(comm, &section);CHKERRQ(ierr);
54578569bb4SMatthew G. Knepley   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
546*1e50132fSMatthew G. Knepley   if (num_cs) {
54778569bb4SMatthew G. Knepley     for (d = 0; d < depth; ++d) {
54878569bb4SMatthew G. Knepley       ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
54978569bb4SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
55078569bb4SMatthew G. Knepley         ierr = PetscSectionSetDof(section, p, nodes[0][d] > 0);CHKERRQ(ierr);
55178569bb4SMatthew G. Knepley       }
55278569bb4SMatthew G. Knepley     }
553*1e50132fSMatthew G. Knepley   }
55478569bb4SMatthew G. Knepley   for (cs = 0; cs < num_cs; ++cs) {
55578569bb4SMatthew G. Knepley     IS              stratumIS;
55678569bb4SMatthew G. Knepley     const PetscInt *cells;
55778569bb4SMatthew G. Knepley     PetscInt        csSize, c;
55878569bb4SMatthew G. Knepley 
55978569bb4SMatthew G. Knepley     ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
56078569bb4SMatthew G. Knepley     ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
56178569bb4SMatthew G. Knepley     ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
56278569bb4SMatthew G. Knepley     for (c = 0; c < csSize; ++c) {
56378569bb4SMatthew G. Knepley       ierr = PetscSectionSetDof(section, cells[c], nodes[cs][3] > 0);CHKERRQ(ierr);
56478569bb4SMatthew G. Knepley     }
56578569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
56678569bb4SMatthew G. Knepley     ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
56778569bb4SMatthew G. Knepley   }
56878569bb4SMatthew G. Knepley   if (num_cs > 0) {
56978569bb4SMatthew G. Knepley     ierr = ISRestoreIndices(csIS, &csIdx);CHKERRQ(ierr);
57078569bb4SMatthew G. Knepley     ierr = ISDestroy(&csIS);CHKERRQ(ierr);
57178569bb4SMatthew G. Knepley   }
57278569bb4SMatthew G. Knepley   ierr = PetscFree(nodes);CHKERRQ(ierr);
57378569bb4SMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
57478569bb4SMatthew G. Knepley   if (numNodes > 0) {
57578569bb4SMatthew G. Knepley     const char  *coordNames[3] = {"x", "y", "z"};
57678569bb4SMatthew G. Knepley     PetscScalar *coords, *closure;
57778569bb4SMatthew G. Knepley     PetscReal   *cval;
57878569bb4SMatthew G. Knepley     PetscInt     hasDof, n = 0;
57978569bb4SMatthew G. Knepley 
58078569bb4SMatthew G. Knepley     /* There can't be more than 24 values in the closure of a point for the coord section */
5816de834b4SMatthew G. Knepley     ierr = PetscCalloc3(numNodes*3, &coords, dim, &cval, 24, &closure);CHKERRQ(ierr);
58278569bb4SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr);
58378569bb4SMatthew G. Knepley     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
58478569bb4SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
58578569bb4SMatthew G. Knepley       ierr = PetscSectionGetDof(section, p, &hasDof);
58678569bb4SMatthew G. Knepley       if (hasDof) {
58778569bb4SMatthew G. Knepley         PetscInt closureSize = 24, j;
58878569bb4SMatthew G. Knepley 
58978569bb4SMatthew G. Knepley         ierr = DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure);CHKERRQ(ierr);
59078569bb4SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
59178569bb4SMatthew G. Knepley           cval[d] = 0.0;
59278569bb4SMatthew G. Knepley           for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d];
59378569bb4SMatthew G. Knepley           coords[d*numNodes+n] = cval[d] * dim / closureSize;
59478569bb4SMatthew G. Knepley         }
59578569bb4SMatthew G. Knepley         ++n;
59678569bb4SMatthew G. Knepley       }
59778569bb4SMatthew G. Knepley     }
5986de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_put_coord,(exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes]));
59978569bb4SMatthew G. Knepley     ierr = PetscFree3(coords, cval, closure);CHKERRQ(ierr);
6006de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_put_coord_names,(exoid, (char **) coordNames));
60178569bb4SMatthew G. Knepley   }
60278569bb4SMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
603fe209ef9SBlaise Bourdin   /* --- Node Sets/Vertex Sets --- */
604fe209ef9SBlaise Bourdin   DMHasLabel(dm, "Vertex Sets", &hasLabel);
605fe209ef9SBlaise Bourdin   if (hasLabel) {
606fe209ef9SBlaise Bourdin     PetscInt        i, vs, vsSize;
607fe209ef9SBlaise Bourdin     const PetscInt *vsIdx, *vertices;
608fe209ef9SBlaise Bourdin     PetscInt       *nodeList;
609fe209ef9SBlaise Bourdin     IS              vsIS, stratumIS;
610fe209ef9SBlaise Bourdin     DMLabel         vsLabel;
611fe209ef9SBlaise Bourdin     ierr = DMGetLabel(dm, "Vertex Sets", &vsLabel);CHKERRQ(ierr);
612fe209ef9SBlaise Bourdin     ierr = DMLabelGetValueIS(vsLabel, &vsIS);CHKERRQ(ierr);
613fe209ef9SBlaise Bourdin     ierr = ISGetIndices(vsIS, &vsIdx);CHKERRQ(ierr);
614fe209ef9SBlaise Bourdin     for (vs=0; vs<num_vs; ++vs) {
615fe209ef9SBlaise Bourdin       ierr = DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS);CHKERRQ(ierr);
616fe209ef9SBlaise Bourdin       ierr = ISGetIndices(stratumIS, &vertices);CHKERRQ(ierr);
617fe209ef9SBlaise Bourdin       ierr = ISGetSize(stratumIS, &vsSize);CHKERRQ(ierr);
618fe209ef9SBlaise Bourdin       ierr = PetscMalloc1(vsSize, &nodeList);
619fe209ef9SBlaise Bourdin       for (i=0; i<vsSize; ++i) {
620fe209ef9SBlaise Bourdin         nodeList[i] = vertices[i] - skipCells + 1;
621fe209ef9SBlaise Bourdin       }
622fe209ef9SBlaise Bourdin       PetscStackCallStandard(ex_put_set_param,(exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0));
623fe209ef9SBlaise Bourdin       PetscStackCallStandard(ex_put_set,(exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL));
624fe209ef9SBlaise Bourdin       ierr = ISRestoreIndices(stratumIS, &vertices);CHKERRQ(ierr);
625fe209ef9SBlaise Bourdin       ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
626fe209ef9SBlaise Bourdin       ierr = PetscFree(nodeList);CHKERRQ(ierr);
627fe209ef9SBlaise Bourdin     }
628fe209ef9SBlaise Bourdin     ierr = ISRestoreIndices(vsIS, &vsIdx);CHKERRQ(ierr);
629fe209ef9SBlaise Bourdin     ierr = ISDestroy(&vsIS);CHKERRQ(ierr);
630fe209ef9SBlaise Bourdin   }
631fe209ef9SBlaise Bourdin   /* --- Side Sets/Face Sets --- */
632fe209ef9SBlaise Bourdin   ierr = DMHasLabel(dm, "Face Sets", &hasLabel);CHKERRQ(ierr);
633fe209ef9SBlaise Bourdin   if (hasLabel) {
634fe209ef9SBlaise Bourdin     PetscInt        i, j, fs, fsSize;
635fe209ef9SBlaise Bourdin     const PetscInt *fsIdx, *faces;
636fe209ef9SBlaise Bourdin     IS              fsIS, stratumIS;
637fe209ef9SBlaise Bourdin     DMLabel         fsLabel;
638fe209ef9SBlaise Bourdin     PetscInt        numPoints, *points;
639fe209ef9SBlaise Bourdin     PetscInt        elem_list_size = 0;
640fe209ef9SBlaise Bourdin     PetscInt       *elem_list, *elem_ind, *side_list;
641fe209ef9SBlaise Bourdin 
642fe209ef9SBlaise Bourdin     ierr = DMGetLabel(dm, "Face Sets", &fsLabel);CHKERRQ(ierr);
643fe209ef9SBlaise Bourdin     /* Compute size of Node List and Element List */
644fe209ef9SBlaise Bourdin     ierr = DMLabelGetValueIS(fsLabel, &fsIS);CHKERRQ(ierr);
645fe209ef9SBlaise Bourdin     ierr = ISGetIndices(fsIS, &fsIdx);CHKERRQ(ierr);
646fe209ef9SBlaise Bourdin     for (fs=0; fs<num_fs; ++fs) {
647feb237baSPierre Jolivet       ierr = DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS);CHKERRQ(ierr);
648fe209ef9SBlaise Bourdin       ierr = ISGetSize(stratumIS, &fsSize);CHKERRQ(ierr);
649fe209ef9SBlaise Bourdin       elem_list_size += fsSize;
650fe209ef9SBlaise Bourdin       ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
651fe209ef9SBlaise Bourdin     }
652fe209ef9SBlaise Bourdin     ierr = PetscMalloc3(num_fs, &elem_ind,
653fe209ef9SBlaise Bourdin                         elem_list_size, &elem_list,
654fe209ef9SBlaise Bourdin                         elem_list_size, &side_list);CHKERRQ(ierr);
655fe209ef9SBlaise Bourdin     elem_ind[0] = 0;
656fe209ef9SBlaise Bourdin     for (fs=0; fs<num_fs; ++fs) {
657fe209ef9SBlaise Bourdin       ierr = DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS);CHKERRQ(ierr);
658fe209ef9SBlaise Bourdin       ierr = ISGetIndices(stratumIS, &faces);CHKERRQ(ierr);
659fe209ef9SBlaise Bourdin       ierr = ISGetSize(stratumIS, &fsSize);CHKERRQ(ierr);
660fe209ef9SBlaise Bourdin       /* Set Parameters */
661fe209ef9SBlaise Bourdin       PetscStackCallStandard(ex_put_set_param,(exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0));
662fe209ef9SBlaise Bourdin       /* Indices */
663fe209ef9SBlaise Bourdin       if (fs<num_fs-1) {
664fe209ef9SBlaise Bourdin         elem_ind[fs+1] = elem_ind[fs] + fsSize;
665fe209ef9SBlaise Bourdin       }
666fe209ef9SBlaise Bourdin 
667fe209ef9SBlaise Bourdin       for (i=0; i<fsSize; ++i) {
668fe209ef9SBlaise Bourdin         /* Element List */
669fe209ef9SBlaise Bourdin         points = NULL;
670fe209ef9SBlaise Bourdin         ierr = DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points);CHKERRQ(ierr);
671fe209ef9SBlaise Bourdin         elem_list[elem_ind[fs] + i] = points[2] +1;
672fe209ef9SBlaise Bourdin         ierr = DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points);CHKERRQ(ierr);
673fe209ef9SBlaise Bourdin 
674fe209ef9SBlaise Bourdin         /* Side List */
675fe209ef9SBlaise Bourdin         points = NULL;
676fe209ef9SBlaise Bourdin         ierr = DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
677fe209ef9SBlaise Bourdin         for (j=1; j<numPoints; ++j) {
678fe209ef9SBlaise Bourdin           if (points[j*2]==faces[i]) {break;}
679fe209ef9SBlaise Bourdin         }
680fe209ef9SBlaise Bourdin         /* Convert HEX sides */
681fe209ef9SBlaise Bourdin         if (numPoints == 27) {
682fe209ef9SBlaise Bourdin           if      (j == 1) {j = 5;}
683fe209ef9SBlaise Bourdin           else if (j == 2) {j = 6;}
684fe209ef9SBlaise Bourdin           else if (j == 3) {j = 1;}
685fe209ef9SBlaise Bourdin           else if (j == 4) {j = 3;}
686fe209ef9SBlaise Bourdin           else if (j == 5) {j = 2;}
687fe209ef9SBlaise Bourdin           else if (j == 6) {j = 4;}
688fe209ef9SBlaise Bourdin         }
689fe209ef9SBlaise Bourdin         /* Convert TET sides */
690fe209ef9SBlaise Bourdin         if (numPoints == 15) {
691fe209ef9SBlaise Bourdin           --j;
692fe209ef9SBlaise Bourdin           if (j == 0) {j = 4;}
693fe209ef9SBlaise Bourdin         }
694fe209ef9SBlaise Bourdin         side_list[elem_ind[fs] + i] = j;
695fe209ef9SBlaise Bourdin         ierr = DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
696fe209ef9SBlaise Bourdin 
697fe209ef9SBlaise Bourdin       }
698fe209ef9SBlaise Bourdin       ierr = ISRestoreIndices(stratumIS, &faces);CHKERRQ(ierr);
699fe209ef9SBlaise Bourdin       ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
700fe209ef9SBlaise Bourdin     }
701fe209ef9SBlaise Bourdin     ierr = ISRestoreIndices(fsIS, &fsIdx);CHKERRQ(ierr);
702fe209ef9SBlaise Bourdin     ierr = ISDestroy(&fsIS);CHKERRQ(ierr);
703fe209ef9SBlaise Bourdin 
704fe209ef9SBlaise Bourdin     /* Put side sets */
705fe209ef9SBlaise Bourdin     for (fs=0; fs<num_fs; ++fs) {
706fe209ef9SBlaise Bourdin       PetscStackCallStandard(ex_put_set,(exoid, EX_SIDE_SET, fsIdx[fs], &elem_list[elem_ind[fs]], &side_list[elem_ind[fs]]));
707fe209ef9SBlaise Bourdin     }
708fe209ef9SBlaise Bourdin     ierr = PetscFree3(elem_ind,elem_list,side_list);CHKERRQ(ierr);
709fe209ef9SBlaise Bourdin   }
71078569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
71178569bb4SMatthew G. Knepley }
71278569bb4SMatthew G. Knepley 
71378569bb4SMatthew G. Knepley /*
71478569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file
71578569bb4SMatthew G. Knepley 
71678569bb4SMatthew G. Knepley   Collective on v
71778569bb4SMatthew G. Knepley 
71878569bb4SMatthew G. Knepley   Input Parameters:
71978569bb4SMatthew G. Knepley + v  - The vector to be written
72078569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
72178569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1)
72278569bb4SMatthew G. Knepley 
72378569bb4SMatthew G. Knepley   Notes:
72478569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
72578569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
72678569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
72778569bb4SMatthew G. Knepley   amongst all nodal variables.
72878569bb4SMatthew G. Knepley 
72978569bb4SMatthew G. Knepley   Level: beginner
73078569bb4SMatthew G. Knepley 
73178569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecLoadNodal_PlexEXO(),VecViewZonal_PlexEXO(),VecLoadZonal_PlexEXO()
73278569bb4SMatthew G. Knepley @*/
73378569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step)
73478569bb4SMatthew G. Knepley {
73578569bb4SMatthew G. Knepley   MPI_Comm         comm;
73678569bb4SMatthew G. Knepley   PetscMPIInt      size;
73778569bb4SMatthew G. Knepley   DM               dm;
73878569bb4SMatthew G. Knepley   Vec              vNatural, vComp;
73922a7544eSVaclav Hapla   const PetscScalar *varray;
74078569bb4SMatthew G. Knepley   const char      *vecname;
74178569bb4SMatthew G. Knepley   PetscInt         xs, xe, bs;
74278569bb4SMatthew G. Knepley   PetscBool        useNatural;
74378569bb4SMatthew G. Knepley   int              offset;
74478569bb4SMatthew G. Knepley   PetscErrorCode   ierr;
74578569bb4SMatthew G. Knepley 
74678569bb4SMatthew G. Knepley   PetscFunctionBegin;
74778569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
74878569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
74978569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
75078569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
75178569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
75278569bb4SMatthew G. Knepley   if (useNatural) {
75378569bb4SMatthew G. Knepley     ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr);
75478569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr);
75578569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr);
75678569bb4SMatthew G. Knepley   } else {
75778569bb4SMatthew G. Knepley     vNatural = v;
75878569bb4SMatthew G. Knepley   }
75978569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
76078569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
76178569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr);
76278569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file\n", vecname);
76378569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
76478569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
76578569bb4SMatthew G. Knepley      We assume that they are stored sequentially */
76678569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
76778569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
76878569bb4SMatthew G. Knepley   if (bs == 1) {
76978569bb4SMatthew G. Knepley     ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr);
7706de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray));
77178569bb4SMatthew G. Knepley     ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr);
77278569bb4SMatthew G. Knepley   } else {
77378569bb4SMatthew G. Knepley     IS       compIS;
77478569bb4SMatthew G. Knepley     PetscInt c;
77578569bb4SMatthew G. Knepley 
77678569bb4SMatthew G. Knepley     ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);
77778569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
77878569bb4SMatthew G. Knepley       ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
77978569bb4SMatthew G. Knepley       ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
78078569bb4SMatthew G. Knepley       ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr);
7816de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray));
78278569bb4SMatthew G. Knepley       ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr);
78378569bb4SMatthew G. Knepley       ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
78478569bb4SMatthew G. Knepley     }
78578569bb4SMatthew G. Knepley     ierr = ISDestroy(&compIS);CHKERRQ(ierr);
78678569bb4SMatthew G. Knepley   }
78778569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);}
78878569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
78978569bb4SMatthew G. Knepley }
79078569bb4SMatthew G. Knepley 
79178569bb4SMatthew G. Knepley /*
79278569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file
79378569bb4SMatthew G. Knepley 
79478569bb4SMatthew G. Knepley   Collective on v
79578569bb4SMatthew G. Knepley 
79678569bb4SMatthew G. Knepley   Input Parameters:
79778569bb4SMatthew G. Knepley + v  - The vector to be written
79878569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
79978569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1)
80078569bb4SMatthew G. Knepley 
80178569bb4SMatthew G. Knepley   Notes:
80278569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
80378569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
80478569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
80578569bb4SMatthew G. Knepley   amongst all nodal variables.
80678569bb4SMatthew G. Knepley 
80778569bb4SMatthew G. Knepley   Level: beginner
80878569bb4SMatthew G. Knepley 
80978569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO()
81078569bb4SMatthew G. Knepley */
81178569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step)
81278569bb4SMatthew G. Knepley {
81378569bb4SMatthew G. Knepley   MPI_Comm       comm;
81478569bb4SMatthew G. Knepley   PetscMPIInt    size;
81578569bb4SMatthew G. Knepley   DM             dm;
81678569bb4SMatthew G. Knepley   Vec            vNatural, vComp;
81722a7544eSVaclav Hapla   PetscScalar   *varray;
81878569bb4SMatthew G. Knepley   const char    *vecname;
81978569bb4SMatthew G. Knepley   PetscInt       xs, xe, bs;
82078569bb4SMatthew G. Knepley   PetscBool      useNatural;
82178569bb4SMatthew G. Knepley   int            offset;
82278569bb4SMatthew G. Knepley   PetscErrorCode ierr;
82378569bb4SMatthew G. Knepley 
82478569bb4SMatthew G. Knepley   PetscFunctionBegin;
82578569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
82678569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
82778569bb4SMatthew G. Knepley   ierr = VecGetDM(v,&dm);CHKERRQ(ierr);
82878569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
82978569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
83078569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
83178569bb4SMatthew G. Knepley   else            {vNatural = v;}
83278569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
83378569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
83478569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr);
83578569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file", vecname);
83678569bb4SMatthew G. Knepley   /* Read local chunk from the file */
83778569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
83878569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
83978569bb4SMatthew G. Knepley   if (bs == 1) {
84078569bb4SMatthew G. Knepley     ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr);
8416de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray));
84278569bb4SMatthew G. Knepley     ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr);
84378569bb4SMatthew G. Knepley   } else {
84478569bb4SMatthew G. Knepley     IS       compIS;
84578569bb4SMatthew G. Knepley     PetscInt c;
84678569bb4SMatthew G. Knepley 
84778569bb4SMatthew G. Knepley     ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);
84878569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
84978569bb4SMatthew G. Knepley       ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
85078569bb4SMatthew G. Knepley       ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
85178569bb4SMatthew G. Knepley       ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr);
8526de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray));
85378569bb4SMatthew G. Knepley       ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr);
85478569bb4SMatthew G. Knepley       ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
85578569bb4SMatthew G. Knepley     }
85678569bb4SMatthew G. Knepley     ierr = ISDestroy(&compIS);CHKERRQ(ierr);
85778569bb4SMatthew G. Knepley   }
85878569bb4SMatthew G. Knepley   if (useNatural) {
85978569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr);
86078569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr);
86178569bb4SMatthew G. Knepley     ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);
86278569bb4SMatthew G. Knepley   }
86378569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
86478569bb4SMatthew G. Knepley }
86578569bb4SMatthew G. Knepley 
86678569bb4SMatthew G. Knepley /*
86778569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file
86878569bb4SMatthew G. Knepley 
86978569bb4SMatthew G. Knepley   Collective on v
87078569bb4SMatthew G. Knepley 
87178569bb4SMatthew G. Knepley   Input Parameters:
87278569bb4SMatthew G. Knepley + v  - The vector to be written
87378569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
87478569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1)
87578569bb4SMatthew G. Knepley 
87678569bb4SMatthew G. Knepley   Notes:
87778569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
87878569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
87978569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
88078569bb4SMatthew G. Knepley   amongst all zonal variables.
88178569bb4SMatthew G. Knepley 
88278569bb4SMatthew G. Knepley   Level: beginner
88378569bb4SMatthew G. Knepley 
88478569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Zonal_Internal()
88578569bb4SMatthew G. Knepley */
88678569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step)
88778569bb4SMatthew G. Knepley {
88878569bb4SMatthew G. Knepley   MPI_Comm          comm;
88978569bb4SMatthew G. Knepley   PetscMPIInt       size;
89078569bb4SMatthew G. Knepley   DM                dm;
89178569bb4SMatthew G. Knepley   Vec               vNatural, vComp;
89222a7544eSVaclav Hapla   const PetscScalar *varray;
89378569bb4SMatthew G. Knepley   const char       *vecname;
89478569bb4SMatthew G. Knepley   PetscInt          xs, xe, bs;
89578569bb4SMatthew G. Knepley   PetscBool         useNatural;
89678569bb4SMatthew G. Knepley   int               offset;
89778569bb4SMatthew G. Knepley   IS                compIS;
89878569bb4SMatthew G. Knepley   PetscInt         *csSize, *csID;
89978569bb4SMatthew G. Knepley   PetscInt          numCS, set, csxs = 0;
90078569bb4SMatthew G. Knepley   PetscErrorCode    ierr;
90178569bb4SMatthew G. Knepley 
90278569bb4SMatthew G. Knepley   PetscFunctionBegin;
90378569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)v, &comm);CHKERRQ(ierr);
90478569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
90578569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
90678569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
90778569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
90878569bb4SMatthew G. Knepley   if (useNatural) {
90978569bb4SMatthew G. Knepley     ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr);
91078569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr);
91178569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr);
91278569bb4SMatthew G. Knepley   } else {
91378569bb4SMatthew G. Knepley     vNatural = v;
91478569bb4SMatthew G. Knepley   }
91578569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
91678569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
91778569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr);
91878569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname);
91978569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
92078569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
92178569bb4SMatthew G. Knepley      We assume that they are stored sequentially
92278569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
92378569bb4SMatthew G. Knepley      but once the vector has been reordered to natural size, we cannot use the label informations
92478569bb4SMatthew G. Knepley      to figure out what to save where. */
92578569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
92678569bb4SMatthew G. Knepley   ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr);
9276de834b4SMatthew G. Knepley   PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID));
92878569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
92978569bb4SMatthew G. Knepley     ex_block block;
93078569bb4SMatthew G. Knepley 
93178569bb4SMatthew G. Knepley     block.id   = csID[set];
93278569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
9336de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_block_param,(exoid, &block));
93478569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
93578569bb4SMatthew G. Knepley   }
93678569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
93778569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
93878569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);}
93978569bb4SMatthew G. Knepley   for (set = 0; set < numCS; set++) {
94078569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
94178569bb4SMatthew G. Knepley 
94278569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
94378569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
94478569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
94578569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
94678569bb4SMatthew G. Knepley     if (bs == 1) {
94778569bb4SMatthew G. Knepley       ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr);
9486de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]));
94978569bb4SMatthew G. Knepley       ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr);
95078569bb4SMatthew G. Knepley     } else {
95178569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
95278569bb4SMatthew G. Knepley         ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
95378569bb4SMatthew G. Knepley         ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
95478569bb4SMatthew G. Knepley         ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr);
9556de834b4SMatthew G. Knepley         PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]));
95678569bb4SMatthew G. Knepley         ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr);
95778569bb4SMatthew G. Knepley         ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
95878569bb4SMatthew G. Knepley       }
95978569bb4SMatthew G. Knepley     }
96078569bb4SMatthew G. Knepley     csxs += csSize[set];
96178569bb4SMatthew G. Knepley   }
96278569bb4SMatthew G. Knepley   ierr = PetscFree2(csID, csSize);CHKERRQ(ierr);
96378569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);}
96478569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMRestoreGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
96578569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
96678569bb4SMatthew G. Knepley }
96778569bb4SMatthew G. Knepley 
96878569bb4SMatthew G. Knepley /*
96978569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file
97078569bb4SMatthew G. Knepley 
97178569bb4SMatthew G. Knepley   Collective on v
97278569bb4SMatthew G. Knepley 
97378569bb4SMatthew G. Knepley   Input Parameters:
97478569bb4SMatthew G. Knepley + v  - The vector to be written
97578569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
97678569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1)
97778569bb4SMatthew G. Knepley 
97878569bb4SMatthew G. Knepley   Notes:
97978569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
98078569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
98178569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
98278569bb4SMatthew G. Knepley   amongst all zonal variables.
98378569bb4SMatthew G. Knepley 
98478569bb4SMatthew G. Knepley   Level: beginner
98578569bb4SMatthew G. Knepley 
98678569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Zonal_Internal()
98778569bb4SMatthew G. Knepley */
98878569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step)
98978569bb4SMatthew G. Knepley {
99078569bb4SMatthew G. Knepley   MPI_Comm          comm;
99178569bb4SMatthew G. Knepley   PetscMPIInt       size;
99278569bb4SMatthew G. Knepley   DM                dm;
99378569bb4SMatthew G. Knepley   Vec               vNatural, vComp;
99422a7544eSVaclav Hapla   PetscScalar      *varray;
99578569bb4SMatthew G. Knepley   const char       *vecname;
99678569bb4SMatthew G. Knepley   PetscInt          xs, xe, bs;
99778569bb4SMatthew G. Knepley   PetscBool         useNatural;
99878569bb4SMatthew G. Knepley   int               offset;
99978569bb4SMatthew G. Knepley   IS                compIS;
100078569bb4SMatthew G. Knepley   PetscInt         *csSize, *csID;
100178569bb4SMatthew G. Knepley   PetscInt          numCS, set, csxs = 0;
100278569bb4SMatthew G. Knepley   PetscErrorCode    ierr;
100378569bb4SMatthew G. Knepley 
100478569bb4SMatthew G. Knepley   PetscFunctionBegin;
100578569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
100678569bb4SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
100778569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
100878569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
100978569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
101078569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
101178569bb4SMatthew G. Knepley   else            {vNatural = v;}
101278569bb4SMatthew G. Knepley   /* Get the location of the variable in the exodus file */
101378569bb4SMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
101478569bb4SMatthew G. Knepley   ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr);
101578569bb4SMatthew G. Knepley   if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname);
101678569bb4SMatthew G. Knepley   /* Read local chunk of the result in the exodus file
101778569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
101878569bb4SMatthew G. Knepley      We assume that they are stored sequentially
101978569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
102078569bb4SMatthew G. Knepley      but once the vector has been reordered to natural size, we cannot use the label informations
102178569bb4SMatthew G. Knepley      to figure out what to save where. */
102278569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
102378569bb4SMatthew G. Knepley   ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr);
10246de834b4SMatthew G. Knepley   PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID));
102578569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
102678569bb4SMatthew G. Knepley     ex_block block;
102778569bb4SMatthew G. Knepley 
102878569bb4SMatthew G. Knepley     block.id   = csID[set];
102978569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
10306de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_block_param,(exoid, &block));
103178569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
103278569bb4SMatthew G. Knepley   }
103378569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
103478569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
103578569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);}
103678569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
103778569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
103878569bb4SMatthew G. Knepley 
103978569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
104078569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
104178569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
104278569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
104378569bb4SMatthew G. Knepley     if (bs == 1) {
104478569bb4SMatthew G. Knepley       ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr);
10456de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]));
104678569bb4SMatthew G. Knepley       ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr);
104778569bb4SMatthew G. Knepley     } else {
104878569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
104978569bb4SMatthew G. Knepley         ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
105078569bb4SMatthew G. Knepley         ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
105178569bb4SMatthew G. Knepley         ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr);
10526de834b4SMatthew G. Knepley         PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]));
105378569bb4SMatthew G. Knepley         ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr);
105478569bb4SMatthew G. Knepley         ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
105578569bb4SMatthew G. Knepley       }
105678569bb4SMatthew G. Knepley     }
105778569bb4SMatthew G. Knepley     csxs += csSize[set];
105878569bb4SMatthew G. Knepley   }
105978569bb4SMatthew G. Knepley   ierr = PetscFree2(csID, csSize);CHKERRQ(ierr);
106078569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);}
106178569bb4SMatthew G. Knepley   if (useNatural) {
106278569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr);
106378569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr);
106478569bb4SMatthew G. Knepley     ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);
106578569bb4SMatthew G. Knepley   }
106678569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
106778569bb4SMatthew G. Knepley }
1068b53c8456SSatish Balay #endif
106978569bb4SMatthew G. Knepley 
1070*1e50132fSMatthew G. Knepley /*@
1071*1e50132fSMatthew G. Knepley   PetscViewerExodusIIGetId - Get the file id of the ExodusII file
1072*1e50132fSMatthew G. Knepley 
1073*1e50132fSMatthew G. Knepley   Logically Collective on PetscViewer
1074*1e50132fSMatthew G. Knepley 
1075*1e50132fSMatthew G. Knepley   Input Parameter:
1076*1e50132fSMatthew G. Knepley .  viewer - the PetscViewer
1077*1e50132fSMatthew G. Knepley 
1078*1e50132fSMatthew G. Knepley   Output Parameter:
1079*1e50132fSMatthew G. Knepley -  exoid - The ExodusII file id
1080*1e50132fSMatthew G. Knepley 
1081*1e50132fSMatthew G. Knepley   Level: intermediate
1082*1e50132fSMatthew G. Knepley 
1083*1e50132fSMatthew G. Knepley .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
1084*1e50132fSMatthew G. Knepley @*/
1085*1e50132fSMatthew G. Knepley PetscErrorCode PetscViewerExodusIIGetId(PetscViewer viewer, int *exoid)
1086*1e50132fSMatthew G. Knepley {
1087*1e50132fSMatthew G. Knepley   PetscErrorCode ierr;
1088*1e50132fSMatthew G. Knepley 
1089*1e50132fSMatthew G. Knepley   PetscFunctionBegin;
1090*1e50132fSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1091*1e50132fSMatthew G. Knepley   ierr = PetscTryMethod(viewer, "PetscViewerExodusIIGetId_C",(PetscViewer,int*),(viewer,exoid));CHKERRQ(ierr);
1092*1e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1093*1e50132fSMatthew G. Knepley }
1094*1e50132fSMatthew G. Knepley 
109533751fbdSMichael Lange /*@C
1096eaf898f9SPatrick Sanan   DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file.
109733751fbdSMichael Lange 
1098d083f849SBarry Smith   Collective
109933751fbdSMichael Lange 
110033751fbdSMichael Lange   Input Parameters:
110133751fbdSMichael Lange + comm  - The MPI communicator
110233751fbdSMichael Lange . filename - The name of the ExodusII file
110333751fbdSMichael Lange - interpolate - Create faces and edges in the mesh
110433751fbdSMichael Lange 
110533751fbdSMichael Lange   Output Parameter:
110633751fbdSMichael Lange . dm  - The DM object representing the mesh
110733751fbdSMichael Lange 
110833751fbdSMichael Lange   Level: beginner
110933751fbdSMichael Lange 
111033751fbdSMichael Lange .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus()
111133751fbdSMichael Lange @*/
111233751fbdSMichael Lange PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
111333751fbdSMichael Lange {
111433751fbdSMichael Lange   PetscMPIInt    rank;
111533751fbdSMichael Lange   PetscErrorCode ierr;
111633751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
1117e84f7031SBlaise Bourdin   int   CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1;
111833751fbdSMichael Lange   float version;
111933751fbdSMichael Lange #endif
112033751fbdSMichael Lange 
112133751fbdSMichael Lange   PetscFunctionBegin;
112233751fbdSMichael Lange   PetscValidCharPointer(filename, 2);
112333751fbdSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
112433751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
112533751fbdSMichael Lange   if (!rank) {
112633751fbdSMichael Lange     exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version);
112733751fbdSMichael Lange     if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename);
112833751fbdSMichael Lange   }
112933751fbdSMichael Lange   ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr);
11306de834b4SMatthew G. Knepley   if (!rank) {PetscStackCallStandard(ex_close,(exoid));}
113133751fbdSMichael Lange #else
113233751fbdSMichael Lange   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
113333751fbdSMichael Lange #endif
113433751fbdSMichael Lange   PetscFunctionReturn(0);
113533751fbdSMichael Lange }
113633751fbdSMichael Lange 
1137552f7358SJed Brown /*@
113833751fbdSMichael Lange   DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID.
1139552f7358SJed Brown 
1140d083f849SBarry Smith   Collective
1141552f7358SJed Brown 
1142552f7358SJed Brown   Input Parameters:
1143552f7358SJed Brown + comm  - The MPI communicator
1144552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
1145552f7358SJed Brown - interpolate - Create faces and edges in the mesh
1146552f7358SJed Brown 
1147552f7358SJed Brown   Output Parameter:
1148552f7358SJed Brown . dm  - The DM object representing the mesh
1149552f7358SJed Brown 
1150552f7358SJed Brown   Level: beginner
1151552f7358SJed Brown 
1152e84e7769SJed Brown .seealso: DMPLEX, DMCreate()
1153552f7358SJed Brown @*/
1154552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm)
1155552f7358SJed Brown {
1156552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
1157552f7358SJed Brown   PetscMPIInt    num_proc, rank;
1158552f7358SJed Brown   PetscSection   coordSection;
1159552f7358SJed Brown   Vec            coordinates;
1160552f7358SJed Brown   PetscScalar    *coords;
1161552f7358SJed Brown   PetscInt       coordSize, v;
1162552f7358SJed Brown   PetscErrorCode ierr;
1163552f7358SJed Brown   /* Read from ex_get_init() */
1164552f7358SJed Brown   char title[PETSC_MAX_PATH_LEN+1];
1165e8f6893fSMatthew G. Knepley   int  dim    = 0, numVertices = 0, numCells = 0, numHybridCells = 0;
1166552f7358SJed Brown   int  num_cs = 0, num_vs = 0, num_fs = 0;
1167552f7358SJed Brown #endif
1168552f7358SJed Brown 
1169552f7358SJed Brown   PetscFunctionBegin;
1170552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
1171552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1172552f7358SJed Brown   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
1173552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
1174552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1175552f7358SJed Brown   /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */
1176552f7358SJed Brown   if (!rank) {
1177580bdb30SBarry Smith     ierr = PetscMemzero(title,PETSC_MAX_PATH_LEN+1);CHKERRQ(ierr);
11786de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_init,(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs));
1179552f7358SJed Brown     if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n");
1180552f7358SJed Brown   }
1181552f7358SJed Brown   ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
1182552f7358SJed Brown   ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
1183552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr);
1184c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
1185552f7358SJed Brown   ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
1186552f7358SJed Brown 
1187552f7358SJed Brown   /* Read cell sets information */
1188552f7358SJed Brown   if (!rank) {
1189552f7358SJed Brown     PetscInt *cone;
1190e8f6893fSMatthew G. Knepley     int      c, cs, ncs, c_loc, v, v_loc;
1191552f7358SJed Brown     /* Read from ex_get_elem_blk_ids() */
1192e8f6893fSMatthew G. Knepley     int *cs_id, *cs_order;
1193552f7358SJed Brown     /* Read from ex_get_elem_block() */
1194552f7358SJed Brown     char buffer[PETSC_MAX_PATH_LEN+1];
1195e8f6893fSMatthew G. Knepley     int  num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr;
1196552f7358SJed Brown     /* Read from ex_get_elem_conn() */
1197552f7358SJed Brown     int *cs_connect;
1198552f7358SJed Brown 
1199552f7358SJed Brown     /* Get cell sets IDs */
1200e8f6893fSMatthew G. Knepley     ierr = PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order);CHKERRQ(ierr);
12016de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, cs_id));
1202552f7358SJed Brown     /* Read the cell set connectivity table and build mesh topology
1203552f7358SJed Brown        EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
1204e8f6893fSMatthew G. Knepley     /* Check for a hybrid mesh */
1205e8f6893fSMatthew G. Knepley     for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) {
1206e8f6893fSMatthew G. Knepley       PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr));
1207e8f6893fSMatthew G. Knepley       switch (dim) {
1208e8f6893fSMatthew G. Knepley         case 3:
1209e8f6893fSMatthew G. Knepley         switch (num_vertex_per_cell) {
1210e8f6893fSMatthew G. Knepley           case 6:
1211e8f6893fSMatthew G. Knepley             cs_order[cs] = cs;
1212e8f6893fSMatthew G. Knepley             numHybridCells += num_cell_in_set;
1213e8f6893fSMatthew G. Knepley             ++num_hybrid;
1214e8f6893fSMatthew G. Knepley           break;
1215e8f6893fSMatthew G. Knepley           default:
1216e8f6893fSMatthew G. Knepley             for (c = cs; c > cs-num_hybrid; --c) cs_order[c] = cs_order[c-1];
1217e8f6893fSMatthew G. Knepley             cs_order[cs-num_hybrid] = cs;
1218e8f6893fSMatthew G. Knepley         }
1219e8f6893fSMatthew G. Knepley         break;
1220e8f6893fSMatthew G. Knepley       default:
1221e8f6893fSMatthew G. Knepley         for (c = cs; c > cs-num_hybrid; --c) cs_order[c] = cs_order[c-1];
1222e8f6893fSMatthew G. Knepley         cs_order[cs-num_hybrid] = cs;
1223e8f6893fSMatthew G. Knepley       }
1224e8f6893fSMatthew G. Knepley     }
1225e8f6893fSMatthew G. Knepley     if (num_hybrid) {ierr = DMPlexSetHybridBounds(*dm, numCells-numHybridCells, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);}
1226552f7358SJed Brown     /* First set sizes */
1227e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
1228e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
12296de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr));
1230552f7358SJed Brown       for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
1231552f7358SJed Brown         ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr);
1232552f7358SJed Brown       }
1233552f7358SJed Brown     }
1234552f7358SJed Brown     ierr = DMSetUp(*dm);CHKERRQ(ierr);
1235e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
1236e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
12376de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr));
1238dcca6d9dSJed Brown       ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr);
12396de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_conn,(exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL));
1240eaf898f9SPatrick Sanan       /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
1241552f7358SJed Brown       for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
1242552f7358SJed Brown         for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) {
1243552f7358SJed Brown           cone[v_loc] = cs_connect[v]+numCells-1;
1244552f7358SJed Brown         }
1245731dcdf4SMatthew G. Knepley         if (dim == 3) {
12462e1b13c2SMatthew G. Knepley           /* Tetrahedra are inverted */
12472e1b13c2SMatthew G. Knepley           if (num_vertex_per_cell == 4) {
12482e1b13c2SMatthew G. Knepley             PetscInt tmp = cone[0];
12492e1b13c2SMatthew G. Knepley             cone[0] = cone[1];
12502e1b13c2SMatthew G. Knepley             cone[1] = tmp;
12512e1b13c2SMatthew G. Knepley           }
12522e1b13c2SMatthew G. Knepley           /* Hexahedra are inverted */
12532e1b13c2SMatthew G. Knepley           if (num_vertex_per_cell == 8) {
12542e1b13c2SMatthew G. Knepley             PetscInt tmp = cone[1];
12552e1b13c2SMatthew G. Knepley             cone[1] = cone[3];
12562e1b13c2SMatthew G. Knepley             cone[3] = tmp;
12572e1b13c2SMatthew G. Knepley           }
1258731dcdf4SMatthew G. Knepley         }
1259552f7358SJed Brown         ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
1260c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr);
1261552f7358SJed Brown       }
1262552f7358SJed Brown       ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr);
1263552f7358SJed Brown     }
1264e8f6893fSMatthew G. Knepley     ierr = PetscFree2(cs_id, cs_order);CHKERRQ(ierr);
1265552f7358SJed Brown   }
1266552f7358SJed Brown   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
1267552f7358SJed Brown   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
1268552f7358SJed Brown   if (interpolate) {
12695fd9971aSMatthew G. Knepley     DM idm;
1270552f7358SJed Brown 
12719f074e33SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
1272552f7358SJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
1273552f7358SJed Brown     *dm  = idm;
1274552f7358SJed Brown   }
1275552f7358SJed Brown 
1276552f7358SJed Brown   /* Create vertex set label */
1277552f7358SJed Brown   if (!rank && (num_vs > 0)) {
1278552f7358SJed Brown     int vs, v;
1279552f7358SJed Brown     /* Read from ex_get_node_set_ids() */
1280552f7358SJed Brown     int *vs_id;
1281552f7358SJed Brown     /* Read from ex_get_node_set_param() */
1282f958083aSBlaise Bourdin     int num_vertex_in_set;
1283552f7358SJed Brown     /* Read from ex_get_node_set() */
1284552f7358SJed Brown     int *vs_vertex_list;
1285552f7358SJed Brown 
1286552f7358SJed Brown     /* Get vertex set ids */
1287785e854fSJed Brown     ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr);
12886de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_ids,(exoid, EX_NODE_SET, vs_id));
1289552f7358SJed Brown     for (vs = 0; vs < num_vs; ++vs) {
12906de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_set_param,(exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL));
1291785e854fSJed Brown       ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr);
12926de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_set,(exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL));
1293552f7358SJed Brown       for (v = 0; v < num_vertex_in_set; ++v) {
1294c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr);
1295552f7358SJed Brown       }
1296552f7358SJed Brown       ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr);
1297552f7358SJed Brown     }
1298552f7358SJed Brown     ierr = PetscFree(vs_id);CHKERRQ(ierr);
1299552f7358SJed Brown   }
1300552f7358SJed Brown   /* Read coordinates */
130169d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
1302552f7358SJed Brown   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
1303552f7358SJed Brown   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
1304552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
1305552f7358SJed Brown   for (v = numCells; v < numCells+numVertices; ++v) {
1306552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
1307552f7358SJed Brown     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
1308552f7358SJed Brown   }
1309552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1310552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
13118b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
1312552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1313552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
13144e90ef8eSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
13152eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
1316552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
13170aba08f6SKarl Rupp   if (!rank) {
1318e84f7031SBlaise Bourdin     PetscReal *x, *y, *z;
1319552f7358SJed Brown 
1320dcca6d9dSJed Brown     ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr);
13216de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_coord,(exoid, x, y, z));
13220d644c17SKarl Rupp     if (dim > 0) {
13230d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v];
13240d644c17SKarl Rupp     }
13250d644c17SKarl Rupp     if (dim > 1) {
13260d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v];
13270d644c17SKarl Rupp     }
13280d644c17SKarl Rupp     if (dim > 2) {
13290d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v];
13300d644c17SKarl Rupp     }
1331552f7358SJed Brown     ierr = PetscFree3(x,y,z);CHKERRQ(ierr);
1332552f7358SJed Brown   }
1333552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1334552f7358SJed Brown   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
1335552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1336552f7358SJed Brown 
1337552f7358SJed Brown   /* Create side set label */
1338552f7358SJed Brown   if (!rank && interpolate && (num_fs > 0)) {
1339552f7358SJed Brown     int fs, f, voff;
1340552f7358SJed Brown     /* Read from ex_get_side_set_ids() */
1341552f7358SJed Brown     int *fs_id;
1342552f7358SJed Brown     /* Read from ex_get_side_set_param() */
1343f958083aSBlaise Bourdin     int num_side_in_set;
1344552f7358SJed Brown     /* Read from ex_get_side_set_node_list() */
1345552f7358SJed Brown     int *fs_vertex_count_list, *fs_vertex_list;
1346ef073283Smichael_afanasiev     /* Read side set labels */
1347ef073283Smichael_afanasiev     char fs_name[MAX_STR_LENGTH+1];
1348552f7358SJed Brown 
1349552f7358SJed Brown     /* Get side set ids */
1350785e854fSJed Brown     ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr);
13516de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_ids,(exoid, EX_SIDE_SET, fs_id));
1352552f7358SJed Brown     for (fs = 0; fs < num_fs; ++fs) {
13536de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_set_param,(exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL));
1354dcca6d9dSJed Brown       ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr);
13556de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_side_set_node_list,(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list));
1356ef073283Smichael_afanasiev       /* Get the specific name associated with this side set ID. */
1357ef073283Smichael_afanasiev       int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name);
1358552f7358SJed Brown       for (f = 0, voff = 0; f < num_side_in_set; ++f) {
13590298fd71SBarry Smith         const PetscInt *faces   = NULL;
1360552f7358SJed Brown         PetscInt       faceSize = fs_vertex_count_list[f], numFaces;
1361552f7358SJed Brown         PetscInt       faceVertices[4], v;
1362552f7358SJed Brown 
1363552f7358SJed Brown         if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize);
1364552f7358SJed Brown         for (v = 0; v < faceSize; ++v, ++voff) {
1365552f7358SJed Brown           faceVertices[v] = fs_vertex_list[voff]+numCells-1;
1366552f7358SJed Brown         }
1367552f7358SJed Brown         ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
1368552f7358SJed Brown         if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces);
1369c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr);
1370ef073283Smichael_afanasiev         /* Only add the label if one has been detected for this side set. */
1371ef073283Smichael_afanasiev         if (!fs_name_err) {
1372ef073283Smichael_afanasiev           ierr = DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]);CHKERRQ(ierr);
1373ef073283Smichael_afanasiev         }
1374552f7358SJed Brown         ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
1375552f7358SJed Brown       }
1376552f7358SJed Brown       ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr);
1377552f7358SJed Brown     }
1378552f7358SJed Brown     ierr = PetscFree(fs_id);CHKERRQ(ierr);
1379552f7358SJed Brown   }
1380552f7358SJed Brown #else
1381552f7358SJed Brown   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
1382552f7358SJed Brown #endif
1383552f7358SJed Brown   PetscFunctionReturn(0);
1384552f7358SJed Brown }
1385