Actual source code: meshexodus.c
1: #include<petscmesh_formats.hh> /*I "petscmesh.h" I*/
3: #ifdef PETSC_HAVE_EXODUSII
5: // This is what I needed in my petscvariables:
6: //
7: // EXODUSII_INCLUDE = -I/PETSc3/mesh/exodusii-4.71/cbind/include
8: // EXODUSII_LIB = -L/PETSc3/mesh/exodusii-4.71/forbind/src -lexoIIv2for -L/PETSc3/mesh/exodusii-4.71/cbind/src -lexoIIv2c -lnetcdf
10: #include<netcdf.h>
11: #include<exodusII.h>
15: PetscErrorCode PetscReadExodusII(MPI_Comm comm, const char filename[], ALE::Obj<PETSC_MESH_TYPE>& mesh)
16: {
17: int exoid;
18: int CPU_word_size = 0, IO_word_size = 0;
19: const PetscMPIInt rank = mesh->commRank();
20: float version;
21: char title[MAX_LINE_LENGTH+1], elem_type[MAX_STR_LENGTH+1];
22: int num_dim, num_nodes, num_elem, num_elem_blk, num_node_sets, num_side_sets;
23: int ierr;
26: // Open EXODUS II file
27: exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version);CHKERRQ(!exoid);
28: // Read database parameters
29: ex_get_init(exoid, title, &num_dim, &num_nodes, &num_elem, &num_elem_blk, &num_node_sets, &num_side_sets);
31: // Read vertex coordinates
32: float *x, *y, *z;
33: PetscMalloc3(num_nodes,float,&x,num_nodes,float,&y,num_nodes,float,&z);
34: ex_get_coord(exoid, x, y, z);
36: // Read element connectivity
37: int *eb_ids, *num_elem_in_block, *num_nodes_per_elem, *num_attr;
38: int **connect;
39: char **block_names;
40: if (num_elem_blk > 0) {
41: PetscMalloc5(num_elem_blk,int,&eb_ids,num_elem_blk,int,&num_elem_in_block,num_elem_blk,int,&num_nodes_per_elem,num_elem_blk,int,&num_attr,num_elem_blk,char*,&block_names);
42: ex_get_elem_blk_ids(exoid, eb_ids);
43: for(int eb = 0; eb < num_elem_blk; ++eb) {
44: PetscMalloc((MAX_STR_LENGTH+1) * sizeof(char), &block_names[eb]);
45: }
46: ex_get_names(exoid, EX_ELEM_BLOCK, block_names);
47: for(int eb = 0; eb < num_elem_blk; ++eb) {
48: ex_get_elem_block(exoid, eb_ids[eb], elem_type, &num_elem_in_block[eb], &num_nodes_per_elem[eb], &num_attr[eb]);
49: PetscFree(block_names[eb]);
50: }
51: PetscMalloc(num_elem_blk * sizeof(int*),&connect);
52: for(int eb = 0; eb < num_elem_blk; ++eb) {
53: if (num_elem_in_block[eb] > 0) {
54: PetscMalloc(num_nodes_per_elem[eb]*num_elem_in_block[eb] * sizeof(int),&connect[eb]);
55: ex_get_elem_conn(exoid, eb_ids[eb], connect[eb]);
56: }
57: }
58: }
60: // Read node sets
61: int *ns_ids, *num_nodes_in_set;
62: int **node_list;
63: if (num_node_sets > 0) {
64: PetscMalloc3(num_node_sets,int,&ns_ids,num_node_sets,int,&num_nodes_in_set,num_node_sets,int*,&node_list);
65: ex_get_node_set_ids(exoid, ns_ids);
66: for(int ns = 0; ns < num_node_sets; ++ns) {
67: int num_df_in_set;
68: ex_get_node_set_param (exoid, ns_ids[ns], &num_nodes_in_set[ns], &num_df_in_set);
69: PetscMalloc(num_nodes_in_set[ns] * sizeof(int), &node_list[ns]);
70: ex_get_node_set(exoid, ns_ids[ns], node_list[ns]);
71: }
72: }
73: ex_close(exoid);
75: // Build mesh topology
76: int numCorners = num_nodes_per_elem[0];
77: int *cells;
78: mesh->setDimension(num_dim);
79: PetscMalloc(numCorners*num_elem * sizeof(int), &cells);
80: for(int eb = 0, k = 0; eb < num_elem_blk; ++eb) {
81: for(int e = 0; e < num_elem_in_block[eb]; ++e, ++k) {
82: for(int c = 0; c < numCorners; ++c) {
83: cells[k*numCorners+c] = connect[eb][e*numCorners+c];
84: }
85: }
86: PetscFree(connect[eb]);
87: }
88: ALE::Obj<PETSC_MESH_TYPE::sieve_type> sieve = new PETSC_MESH_TYPE::sieve_type(mesh->comm(), mesh->debug());
89: bool interpolate = false;
91: try {
92: mesh->setSieve(sieve);
93: if (0 == rank) {
94: if (!interpolate) {
95: // Create the ISieve
96: sieve->setChart(PETSC_MESH_TYPE::sieve_type::chart_type(0, num_elem+num_nodes));
97: // Set cone and support sizes
98: for (int c = 0; c < num_elem; ++c) {
99: sieve->setConeSize(c, numCorners);
100: }
101: sieve->symmetrizeSizes(num_elem, numCorners, cells, num_elem - 1); /* Notice the -1 for 1-based indexing in cells[] */
102: // Allocate point storage
103: sieve->allocate();
104: // Fill up cones
105: int *cone = new int[numCorners];
106: int *coneO = new int[numCorners];
107: for (int v = 0; v < numCorners; ++v) {
108: coneO[v] = 1;
109: }
110: for (int c = 0; c < num_elem; ++c) {
111: for (int v = 0; v < numCorners; ++v) {
112: cone[v] = cells[c*numCorners+v]+num_elem - 1;
113: }
114: sieve->setCone(cone, c);
115: sieve->setConeOrientation(coneO, c);
116: } // for
117: delete[] cone; cone = 0;
118: delete[] coneO; coneO = 0;
119: // Symmetrize to fill up supports
120: sieve->symmetrize();
121: } else {
122: // Same old thing
123: ALE::Obj<ALE::Mesh::sieve_type> s = new ALE::Mesh::sieve_type(sieve->comm(), sieve->debug());
125: ALE::SieveBuilder<ALE::Mesh>::buildTopology(s, num_dim, num_elem, cells, num_nodes, interpolate, numCorners);
126: std::map<ALE::Mesh::point_type,ALE::Mesh::point_type> renumbering;
127: ALE::ISieveConverter::convertSieve(*s, *sieve, renumbering);
128: }
129: if (!interpolate) {
130: // Optimized stratification
131: const ALE::Obj<PETSC_MESH_TYPE::label_type>& height = mesh->createLabel("height");
132: const ALE::Obj<PETSC_MESH_TYPE::label_type>& depth = mesh->createLabel("depth");
134: for(int c = 0; c < num_elem; ++c) {
135: height->setCone(0, c);
136: depth->setCone(1, c);
137: }
138: for(int v = num_elem; v < num_elem+num_nodes; ++v) {
139: height->setCone(1, v);
140: depth->setCone(0, v);
141: }
142: mesh->setHeight(1);
143: mesh->setDepth(1);
144: } else {
145: mesh->stratify();
146: }
147: } else {
148: mesh->getSieve()->setChart(PETSC_MESH_TYPE::sieve_type::chart_type());
149: mesh->getSieve()->allocate();
150: mesh->stratify();
151: }
152: } catch (ALE::Exception e) {
153: SETERRQ(PETSC_ERR_LIB, e.msg().c_str());
154: }
155: PetscFree(cells);
157: // Build cell blocks
158: const ALE::Obj<PETSC_MESH_TYPE::label_type>& cellBlocks = mesh->createLabel("CellBlocks");
159: if (rank == 0) {
160: for(int eb = 0, k = 0; eb < num_elem_blk; ++eb) {
161: for(int e = 0; e < num_elem_in_block[eb]; ++e, ++k) {
162: mesh->setValue(cellBlocks, k, eb_ids[eb]);
163: }
164: }
165: }
166: if (num_elem_blk > 0) {
167: PetscFree(connect);
168: PetscFree5(eb_ids, num_elem_in_block, num_nodes_per_elem, num_attr, block_names);
169: }
171: // Build coordinates
172: double *coords;
173: PetscMalloc(num_dim*num_nodes * sizeof(double), &coords);
174: if (num_dim > 0) {for(int v = 0; v < num_nodes; ++v) {coords[v*num_dim+0] = x[v];}}
175: if (num_dim > 1) {for(int v = 0; v < num_nodes; ++v) {coords[v*num_dim+1] = y[v];}}
176: if (num_dim > 2) {for(int v = 0; v < num_nodes; ++v) {coords[v*num_dim+2] = z[v];}}
177: ALE::SieveBuilder<PETSC_MESH_TYPE>::buildCoordinates(mesh, num_dim, coords);
178: PetscFree(coords);
179: PetscFree3(x, y, z);
181: // Build vertex sets
182: const ALE::Obj<PETSC_MESH_TYPE::label_type>& vertexSets = mesh->createLabel("VertexSets");
183: if (rank == 0) {
184: for(int ns = 0; ns < num_node_sets; ++ns) {
185: for(int n = 0; n < num_nodes_in_set[ns]; ++n) {
186: mesh->addValue(vertexSets, node_list[ns][n]+num_elem-1, ns_ids[ns]);
187: }
188: PetscFree(node_list[ns]);
189: }
190: }
191: if (num_node_sets > 0) {
192: PetscFree3(ns_ids,num_nodes_in_set,node_list);
193: }
195: //cellBlocks->view("Cell Blocks");
196: //vertexSets->view("Vertex Sets");
198: // Get coords and print in F90
199: // Get connectivity and print in F90
200: // Calculate cost function
201: // Do in parallel
202: // Read in parallel
203: // Distribute
204: // Print out local meshes
205: // Do Blaise's assembly loop in parallel
206: // Assemble function into Section
207: // Assemble jacobian into Mat
208: // Assemble in parallel
209: // Convert Section to Vec
210: return(0);
211: }
213: #endif // PETSC_HAVE_EXODUSII
217: /*@C
218: MeshCreateExodus - Create a Mesh from an ExodusII file.
220: Not Collective
222: Input Parameters:
223: + comm - The MPI communicator
224: - filename - The ExodusII filename
226: Output Parameter:
227: . mesh - The Mesh object
229: Level: beginner
231: .keywords: mesh, ExodusII
232: .seealso: MeshCreate()
233: @*/
234: PetscErrorCode MeshCreateExodus(MPI_Comm comm, const char filename[], Mesh *mesh)
235: {
236: PetscInt debug = 0;
237: PetscTruth flag;
241: MeshCreate(comm, mesh);
242: PetscOptionsGetInt(PETSC_NULL, "-debug", &debug, &flag);
243: ALE::Obj<PETSC_MESH_TYPE> m = new PETSC_MESH_TYPE(comm, -1, debug);
244: #ifdef PETSC_HAVE_EXODUSII
245: PetscReadExodusII(comm, filename, m);
246: #else
247: SETERRQ(PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --with-exodus-dir=/path/to/exodus");
248: #endif
249: if (debug) {m->view("Mesh");}
250: MeshSetMesh(*mesh, m);
251: return(0);
252: }
256: /*@C
257: MeshExodusGetInfo - Get information about an ExodusII Mesh.
259: Not Collective
261: Input Parameter:
262: . mesh - The Mesh object
264: Output Parameters:
265: + dim - The mesh dimension
266: . numVertices - The number of vertices in the mesh
267: . numCells - The number of cells in the mesh
268: . numCellBlocks - The number of cell blocks in the mesh
269: - numVertexSets - The number of vertex sets in the mesh
271: Level: beginner
273: .keywords: mesh, ExodusII
274: .seealso: MeshCreateExodus()
275: @*/
276: PetscErrorCode MeshExodusGetInfo(Mesh mesh, PetscInt *dim, PetscInt *numVertices, PetscInt *numCells, PetscInt *numCellBlocks, PetscInt *numVertexSets)
277: {
278: ALE::Obj<PETSC_MESH_TYPE> m;
279: PetscErrorCode ierr;
282: MeshGetMesh(mesh, m);
283: *dim = m->getDimension();
284: *numVertices = m->depthStratum(0)->size();
285: *numCells = m->heightStratum(0)->size();
286: *numCellBlocks = m->getLabel("CellBlocks")->getCapSize();
287: *numVertexSets = m->getLabel("VertexSets")->getCapSize();
288: return(0);
289: }
293: /*@C
294: MeshGetLabelSize - Get the number of different integer ids in a Label
296: Not Collective
298: Input Parameters:
299: + mesh - The Mesh object
300: - name - The label name
302: Output Parameter:
303: . size - The label size (number of different integer ids)
305: Level: beginner
307: .keywords: mesh, ExodusII
308: .seealso: MeshCreateExodus()
309: @*/
310: PetscErrorCode MeshGetLabelSize(Mesh mesh, const char name[], PetscInt *size)
311: {
312: ALE::Obj<PETSC_MESH_TYPE> m;
313: PetscErrorCode ierr;
316: MeshGetMesh(mesh, m);
317: *size = m->getLabel(name)->getCapSize();
318: return(0);
319: }
323: /*@C
324: MeshGetLabelIds - Get the integer ids in a label
326: Not Collective
328: Input Parameters:
329: + mesh - The Mesh object
330: . name - The label name
331: - ids - The id storage array
333: Output Parameter:
334: . ids - The integer ids
336: Level: beginner
338: .keywords: mesh, ExodusII
339: .seealso: MeshCreateExodus()
340: @*/
341: PetscErrorCode MeshGetLabelIds(Mesh mesh, const char name[], PetscInt *ids)
342: {
343: ALE::Obj<PETSC_MESH_TYPE> m;
344: PetscErrorCode ierr;
347: MeshGetMesh(mesh, m);
348: const ALE::Obj<PETSC_MESH_TYPE::label_type::capSequence>& labelIds = m->getLabel(name)->cap();
349: const PETSC_MESH_TYPE::label_type::capSequence::const_iterator iEnd = labelIds->end();
350: PetscInt i = 0;
352: for(PETSC_MESH_TYPE::label_type::capSequence::const_iterator i_iter = labelIds->begin(); i_iter != iEnd; ++i_iter, ++i) {
353: ids[i] = *i_iter;
354: }
355: return(0);
356: }
360: /*@C
361: MeshGetStratumSize - Get the number of points in a label stratum
363: Not Collective
365: Input Parameters:
366: + mesh - The Mesh object
367: . name - The label name
368: - value - The stratum value
370: Output Parameter:
371: . size - The stratum size
373: Level: beginner
375: .keywords: mesh, ExodusII
376: .seealso: MeshCreateExodus()
377: @*/
378: PetscErrorCode MeshGetStratumSize(Mesh mesh, const char name[], PetscInt value, PetscInt *size)
379: {
380: ALE::Obj<PETSC_MESH_TYPE> m;
381: PetscErrorCode ierr;
384: MeshGetMesh(mesh, m);
385: *size = m->getLabelStratum(name, value)->size();
386: return(0);
387: }
391: /*@C
392: MeshGetStratum - Get the points in a label stratum
394: Not Collective
396: Input Parameters:
397: + mesh - The Mesh object
398: . name - The label name
399: . value - The stratum value
400: - points - The stratum points storage array
402: Output Parameter:
403: . points - The stratum points
405: Level: beginner
407: .keywords: mesh, ExodusII
408: .seealso: MeshCreateExodus()
409: @*/
410: PetscErrorCode MeshGetStratum(Mesh mesh, const char name[], PetscInt value, PetscInt *points)
411: {
412: ALE::Obj<PETSC_MESH_TYPE> m;
413: PetscErrorCode ierr;
416: MeshGetMesh(mesh, m);
417: const ALE::Obj<PETSC_MESH_TYPE::label_sequence>& stratum = m->getLabelStratum(name, value);
418: const PETSC_MESH_TYPE::label_sequence::iterator sEnd = stratum->end();
419: PetscInt s = 0;
421: for(PETSC_MESH_TYPE::label_sequence::iterator s_iter = stratum->begin(); s_iter != sEnd; ++s_iter, ++s) {
422: points[s] = *s_iter;
423: }
424: return(0);
425: }