Actual source code: meshexodus.c

  1: #include<petscmesh_formats.hh>   /*I      "petscmesh.h"   I*/

  3: #ifdef PETSC_HAVE_EXODUS

  5: // This is what I needed in my petscvariables:
  6: //
  7: // EXODUS_INCLUDE = -I/PETSc3/mesh/exodusii-4.71/cbind/include
  8: // EXODUS_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:   ALE::Obj<ALE::Mesh::sieve_type>       s     = new ALE::Mesh::sieve_type(mesh->comm(), mesh->debug());
 90:   ALE::SieveBuilder<ALE::Mesh>::buildTopology(s, num_dim, num_elem, cells, num_nodes, false, numCorners);
 91:   std::map<PETSC_MESH_TYPE::point_type,PETSC_MESH_TYPE::point_type> renumbering;
 92:   ALE::ISieveConverter::convertSieve(*s, *sieve, renumbering);
 93:   mesh->setSieve(sieve);
 94:   mesh->stratify();
 95:   PetscFree(cells);

 97:   // Build cell blocks
 98:   const ALE::Obj<PETSC_MESH_TYPE::label_type>& cellBlocks = mesh->createLabel("CellBlocks");
 99:   if (rank == 0) {
100:     for(int eb = 0, k = 0; eb < num_elem_blk; ++eb) {
101:       for(int e = 0; e < num_elem_in_block[eb]; ++e, ++k) {
102:         mesh->setValue(cellBlocks, k, eb);
103:       }
104:     }
105:   }
106:   if (num_elem_blk > 0) {
107:     PetscFree(connect);
108:     PetscFree5(eb_ids, num_elem_in_block, num_nodes_per_elem, num_attr, block_names);
109:   }

111:   // Build coordinates
112:   double *coords;
113:   PetscMalloc(num_dim*num_nodes * sizeof(double), &coords);
114:   if (num_dim > 0) {for(int v = 0; v < num_nodes; ++v) {coords[v*num_dim+0] = x[v];}}
115:   if (num_dim > 1) {for(int v = 0; v < num_nodes; ++v) {coords[v*num_dim+1] = y[v];}}
116:   if (num_dim > 2) {for(int v = 0; v < num_nodes; ++v) {coords[v*num_dim+2] = z[v];}}
117:   ALE::SieveBuilder<PETSC_MESH_TYPE>::buildCoordinates(mesh, num_dim, coords);
118:   PetscFree(coords);
119:   PetscFree3(x, y, z);

121:   // Build vertex sets
122:   const ALE::Obj<PETSC_MESH_TYPE::label_type>& vertexSets = mesh->createLabel("VertexSets");
123:   if (rank == 0) {
124:     for(int ns = 0; ns < num_node_sets; ++ns) {
125:       for(int n = 0; n < num_nodes_in_set[ns]; ++n) {
126:         mesh->setValue(vertexSets, node_list[ns][n]+num_elem, ns);
127:       }
128:       PetscFree(node_list[ns]);
129:     }
130:   }
131:   if (num_node_sets > 0) {
132:     PetscFree3(ns_ids,num_nodes_in_set,node_list);
133:   }

135:   //cellBlocks->view("Cell Blocks");
136:   //vertexSets->view("Vertex Sets");

138:   // Get coords and print in F90
139:   // Get connectivity and print in F90
140:   // Calculate cost function
141:   // Do in parallel
142:   //   Read in parallel
143:   //   Distribute
144:   //   Print out local meshes
145:   //   Do Blaise's assembly loop in parallel
146:   // Assemble function into Section
147:   // Assemble jacobian into Mat
148:   // Assemble in parallel
149:   // Convert Section to Vec
150:   return(0);
151: }

153: #endif // PETSC_HAVE_EXODUS

157: /*@C
158:   MeshCreateExodus - Create a Mesh from an ExodusII file.

160:   Not Collective

162:   Input Parameters:
163: + comm - The MPI communicator
164: - filename - The ExodusII filename

166:   Output Parameter:
167: . mesh - The Mesh object

169:   Level: beginner

171: .keywords: mesh, ExodusII
172: .seealso: MeshCreate()
173: @*/
174: PetscErrorCode MeshCreateExodus(MPI_Comm comm, const char filename[], Mesh *mesh)
175: {
176:   PetscInt       debug = 0;
177:   PetscTruth     flag;

181:   MeshCreate(comm, mesh);
182:   PetscOptionsGetInt(PETSC_NULL, "-debug", &debug, &flag);
183:   ALE::Obj<PETSC_MESH_TYPE> m = new PETSC_MESH_TYPE(comm, -1, debug);
184: #ifdef PETSC_HAVE_EXODUS
185:   PetscReadExodusII(comm, filename, m);
186: #else
187:   SETERRQ(PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --with-exodus-dir=/path/to/exodus");
188: #endif
189:   if (debug) {m->view("Mesh");}
190:   MeshSetMesh(*mesh, m);
191:   return(0);
192: }

196: /*@C
197:   MeshExodusGetInfo - Get information about an ExodusII Mesh.

199:   Not Collective

201:   Input Parameter:
202: . mesh - The Mesh object

204:   Output Parameters:
205: + dim - The mesh dimension
206: . numVertices - The number of vertices in the mesh
207: . numCells - The number of cells in the mesh
208: . numCellBlocks - The number of cell blocks in the mesh
209: - numVertexSets - The number of vertex sets in the mesh

211:   Level: beginner

213: .keywords: mesh, ExodusII
214: .seealso: MeshCreateExodus()
215: @*/
216: PetscErrorCode MeshExodusGetInfo(Mesh mesh, PetscInt *dim, PetscInt *numVertices, PetscInt *numCells, PetscInt *numCellBlocks, PetscInt *numVertexSets)
217: {
218:   ALE::Obj<PETSC_MESH_TYPE> m;
219:   PetscErrorCode            ierr;

222:   MeshGetMesh(mesh, m);
223:   *dim           = m->getDimension();
224:   *numVertices   = m->depthStratum(0)->size();
225:   *numCells      = m->heightStratum(0)->size();
226:   *numCellBlocks = m->getLabel("CellBlocks")->getCapSize();
227:   *numVertexSets = m->getLabel("VertexSets")->getCapSize();
228:   return(0);
229: }

233: /*@C
234:   MeshGetLabelSize - Get the number of different integer ids in a Label

236:   Not Collective

238:   Input Parameters:
239: + mesh - The Mesh object
240: - name - The label name

242:   Output Parameter:
243: . size - The label size (number of different integer ids)

245:   Level: beginner

247: .keywords: mesh, ExodusII
248: .seealso: MeshCreateExodus()
249: @*/
250: PetscErrorCode MeshGetLabelSize(Mesh mesh, const char name[], PetscInt *size)
251: {
252:   ALE::Obj<PETSC_MESH_TYPE> m;
253:   PetscErrorCode            ierr;

256:   MeshGetMesh(mesh, m);
257:   *size = m->getLabel(name)->getCapSize();
258:   return(0);
259: }

263: /*@C
264:   MeshGetLabelIds - Get the integer ids in a label

266:   Not Collective

268:   Input Parameters:
269: + mesh - The Mesh object
270: . name - The label name
271: - ids - The id storage array

273:   Output Parameter:
274: . ids - The integer ids

276:   Level: beginner

278: .keywords: mesh, ExodusII
279: .seealso: MeshCreateExodus()
280: @*/
281: PetscErrorCode MeshGetLabelIds(Mesh mesh, const char name[], PetscInt *ids)
282: {
283:   ALE::Obj<PETSC_MESH_TYPE> m;
284:   PetscErrorCode            ierr;

287:   MeshGetMesh(mesh, m);
288:   const ALE::Obj<PETSC_MESH_TYPE::label_type::capSequence>&      labelIds = m->getLabel(name)->cap();
289:   const PETSC_MESH_TYPE::label_type::capSequence::const_iterator iEnd     = labelIds->end();
290:   PetscInt                                                       i        = 0;

292:   for(PETSC_MESH_TYPE::label_type::capSequence::const_iterator i_iter = labelIds->begin(); i_iter != iEnd; ++i_iter, ++i) {
293:     ids[i] = *i_iter;
294:   }
295:   return(0);
296: }

300: /*@C
301:   MeshGetStratumSize - Get the number of points in a label stratum

303:   Not Collective

305:   Input Parameters:
306: + mesh - The Mesh object
307: . name - The label name
308: - value - The stratum value

310:   Output Parameter:
311: . size - The stratum size

313:   Level: beginner

315: .keywords: mesh, ExodusII
316: .seealso: MeshCreateExodus()
317: @*/
318: PetscErrorCode MeshGetStratumSize(Mesh mesh, const char name[], PetscInt value, PetscInt *size)
319: {
320:   ALE::Obj<PETSC_MESH_TYPE> m;
321:   PetscErrorCode            ierr;

324:   MeshGetMesh(mesh, m);
325:   *size = m->getLabelStratum(name, value)->size();
326:   return(0);
327: }

331: /*@C
332:   MeshGetStratum - Get the points in a label stratum

334:   Not Collective

336:   Input Parameters:
337: + mesh - The Mesh object
338: . name - The label name
339: . value - The stratum value
340: - points - The stratum points storage array

342:   Output Parameter:
343: . points - The stratum points

345:   Level: beginner

347: .keywords: mesh, ExodusII
348: .seealso: MeshCreateExodus()
349: @*/
350: PetscErrorCode MeshGetStratum(Mesh mesh, const char name[], PetscInt value, PetscInt *points)
351: {
352:   ALE::Obj<PETSC_MESH_TYPE> m;
353:   PetscErrorCode            ierr;

356:   MeshGetMesh(mesh, m);
357:   const ALE::Obj<PETSC_MESH_TYPE::label_sequence>& stratum = m->getLabelStratum(name, value);
358:   const PETSC_MESH_TYPE::label_sequence::iterator  sEnd    = stratum->end();
359:   PetscInt                                         s       = 0;

361:   for(PETSC_MESH_TYPE::label_sequence::iterator s_iter = stratum->begin(); s_iter != sEnd; ++s_iter, ++s) {
362:     points[s] = *s_iter;
363:   }
364:   return(0);
365: }