xref: /petsc/src/dm/impls/plex/plexexodusii.c (revision 33751fbdd5231c54fff9990cabc9db80831a50fa)
1552f7358SJed Brown #define PETSCDM_DLL
234541f0dSBarry 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 
9552f7358SJed Brown #undef __FUNCT__
10*33751fbdSMichael Lange #define __FUNCT__ "DMPlexCreateExodusFromFile"
11*33751fbdSMichael Lange /*@C
12*33751fbdSMichael Lange   DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file.
13*33751fbdSMichael Lange 
14*33751fbdSMichael Lange   Collective on comm
15*33751fbdSMichael Lange 
16*33751fbdSMichael Lange   Input Parameters:
17*33751fbdSMichael Lange + comm  - The MPI communicator
18*33751fbdSMichael Lange . filename - The name of the ExodusII file
19*33751fbdSMichael Lange - interpolate - Create faces and edges in the mesh
20*33751fbdSMichael Lange 
21*33751fbdSMichael Lange   Output Parameter:
22*33751fbdSMichael Lange . dm  - The DM object representing the mesh
23*33751fbdSMichael Lange 
24*33751fbdSMichael Lange   Level: beginner
25*33751fbdSMichael Lange 
26*33751fbdSMichael Lange .keywords: mesh,ExodusII
27*33751fbdSMichael Lange .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus()
28*33751fbdSMichael Lange @*/
29*33751fbdSMichael Lange PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
30*33751fbdSMichael Lange {
31*33751fbdSMichael Lange   PetscMPIInt    rank;
32*33751fbdSMichael Lange   PetscErrorCode ierr;
33*33751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
34*33751fbdSMichael Lange   int   CPU_word_size = 0, IO_word_size = 0, exoid = -1;
35*33751fbdSMichael Lange   float version;
36*33751fbdSMichael Lange #endif
37*33751fbdSMichael Lange 
38*33751fbdSMichael Lange   PetscFunctionBegin;
39*33751fbdSMichael Lange   PetscValidCharPointer(filename, 2);
40*33751fbdSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
41*33751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
42*33751fbdSMichael Lange   if (!rank) {
43*33751fbdSMichael Lange     exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version);
44*33751fbdSMichael Lange     if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename);
45*33751fbdSMichael Lange   }
46*33751fbdSMichael Lange   ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr);
47*33751fbdSMichael Lange   if (!rank) {ierr = ex_close(exoid);CHKERRQ(ierr);}
48*33751fbdSMichael Lange #else
49*33751fbdSMichael Lange   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
50*33751fbdSMichael Lange #endif
51*33751fbdSMichael Lange   PetscFunctionReturn(0);
52*33751fbdSMichael Lange }
53*33751fbdSMichael Lange 
54*33751fbdSMichael Lange #undef __FUNCT__
55552f7358SJed Brown #define __FUNCT__ "DMPlexCreateExodus"
56552f7358SJed Brown /*@
57*33751fbdSMichael Lange   DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID.
58552f7358SJed Brown 
59552f7358SJed Brown   Collective on comm
60552f7358SJed Brown 
61552f7358SJed Brown   Input Parameters:
62552f7358SJed Brown + comm  - The MPI communicator
63552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
64552f7358SJed Brown - interpolate - Create faces and edges in the mesh
65552f7358SJed Brown 
66552f7358SJed Brown   Output Parameter:
67552f7358SJed Brown . dm  - The DM object representing the mesh
68552f7358SJed Brown 
69552f7358SJed Brown   Level: beginner
70552f7358SJed Brown 
71552f7358SJed Brown .keywords: mesh,ExodusII
72e84e7769SJed Brown .seealso: DMPLEX, DMCreate()
73552f7358SJed Brown @*/
74552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm)
75552f7358SJed Brown {
76552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
77552f7358SJed Brown   PetscMPIInt    num_proc, rank;
78552f7358SJed Brown   PetscSection   coordSection;
79552f7358SJed Brown   Vec            coordinates;
80552f7358SJed Brown   PetscScalar    *coords;
81552f7358SJed Brown   PetscInt       coordSize, v;
82552f7358SJed Brown   PetscErrorCode ierr;
83552f7358SJed Brown   /* Read from ex_get_init() */
84552f7358SJed Brown   char title[PETSC_MAX_PATH_LEN+1];
85552f7358SJed Brown   int  dim    = 0, numVertices = 0, numCells = 0;
86552f7358SJed Brown   int  num_cs = 0, num_vs = 0, num_fs = 0;
87552f7358SJed Brown #endif
88552f7358SJed Brown 
89552f7358SJed Brown   PetscFunctionBegin;
90552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
91552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
92552f7358SJed Brown   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
93552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
94552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
95552f7358SJed Brown   /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */
96552f7358SJed Brown   if (!rank) {
9739bba695SBarry Smith     ierr = PetscMemzero(title,(PETSC_MAX_PATH_LEN+1)*sizeof(char));CHKERRQ(ierr);
98552f7358SJed Brown     ierr = ex_get_init(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs);
99552f7358SJed Brown     if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ExodusII ex_get_init() failed with error code %D\n",ierr);
100552f7358SJed Brown     if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n");
101552f7358SJed Brown   }
102552f7358SJed Brown   ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
103552f7358SJed Brown   ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
104552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr);
105552f7358SJed Brown   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
106552f7358SJed Brown   ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
107552f7358SJed Brown 
108552f7358SJed Brown   /* Read cell sets information */
109552f7358SJed Brown   if (!rank) {
110552f7358SJed Brown     PetscInt *cone;
111552f7358SJed Brown     int      c, cs, c_loc, v, v_loc;
112552f7358SJed Brown     /* Read from ex_get_elem_blk_ids() */
113552f7358SJed Brown     int *cs_id;
114552f7358SJed Brown     /* Read from ex_get_elem_block() */
115552f7358SJed Brown     char buffer[PETSC_MAX_PATH_LEN+1];
116552f7358SJed Brown     int  num_cell_in_set, num_vertex_per_cell, num_attr;
117552f7358SJed Brown     /* Read from ex_get_elem_conn() */
118552f7358SJed Brown     int *cs_connect;
119552f7358SJed Brown 
120552f7358SJed Brown     /* Get cell sets IDs */
121785e854fSJed Brown     ierr = PetscMalloc1(num_cs, &cs_id);CHKERRQ(ierr);
122552f7358SJed Brown     ierr = ex_get_elem_blk_ids(exoid, cs_id);CHKERRQ(ierr);
123552f7358SJed Brown     /* Read the cell set connectivity table and build mesh topology
124552f7358SJed Brown        EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
125552f7358SJed Brown     /* First set sizes */
126552f7358SJed Brown     for (cs = 0, c = 0; cs < num_cs; cs++) {
127552f7358SJed Brown       ierr = ex_get_elem_block(exoid, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, &num_attr);CHKERRQ(ierr);
128552f7358SJed Brown       for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
129552f7358SJed Brown         ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr);
130552f7358SJed Brown       }
131552f7358SJed Brown     }
132552f7358SJed Brown     ierr = DMSetUp(*dm);CHKERRQ(ierr);
133552f7358SJed Brown     for (cs = 0, c = 0; cs < num_cs; cs++) {
134552f7358SJed Brown       ierr = ex_get_elem_block(exoid, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, &num_attr);CHKERRQ(ierr);
135dcca6d9dSJed Brown       ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr);
136552f7358SJed Brown       ierr = ex_get_elem_conn(exoid, cs_id[cs], cs_connect);CHKERRQ(ierr);
137552f7358SJed Brown       /* EXO uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */
138552f7358SJed Brown       for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
139552f7358SJed Brown         for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) {
140552f7358SJed Brown           cone[v_loc] = cs_connect[v]+numCells-1;
141552f7358SJed Brown         }
142731dcdf4SMatthew G. Knepley         if (dim == 3) {
1432e1b13c2SMatthew G. Knepley           /* Tetrahedra are inverted */
1442e1b13c2SMatthew G. Knepley           if (num_vertex_per_cell == 4) {
1452e1b13c2SMatthew G. Knepley             PetscInt tmp = cone[0];
1462e1b13c2SMatthew G. Knepley             cone[0] = cone[1];
1472e1b13c2SMatthew G. Knepley             cone[1] = tmp;
1482e1b13c2SMatthew G. Knepley           }
1492e1b13c2SMatthew G. Knepley           /* Hexahedra are inverted */
1502e1b13c2SMatthew G. Knepley           if (num_vertex_per_cell == 8) {
1512e1b13c2SMatthew G. Knepley             PetscInt tmp = cone[1];
1522e1b13c2SMatthew G. Knepley             cone[1] = cone[3];
1532e1b13c2SMatthew G. Knepley             cone[3] = tmp;
1542e1b13c2SMatthew G. Knepley           }
155731dcdf4SMatthew G. Knepley         }
156552f7358SJed Brown         ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
157552f7358SJed Brown         ierr = DMPlexSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr);
158552f7358SJed Brown       }
159552f7358SJed Brown       ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr);
160552f7358SJed Brown     }
161552f7358SJed Brown     ierr = PetscFree(cs_id);CHKERRQ(ierr);
162552f7358SJed Brown   }
163552f7358SJed Brown   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
164552f7358SJed Brown   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
165552f7358SJed Brown   if (interpolate) {
166552f7358SJed Brown     DM idm;
167552f7358SJed Brown 
1689f074e33SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
169552f7358SJed Brown     /* Maintain Cell Sets label */
170552f7358SJed Brown     {
171552f7358SJed Brown       DMLabel label;
172552f7358SJed Brown 
173552f7358SJed Brown       ierr = DMPlexRemoveLabel(*dm, "Cell Sets", &label);CHKERRQ(ierr);
174552f7358SJed Brown       if (label) {ierr = DMPlexAddLabel(idm, label);CHKERRQ(ierr);}
175552f7358SJed Brown     }
176552f7358SJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
177552f7358SJed Brown     *dm  = idm;
178552f7358SJed Brown   }
179552f7358SJed Brown 
180552f7358SJed Brown   /* Create vertex set label */
181552f7358SJed Brown   if (!rank && (num_vs > 0)) {
182552f7358SJed Brown     int vs, v;
183552f7358SJed Brown     /* Read from ex_get_node_set_ids() */
184552f7358SJed Brown     int *vs_id;
185552f7358SJed Brown     /* Read from ex_get_node_set_param() */
186552f7358SJed Brown     int num_vertex_in_set, num_attr;
187552f7358SJed Brown     /* Read from ex_get_node_set() */
188552f7358SJed Brown     int *vs_vertex_list;
189552f7358SJed Brown 
190552f7358SJed Brown     /* Get vertex set ids */
191785e854fSJed Brown     ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr);
192552f7358SJed Brown     ierr = ex_get_node_set_ids(exoid, vs_id);CHKERRQ(ierr);
193552f7358SJed Brown     for (vs = 0; vs < num_vs; ++vs) {
194552f7358SJed Brown       ierr = ex_get_node_set_param(exoid, vs_id[vs], &num_vertex_in_set, &num_attr);CHKERRQ(ierr);
195785e854fSJed Brown       ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr);
196552f7358SJed Brown       ierr = ex_get_node_set(exoid, vs_id[vs], vs_vertex_list);CHKERRQ(ierr);
197552f7358SJed Brown       for (v = 0; v < num_vertex_in_set; ++v) {
198552f7358SJed Brown         ierr = DMPlexSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr);
199552f7358SJed Brown       }
200552f7358SJed Brown       ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr);
201552f7358SJed Brown     }
202552f7358SJed Brown     ierr = PetscFree(vs_id);CHKERRQ(ierr);
203552f7358SJed Brown   }
204552f7358SJed Brown   /* Read coordinates */
20569d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
206552f7358SJed Brown   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
207552f7358SJed Brown   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
208552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
209552f7358SJed Brown   for (v = numCells; v < numCells+numVertices; ++v) {
210552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
211552f7358SJed Brown     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
212552f7358SJed Brown   }
213552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
214552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
215552f7358SJed Brown   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
216552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
217552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
2182eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
219552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
2200aba08f6SKarl Rupp   if (!rank) {
221552f7358SJed Brown     float *x, *y, *z;
222552f7358SJed Brown 
223dcca6d9dSJed Brown     ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr);
224552f7358SJed Brown     ierr = ex_get_coord(exoid, x, y, z);CHKERRQ(ierr);
2250d644c17SKarl Rupp     if (dim > 0) {
2260d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v];
2270d644c17SKarl Rupp     }
2280d644c17SKarl Rupp     if (dim > 1) {
2290d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v];
2300d644c17SKarl Rupp     }
2310d644c17SKarl Rupp     if (dim > 2) {
2320d644c17SKarl Rupp       for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v];
2330d644c17SKarl Rupp     }
234552f7358SJed Brown     ierr = PetscFree3(x,y,z);CHKERRQ(ierr);
235552f7358SJed Brown   }
236552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
237552f7358SJed Brown   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
238552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
239552f7358SJed Brown 
240552f7358SJed Brown   /* Create side set label */
241552f7358SJed Brown   if (!rank && interpolate && (num_fs > 0)) {
242552f7358SJed Brown     int fs, f, voff;
243552f7358SJed Brown     /* Read from ex_get_side_set_ids() */
244552f7358SJed Brown     int *fs_id;
245552f7358SJed Brown     /* Read from ex_get_side_set_param() */
246552f7358SJed Brown     int num_side_in_set, num_dist_fact_in_set;
247552f7358SJed Brown     /* Read from ex_get_side_set_node_list() */
248552f7358SJed Brown     int *fs_vertex_count_list, *fs_vertex_list;
249552f7358SJed Brown 
250552f7358SJed Brown     /* Get side set ids */
251785e854fSJed Brown     ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr);
252552f7358SJed Brown     ierr = ex_get_side_set_ids(exoid, fs_id);CHKERRQ(ierr);
253552f7358SJed Brown     for (fs = 0; fs < num_fs; ++fs) {
254552f7358SJed Brown       ierr = ex_get_side_set_param(exoid, fs_id[fs], &num_side_in_set, &num_dist_fact_in_set);CHKERRQ(ierr);
255dcca6d9dSJed Brown       ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr);
256552f7358SJed Brown       ierr = ex_get_side_set_node_list(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);CHKERRQ(ierr);
257552f7358SJed Brown       for (f = 0, voff = 0; f < num_side_in_set; ++f) {
2580298fd71SBarry Smith         const PetscInt *faces   = NULL;
259552f7358SJed Brown         PetscInt       faceSize = fs_vertex_count_list[f], numFaces;
260552f7358SJed Brown         PetscInt       faceVertices[4], v;
261552f7358SJed Brown 
262552f7358SJed Brown         if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize);
263552f7358SJed Brown         for (v = 0; v < faceSize; ++v, ++voff) {
264552f7358SJed Brown           faceVertices[v] = fs_vertex_list[voff]+numCells-1;
265552f7358SJed Brown         }
266552f7358SJed Brown         ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
267552f7358SJed Brown         if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces);
268552f7358SJed Brown         ierr = DMPlexSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr);
269552f7358SJed Brown         ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
270552f7358SJed Brown       }
271552f7358SJed Brown       ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr);
272552f7358SJed Brown     }
273552f7358SJed Brown     ierr = PetscFree(fs_id);CHKERRQ(ierr);
274552f7358SJed Brown   }
275552f7358SJed Brown #else
276552f7358SJed Brown   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
277552f7358SJed Brown #endif
278552f7358SJed Brown   PetscFunctionReturn(0);
279552f7358SJed Brown }
280