Actual source code: meshexodus.c

  1: #include<petscdmmesh_formats.hh>   /*I "petscdmmesh.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:       typedef ALE::Mesh<PetscInt,PetscScalar> FlexMesh;
124:       ALE::Obj<FlexMesh::sieve_type> s = new FlexMesh::sieve_type(sieve->comm(), sieve->debug());

126:       ALE::SieveBuilder<FlexMesh>::buildTopology(s, num_dim, num_elem, cells, num_nodes, interpolate, numCorners);
127:       std::map<FlexMesh::point_type,FlexMesh::point_type> renumbering;
128:       ALE::ISieveConverter::convertSieve(*s, *sieve, renumbering);
129:     }
130:     if (!interpolate) {
131:       // Optimized stratification
132:       const ALE::Obj<PETSC_MESH_TYPE::label_type>& height = mesh->createLabel("height");
133:       const ALE::Obj<PETSC_MESH_TYPE::label_type>& depth  = mesh->createLabel("depth");

135:       for(int c = 0; c < num_elem; ++c) {
136:         height->setCone(0, c);
137:         depth->setCone(1, c);
138:       }
139:       for(int v = num_elem; v < num_elem+num_nodes; ++v) {
140:         height->setCone(1, v);
141:         depth->setCone(0, v);
142:       }
143:       mesh->setHeight(1);
144:       mesh->setDepth(1);
145:     } else {
146:       mesh->stratify();
147:     }
148:   } else {
149:     mesh->getSieve()->setChart(PETSC_MESH_TYPE::sieve_type::chart_type());
150:     mesh->getSieve()->allocate();
151:     mesh->stratify();
152:   }
153:   } catch (ALE::Exception e) {
154:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB, e.msg().c_str());
155:   }
156:   PetscFree(cells);

158:   // Build cell blocks
159:   const ALE::Obj<PETSC_MESH_TYPE::label_type>& cellBlocks = mesh->createLabel("CellBlocks");
160:   if (rank == 0) {
161:     for(int eb = 0, k = 0; eb < num_elem_blk; ++eb) {
162:       for(int e = 0; e < num_elem_in_block[eb]; ++e, ++k) {
163:         mesh->setValue(cellBlocks, k, eb_ids[eb]);
164:       }
165:     }
166:   }
167:   if (num_elem_blk > 0) {
168:     PetscFree(connect);
169:     PetscFree5(eb_ids, num_elem_in_block, num_nodes_per_elem, num_attr, block_names);
170:   }

172:   // Build coordinates
173:   double *coords;
174:   PetscMalloc(num_dim*num_nodes * sizeof(double), &coords);
175:   if (num_dim > 0) {for(int v = 0; v < num_nodes; ++v) {coords[v*num_dim+0] = x[v];}}
176:   if (num_dim > 1) {for(int v = 0; v < num_nodes; ++v) {coords[v*num_dim+1] = y[v];}}
177:   if (num_dim > 2) {for(int v = 0; v < num_nodes; ++v) {coords[v*num_dim+2] = z[v];}}
178:   ALE::SieveBuilder<PETSC_MESH_TYPE>::buildCoordinates(mesh, num_dim, coords);
179:   PetscFree(coords);
180:   PetscFree3(x, y, z);

182:   // Build vertex sets
183:   const ALE::Obj<PETSC_MESH_TYPE::label_type>& vertexSets = mesh->createLabel("VertexSets");
184:   if (rank == 0) {
185:     for(int ns = 0; ns < num_node_sets; ++ns) {
186:       for(int n = 0; n < num_nodes_in_set[ns]; ++n) {
187:         mesh->addValue(vertexSets, node_list[ns][n]+num_elem-1, ns_ids[ns]);
188:       }
189:       PetscFree(node_list[ns]);
190:     }
191:   }
192:   if (num_node_sets > 0) {
193:     PetscFree3(ns_ids,num_nodes_in_set,node_list);
194:   }

196:   //cellBlocks->view("Cell Blocks");
197:   //vertexSets->view("Vertex Sets");

199:   // Get coords and print in F90
200:   // Get connectivity and print in F90
201:   // Calculate cost function
202:   // Do in parallel
203:   //   Read in parallel
204:   //   Distribute
205:   //   Print out local meshes
206:   //   Do Blaise's assembly loop in parallel
207:   // Assemble function into Section
208:   // Assemble jacobian into Mat
209:   // Assemble in parallel
210:   // Convert Section to Vec
211:   return(0);
212: }

214: #endif // PETSC_HAVE_EXODUSII

218: /*@C
219:   DMMeshCreateExodus - Create a DMMesh from an ExodusII file.

221:   Not Collective

223:   Input Parameters:
224: + comm - The MPI communicator
225: - filename - The ExodusII filename

227:   Output Parameter:
228: . dm - The DMMesh object

230:   Level: beginner

232: .keywords: mesh, ExodusII
233: .seealso: DMMeshCreate()
234: @*/
235: PetscErrorCode DMMeshCreateExodus(MPI_Comm comm, const char filename[], DM *dm)
236: {
237:   PetscInt       debug = 0;
238:   PetscBool      flag;

242:   DMMeshCreate(comm, dm);
243:   PetscOptionsGetInt(PETSC_NULL, "-debug", &debug, &flag);
244:   ALE::Obj<PETSC_MESH_TYPE> m = new PETSC_MESH_TYPE(comm, -1, debug);
245: #ifdef PETSC_HAVE_EXODUSII
246:   PetscReadExodusII(comm, filename, m);
247: #else
248:   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --with-exodus-dir=/path/to/exodus");
249: #endif
250:   if (debug) {m->view("Mesh");}
251:   DMMeshSetMesh(*dm, m);
252:   return(0);
253: }

257: /*@
258:   DMMeshExodusGetInfo - Get information about an ExodusII Mesh.

260:   Not Collective

262:   Input Parameter:
263: . dm - The DMMesh object

265:   Output Parameters:
266: + dim - The mesh dimension
267: . numVertices - The number of vertices in the mesh
268: . numCells - The number of cells in the mesh
269: . numCellBlocks - The number of cell blocks in the mesh
270: - numVertexSets - The number of vertex sets in the mesh

272:   Level: beginner

274: .keywords: mesh, ExodusII
275: .seealso: DMMeshCreateExodus()
276: @*/
277: PetscErrorCode DMMeshExodusGetInfo(DM dm, PetscInt *dim, PetscInt *numVertices, PetscInt *numCells, PetscInt *numCellBlocks, PetscInt *numVertexSets)
278: {
279:   ALE::Obj<PETSC_MESH_TYPE> m;
280:   PetscErrorCode            ierr;

283:   DMMeshGetMesh(dm, m);
284:   *dim           = m->getDimension();
285:   *numVertices   = m->depthStratum(0)->size();
286:   *numCells      = m->heightStratum(0)->size();
287:   *numCellBlocks = m->getLabel("CellBlocks")->getCapSize();
288:   *numVertexSets = m->getLabel("VertexSets")->getCapSize();
289:   return(0);
290: }