libMesh::VTKIO Class Reference
#include <vtk_io.h>

Public Member Functions | |
| VTKIO (MeshBase &mesh, MeshData *mesh_data=NULL) | |
| VTKIO (const MeshBase &mesh, MeshData *mesh_data=NULL) | |
| virtual void | write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &) |
| virtual void | read (const std::string &) |
| virtual void | write (const std::string &) |
| vtkUnstructuredGrid * | get_vtk_grid () |
| void | set_compression (bool b) |
| virtual void | write_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=NULL) |
| unsigned int & | ascii_precision () |
Protected Member Functions | |
| MeshBase & | mesh () |
| void | skip_comment_lines (std::istream &in, const char comment_start) |
| const MeshBase & | mesh () const |
Protected Attributes | |
| std::vector< bool > | elems_of_dimension |
| const bool | _is_parallel_format |
Private Member Functions | |
| vtkIdType | get_elem_type (ElemType type) |
| void | nodes_to_vtk () |
| void | cells_to_vtk () |
| void | system_vectors_to_vtk (const EquationSystems &es, vtkUnstructuredGrid *&grid) |
Private Attributes | |
| vtkUnstructuredGrid * | _vtk_grid |
| MeshData * | _mesh_data |
| bool | _compress |
| std::map< dof_id_type, dof_id_type > | _local_node_map |
Detailed Description
This class implements reading and writing meshes in the VTK format. Format description: cf. VTK home page.
This class will not have any functionality unless VTK is detected during configure and hence LIBMESH_HAVE_VTK is defined.
Definition at line 61 of file vtk_io.h.
Constructor & Destructor Documentation
Constructor. Takes a writeable reference to a mesh object. This is the constructor required to read a mesh.
Definition at line 170 of file vtk_io.h.
References _vtk_grid.
00170 : 00171 MeshInput<MeshBase> (mesh), 00172 MeshOutput<MeshBase>(mesh), 00173 _mesh_data(mesh_data), 00174 _compress(false), 00175 _local_node_map() 00176 { 00177 _vtk_grid = NULL; 00178 libmesh_experimental(); 00179 }
Constructor. Takes a read-only reference to a mesh object. This is the constructor required to write a mesh.
Definition at line 184 of file vtk_io.h.
References _vtk_grid.
00184 : 00185 MeshOutput<MeshBase>(mesh), 00186 _mesh_data(mesh_data), 00187 _compress(false), 00188 _local_node_map() 00189 { 00190 _vtk_grid = NULL; 00191 libmesh_experimental(); 00192 }
Member Function Documentation
| unsigned int& libMesh::MeshOutput< MeshBase >::ascii_precision | ( | ) | [inherited] |
Return/set the precision to use when writing ASCII files.
By default we use numeric_limits<Real>::digits10 + 2, which should be enough to write out to ASCII and get the exact same Real back when reading in.
Referenced by libMesh::TecplotIO::write_ascii(), libMesh::GMVIO::write_ascii_new_impl(), and libMesh::GMVIO::write_ascii_old_impl().
| void libMesh::VTKIO::cells_to_vtk | ( | ) | [private] |
write the cells from the mesh into a vtkUnstructuredGrid
Definition at line 191 of file vtk_io.C.
References _local_node_map, _vtk_grid, libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::Elem::connectivity(), end, get_elem_type(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::n_active_local_elem(), libMesh::Elem::n_nodes(), libMesh::MeshBase::node(), libMesh::Elem::node(), libMesh::Elem::type(), and libMeshEnums::VTK.
00192 { 00193 const MeshBase& mesh = MeshOutput<MeshBase>::mesh(); 00194 00195 vtkSmartPointer<vtkCellArray> cells = vtkCellArray::New(); 00196 vtkSmartPointer<vtkIdList> pts = vtkIdList::New(); 00197 00198 std::vector<int> types(mesh.n_active_local_elem()); 00199 int active_element_counter = 0; 00200 00201 MeshBase::const_element_iterator it = mesh.active_local_elements_begin(); 00202 const MeshBase::const_element_iterator end = mesh.active_local_elements_end(); 00203 for (; it != end; it++, ++active_element_counter) 00204 { 00205 Elem *elem = *it; 00206 00207 pts->SetNumberOfIds(elem->n_nodes()); 00208 00209 // get the connectivity for this element 00210 std::vector<dof_id_type> conn; 00211 elem->connectivity(0, VTK, conn); 00212 00213 for (unsigned int i = 0; i < conn.size(); ++i) 00214 { 00215 if (_local_node_map.find(conn[i]) == _local_node_map.end()) 00216 { 00217 // Ghost node 00218 float pt[LIBMESH_DIM]; 00219 dof_id_type node = elem->node(i); 00220 for (unsigned int d = 0; d < LIBMESH_DIM; ++d) 00221 { 00222 pt[d] = mesh.node(node)(d); 00223 } 00224 vtkIdType local = _vtk_grid->GetPoints()->InsertNextPoint(pt); 00225 _local_node_map[node] = local; 00226 } 00227 pts->InsertId(i, _local_node_map[conn[i]]); 00228 } 00229 00230 cells->InsertNextCell(pts); 00231 types[active_element_counter] = this->get_elem_type(elem->type()); 00232 } // end loop over active elements 00233 pts->Delete(); 00234 _vtk_grid->SetCells(&types[0], cells); 00235 cells->Delete(); 00236 }
| vtkIdType libMesh::VTKIO::get_elem_type | ( | ElemType | type | ) | [private] |
Map libMesh element types to VTK element types
Definition at line 74 of file vtk_io.C.
References libMeshEnums::EDGE2, libMeshEnums::EDGE3, libMeshEnums::EDGE4, libMesh::err, libMeshEnums::HEX20, libMeshEnums::HEX27, libMeshEnums::HEX8, libMeshEnums::INFEDGE2, libMeshEnums::INFHEX16, libMeshEnums::INFHEX18, libMeshEnums::INFHEX8, libMeshEnums::INFPRISM12, libMeshEnums::INFPRISM6, libMeshEnums::INFQUAD4, libMeshEnums::INFQUAD6, libMeshEnums::INVALID_ELEM, libMeshEnums::NODEELEM, libMeshEnums::PRISM15, libMeshEnums::PRISM18, libMeshEnums::PRISM6, libMeshEnums::PYRAMID5, libMeshEnums::QUAD4, libMeshEnums::QUAD8, libMeshEnums::QUAD9, libMeshEnums::TET10, libMeshEnums::TET4, libMeshEnums::TRI3, and libMeshEnums::TRI6.
Referenced by cells_to_vtk().
00074 { 00075 vtkIdType celltype = VTK_EMPTY_CELL; // initialize to something to avoid compiler warning 00076 00077 switch(type) 00078 { 00079 case EDGE2: 00080 celltype = VTK_LINE; 00081 break; 00082 case EDGE3: 00083 celltype = VTK_QUADRATIC_EDGE; 00084 break;// 1 00085 case TRI3: 00086 celltype = VTK_TRIANGLE; 00087 break;// 3 00088 case TRI6: 00089 celltype = VTK_QUADRATIC_TRIANGLE; 00090 break;// 4 00091 case QUAD4: 00092 celltype = VTK_QUAD; 00093 break;// 5 00094 case QUAD8: 00095 celltype = VTK_QUADRATIC_QUAD; 00096 break;// 6 00097 case TET4: 00098 celltype = VTK_TETRA; 00099 break;// 8 00100 case TET10: 00101 celltype = VTK_QUADRATIC_TETRA; 00102 break;// 9 00103 case HEX8: 00104 celltype = VTK_HEXAHEDRON; 00105 break;// 10 00106 case HEX20: 00107 celltype = VTK_QUADRATIC_HEXAHEDRON; 00108 break;// 12 00109 case HEX27: 00110 celltype = VTK_TRIQUADRATIC_HEXAHEDRON; 00111 break; 00112 case PRISM6: 00113 celltype = VTK_WEDGE; 00114 break;// 13 00115 case PRISM15: 00116 celltype = VTK_QUADRATIC_WEDGE; 00117 break;// 14 00118 case PRISM18: 00119 celltype = VTK_BIQUADRATIC_QUADRATIC_WEDGE; 00120 break;// 15 00121 case PYRAMID5: 00122 celltype = VTK_PYRAMID; 00123 break;// 16 00124 #if VTK_MAJOR_VERSION > 5 || (VTK_MAJOR_VERSION == 5 && VTK_MINOR_VERSION > 0) 00125 case QUAD9: 00126 celltype = VTK_BIQUADRATIC_QUAD; 00127 break; 00128 #else 00129 case QUAD9: 00130 #endif 00131 case EDGE4: 00132 case INFEDGE2: 00133 case INFQUAD4: 00134 case INFQUAD6: 00135 case INFHEX8: 00136 case INFHEX16: 00137 case INFHEX18: 00138 case INFPRISM6: 00139 case INFPRISM12: 00140 case NODEELEM: 00141 case INVALID_ELEM: 00142 default: 00143 { 00144 libMesh::err<<"element type "<<type<<" not implemented"<<std::endl; 00145 libmesh_error(); 00146 break; 00147 } 00148 } 00149 return celltype; 00150 }
| vtkUnstructuredGrid* libMesh::VTKIO::get_vtk_grid | ( | ) | [inline] |
| const MeshBase & libMesh::MeshOutput< MeshBase >::mesh | ( | ) | const [protected, inherited] |
Returns the object as a read-only reference.
Referenced by libMesh::PostscriptIO::write(), libMesh::FroIO::write(), libMesh::EnsightIO::write(), libMesh::DivaIO::write(), libMesh::TecplotIO::write_ascii(), libMesh::MEDITIO::write_ascii(), libMesh::TecplotIO::write_binary(), libMesh::EnsightIO::write_geometry_ascii(), libMesh::EnsightIO::write_scalar_ascii(), libMesh::GnuPlotIO::write_solution(), libMesh::DivaIO::write_stream(), and libMesh::EnsightIO::write_vector_ascii().
| MeshBase & libMesh::MeshInput< MeshBase >::mesh | ( | ) | [protected, inherited] |
Returns the object as a writeable reference.
Referenced by libMesh::GMVIO::_read_materials(), libMesh::GMVIO::_read_nodes(), libMesh::GMVIO::_read_one_cell(), libMesh::AbaqusIO::assign_boundary_node_ids(), libMesh::AbaqusIO::assign_sideset_ids(), libMesh::AbaqusIO::assign_subdomain_ids(), cells_to_vtk(), libMesh::GMVIO::copy_nodal_solution(), libMesh::UNVIO::element_in(), libMesh::TetGenIO::element_in(), libMesh::UNVIO::element_out(), libMesh::UNVIO::node_in(), libMesh::TetGenIO::node_in(), libMesh::UNVIO::node_out(), nodes_to_vtk(), libMesh::XdrIO::read(), read(), libMesh::TetGenIO::read(), libMesh::Nemesis_IO::read(), libMesh::GMVIO::read(), libMesh::ExodusII_IO::read(), libMesh::AbaqusIO::read(), libMesh::LegacyXdrIO::read_ascii(), libMesh::AbaqusIO::read_elements(), libMesh::UNVIO::read_implementation(), libMesh::UCDIO::read_implementation(), libMesh::LegacyXdrIO::read_mesh(), libMesh::GmshIO::read_mesh(), libMesh::AbaqusIO::read_nodes(), libMesh::XdrIO::read_serialized_bcs(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::OFFIO::read_stream(), libMesh::MatlabIO::read_stream(), libMesh::XdrIO::write(), libMesh::TetGenIO::write(), libMesh::Nemesis_IO::write(), libMesh::GMVIO::write_ascii_new_impl(), libMesh::GMVIO::write_ascii_old_impl(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::UNVIO::write_implementation(), libMesh::UCDIO::write_implementation(), libMesh::LegacyXdrIO::write_mesh(), libMesh::GmshIO::write_mesh(), libMesh::UCDIO::write_nodal_data(), libMesh::Nemesis_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::XdrIO::write_parallel(), libMesh::GmshIO::write_post(), libMesh::XdrIO::write_serialized_bcs(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::LegacyXdrIO::write_soln().
| void libMesh::VTKIO::nodes_to_vtk | ( | ) | [private] |
write the nodes from the mesh into a vtkUnstructuredGrid
Definition at line 153 of file vtk_io.C.
References _local_node_map, _vtk_grid, libMesh::DofObject::id(), libMesh::MeshBase::local_nodes_begin(), libMesh::MeshBase::local_nodes_end(), libMesh::MeshInput< MeshBase >::mesh(), and libMesh::MeshBase::n_local_nodes().
00154 { 00155 const MeshBase& mesh = MeshOutput<MeshBase>::mesh(); 00156 00157 // containers for points and coordinates of points 00158 vtkSmartPointer<vtkPoints> points = vtkPoints::New(); 00159 vtkSmartPointer<vtkDoubleArray> pcoords = vtkDoubleArray::New(); 00160 pcoords->SetNumberOfComponents(LIBMESH_DIM); 00161 points->SetNumberOfPoints(mesh.n_local_nodes()); // it seems that it needs this to prevent a segfault 00162 00163 unsigned int local_node_counter = 0; 00164 00165 MeshBase::const_node_iterator nd = mesh.local_nodes_begin(); 00166 MeshBase::const_node_iterator nd_end = mesh.local_nodes_end(); 00167 for (; nd != nd_end; nd++, local_node_counter++) 00168 { 00169 Node* node = (*nd); 00170 00171 float pnt[LIBMESH_DIM]; 00172 for (unsigned int i = 0; i < LIBMESH_DIM; ++i) 00173 pnt[i] = (*node)(i); 00174 00175 // Fill mapping between global and local node numbers 00176 _local_node_map[node->id()] = local_node_counter; 00177 00178 // add point 00179 pcoords->InsertNextTuple(pnt); // SetTuple(node->id(),pnt); 00180 } 00181 00182 // add coordinates to points 00183 points->SetData(pcoords); 00184 pcoords->Delete(); 00185 // add points to grid 00186 _vtk_grid->SetPoints(points); 00187 points->Delete(); 00188 }
| void libMesh::VTKIO::read | ( | const std::string & | name | ) | [virtual] |
Overloads writing equation systems, this is done because when overloading write_nodal_data there would be no way to export cell centered data This method implements reading a mesh from a specified file in VTK format.
Implements libMesh::MeshInput< MeshBase >.
Definition at line 313 of file vtk_io.C.
References _mesh_data, _vtk_grid, libMesh::MeshBase::add_elem(), libMesh::MeshData::add_foreign_node_id(), libMesh::MeshBase::add_point(), libMesh::MeshBase::clear(), libMesh::Elem::connectivity(), libMesh::Elem::dim(), libMesh::MeshInput< MeshBase >::elems_of_dimension, libMesh::err, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::Elem::n_nodes(), libMesh::MeshBase::node_ptr(), libMesh::processor_id(), libMesh::DofObject::set_id(), libMesh::MeshBase::set_mesh_dimension(), libMesh::Elem::set_node(), and libMeshEnums::VTK.
00314 { 00315 // This is a serial-only process for now; 00316 // the Mesh should be read on processor 0 and 00317 // broadcast later 00318 libmesh_assert_equal_to (libMesh::processor_id(), 0); 00319 00320 // Keep track of what kinds of elements this file contains 00321 elems_of_dimension.clear(); 00322 elems_of_dimension.resize(4, false); 00323 00324 #ifndef LIBMESH_HAVE_VTK 00325 libMesh::err << "Cannot read VTK file: " << name 00326 << "\nYou must have VTK installed and correctly configured to read VTK meshes." 00327 << std::endl; 00328 libmesh_error(); 00329 00330 #else 00331 // libMesh::out<<"read "<<name <<std::endl; 00332 vtkXMLUnstructuredGridReader *reader = vtkXMLUnstructuredGridReader::New(); 00333 reader->SetFileName( name.c_str() ); 00334 //libMesh::out<<"force read"<<std::endl; 00335 // Force reading 00336 reader->Update(); 00337 00338 // read in the grid 00339 // vtkUnstructuredGrid *grid = reader->GetOutput(); 00340 _vtk_grid = reader->GetOutput(); 00341 _vtk_grid->Update(); 00342 reader->Delete(); 00343 00344 // Get a reference to the mesh 00345 MeshBase& mesh = MeshInput<MeshBase>::mesh(); 00346 00347 // Clear out any pre-existing data from the Mesh 00348 mesh.clear(); 00349 00350 // always numbered nicely??, so we can loop like this 00351 // I'm pretty sure it is numbered nicely 00352 for (unsigned int i=0; i < static_cast<unsigned int>(_vtk_grid->GetNumberOfPoints()); ++i) 00353 { 00354 // add to the id map 00355 // and add the actual point 00356 double * pnt = _vtk_grid->GetPoint(static_cast<vtkIdType>(i)); 00357 Point xyz(pnt[0],pnt[1],pnt[2]); 00358 Node* newnode = mesh.add_point(xyz,i); 00359 00360 // Add node to the nodes vector & 00361 // tell the MeshData object the foreign node id. 00362 if (this->_mesh_data != NULL) 00363 this->_mesh_data->add_foreign_node_id (newnode, i); 00364 } 00365 00366 for (unsigned int i=0; i < static_cast<unsigned int>(_vtk_grid->GetNumberOfCells()); ++i) 00367 { 00368 vtkCell* cell = _vtk_grid->GetCell(i); 00369 Elem* elem = NULL; // Initialize to avoid compiler warning 00370 switch(cell->GetCellType()) 00371 { 00372 // FIXME - we're not supporting 2D VTK input yet!? [RHS] 00373 case VTK_TETRA: 00374 elem = new Tet4(); 00375 break; 00376 case VTK_WEDGE: 00377 elem = new Prism6(); 00378 break; 00379 case VTK_HEXAHEDRON: 00380 elem = new Hex8(); 00381 break; 00382 case VTK_PYRAMID: 00383 elem = new Pyramid5(); 00384 break; 00385 case VTK_QUADRATIC_HEXAHEDRON: 00386 elem = new Hex20(); 00387 break; 00388 case VTK_QUADRATIC_TETRA: 00389 elem = new Tet10(); 00390 break; 00391 default: 00392 libMesh::err << "element type not implemented in vtkinterface " << cell->GetCellType() << std::endl; 00393 libmesh_error(); 00394 break; 00395 } 00396 00397 // get the straightforward numbering from the VTK cells 00398 for(unsigned int j=0;j<elem->n_nodes();++j){ 00399 elem->set_node(j) = mesh.node_ptr(cell->GetPointId(j)); 00400 } 00401 // then get the connectivity 00402 std::vector<dof_id_type> conn; 00403 elem->connectivity(0,VTK,conn); 00404 // then reshuffle the nodes according to the connectivity, this 00405 // two-time-assign would evade the definition of the vtk_mapping 00406 for(unsigned int j=0;j<conn.size();++j){ 00407 elem->set_node(j) = mesh.node_ptr(conn[j]); 00408 } 00409 elem->set_id(i); 00410 00411 elems_of_dimension[elem->dim()] = true; 00412 00413 mesh.add_elem(elem); 00414 } // end loop over VTK cells 00415 _vtk_grid->Delete(); 00416 00417 // Set the mesh dimension to the largest encountered for an element 00418 for (unsigned int i=0; i!=4; ++i) 00419 if (elems_of_dimension[i]) 00420 mesh.set_mesh_dimension(i); 00421 00422 #if LIBMESH_DIM < 3 00423 if (mesh.mesh_dimension() > LIBMESH_DIM) 00424 { 00425 libMesh::err << "Cannot open dimension " << 00426 mesh.mesh_dimension() << 00427 " mesh file when configured without " << 00428 mesh.mesh_dimension() << "D support." << 00429 std::endl; 00430 libmesh_error(); 00431 } 00432 #endif 00433 00434 #endif // LIBMESH_HAVE_VTK 00435 }
| void libMesh::VTKIO::set_compression | ( | bool | b | ) | [inline] |
| void libMesh::MeshInput< MeshBase >::skip_comment_lines | ( | std::istream & | in, | |
| const char | comment_start | |||
| ) | [protected, inherited] |
Reads input from in, skipping all the lines that start with the character comment_start.
Referenced by libMesh::TetGenIO::read(), and libMesh::UCDIO::read_implementation().
| void libMesh::VTKIO::system_vectors_to_vtk | ( | const EquationSystems & | es, | |
| vtkUnstructuredGrid *& | grid | |||
| ) | [private] |
write the system vectors to vtk
Definition at line 242 of file vtk_io.C.
References data, libMesh::err, libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::MeshBase::n_nodes(), libMesh::EquationSystems::n_systems(), libMesh::processor_id(), libMesh::System::vectors_begin(), and libMesh::System::vectors_end().
00242 { 00243 if (libMesh::processor_id() == 0){ 00244 00245 std::map<std::string,std::vector<Number> > vecs; 00246 for(unsigned int i=0;i<es.n_systems();++i){ 00247 const System& sys = es.get_system(i); 00248 System::const_vectors_iterator v_end = sys.vectors_end(); 00249 System::const_vectors_iterator it = sys.vectors_begin(); 00250 for(;it!= v_end;++it){ // for all vectors on this system 00251 std::vector<Number> values; 00252 // libMesh::out<<"it "<<it->first<<std::endl; 00253 00254 it->second->localize_to_one(values,0); 00255 // libMesh::out<<"finish localize"<<std::endl; 00256 vecs[it->first] = values; 00257 } 00258 } 00259 00260 00261 std::map<std::string,std::vector<Number> >::iterator it = vecs.begin(); 00262 00263 for(;it!=vecs.end();++it){ 00264 00265 vtkDoubleArray *data = vtkDoubleArray::New(); 00266 00267 data->SetName(it->first.c_str()); 00268 00269 libmesh_assert_equal_to (it->second.size(), es.get_mesh().n_nodes()); 00270 00271 data->SetNumberOfValues(it->second.size()); 00272 00273 for(unsigned int i=0;i<it->second.size();++i){ 00274 00275 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 00276 libmesh_do_once (libMesh::err << "Only writing the real part for complex numbers!\n" 00277 << "if you need this support contact " << LIBMESH_PACKAGE_BUGREPORT 00278 << std::endl); 00279 data->SetValue(i,it->second[i].real()); 00280 #else 00281 data->SetValue(i,it->second[i]); 00282 #endif 00283 00284 } 00285 00286 grid->GetPointData()->AddArray(data); 00287 data->Delete(); 00288 } 00289 00290 } 00291 00292 }
| void libMesh::VTKIO::write | ( | const std::string & | name | ) | [virtual] |
Output the mesh without solutions to a .pvtu file
Implements libMesh::MeshOutput< MeshBase >.
Definition at line 544 of file vtk_io.C.
References write_nodal_data().
00544 { 00545 std::vector<Number> soln; 00546 std::vector<std::string> names; 00547 this->write_nodal_data(name, soln, names); 00548 }
| virtual void libMesh::MeshOutput< MeshBase >::write_equation_systems | ( | const std::string & | , | |
| const EquationSystems & | , | |||
| const std::set< std::string > * | system_names = NULL | |||
| ) | [virtual, inherited] |
This method implements writing a mesh with data to a specified file where the data is taken from the EquationSystems object.
Referenced by libMesh::Nemesis_IO::write_timestep().
| virtual void libMesh::VTKIO::write_nodal_data | ( | const std::string & | , | |
| const std::vector< Number > & | , | |||
| const std::vector< std::string > & | ||||
| ) | [virtual] |
This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are provided.
Reimplemented from libMesh::MeshOutput< MeshBase >.
Referenced by write().
Member Data Documentation
bool libMesh::VTKIO::_compress [private] |
Flag to indicate whether the output should be compressed
Definition at line 157 of file vtk_io.h.
Referenced by set_compression().
const bool libMesh::MeshOutput< MeshBase >::_is_parallel_format [protected, inherited] |
Flag specifying whether this format is parallel-capable. If this is false (default) I/O is only permitted when the mesh has been serialized.
Definition at line 126 of file mesh_output.h.
Referenced by libMesh::PostscriptIO::write(), libMesh::FroIO::write(), libMesh::EnsightIO::write(), and libMesh::DivaIO::write().
std::map<dof_id_type, dof_id_type> libMesh::VTKIO::_local_node_map [private] |
maps global node id to node id of partition
Definition at line 162 of file vtk_io.h.
Referenced by cells_to_vtk(), and nodes_to_vtk().
MeshData* libMesh::VTKIO::_mesh_data [private] |
vtkUnstructuredGrid* libMesh::VTKIO::_vtk_grid [private] |
pointer to the VTK grid
Definition at line 146 of file vtk_io.h.
Referenced by cells_to_vtk(), get_vtk_grid(), nodes_to_vtk(), read(), and VTKIO().
std::vector<bool> libMesh::MeshInput< MeshBase >::elems_of_dimension [protected, inherited] |
A vector of bools describing what dimension elements have been encountered when reading a mesh.
Definition at line 93 of file mesh_input.h.
Referenced by libMesh::GMVIO::_read_one_cell(), libMesh::UNVIO::element_in(), read(), libMesh::Nemesis_IO::read(), libMesh::GMVIO::read(), libMesh::ExodusII_IO::read(), libMesh::UNVIO::read_implementation(), libMesh::UCDIO::read_implementation(), libMesh::LegacyXdrIO::read_mesh(), and libMesh::XdrIO::read_serialized_connectivity().
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:43 UTC
Hosted By: