$extrastylesheet
libMesh::VTKIO Class Reference

#include <vtk_io.h>

Inheritance diagram for libMesh::VTKIO:

List of all members.

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

MeshBasemesh ()
void set_n_partitions (unsigned int n_parts)
void skip_comment_lines (std::istream &in, const char comment_start)
const MeshBasemesh () 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.

Author:
Wout Ruijter, 2007 (Checked in to LibMesh by J.W. Peterson)

Definition at line 61 of file vtk_io.h.


Constructor & Destructor Documentation

libMesh::VTKIO::VTKIO ( MeshBase mesh,
MeshData mesh_data = NULL 
) [explicit]

Constructor. Takes a writeable reference to a mesh object. This is the constructor required to read a mesh.

Definition at line 358 of file vtk_io.C.

References _vtk_grid.

                                                 :
  MeshInput<MeshBase> (mesh),
  MeshOutput<MeshBase>(mesh),
  _mesh_data(mesh_data),
  _compress(false),
  _local_node_map()
{
  _vtk_grid = NULL;
  libmesh_experimental();
}
libMesh::VTKIO::VTKIO ( const MeshBase mesh,
MeshData mesh_data = NULL 
) [explicit]

Constructor. Takes a read-only reference to a mesh object. This is the constructor required to write a mesh.

Definition at line 372 of file vtk_io.C.

References _vtk_grid.

                                                       :
  MeshOutput<MeshBase>(mesh),
  _mesh_data(mesh_data),
  _compress(false),
  _local_node_map()
{
  _vtk_grid = NULL;
  libmesh_experimental();
}

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 208 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::DofObject::id(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::n_active_local_elem(), libMesh::Elem::n_nodes(), libMesh::Elem::node(), libMesh::MeshBase::node_ptr(), libMesh::Real, libMesh::Elem::subdomain_id(), libMesh::Elem::type(), and libMesh::VTK.

Referenced by write_nodal_data().

{
  const MeshBase& mesh = MeshOutput<MeshBase>::mesh();

  vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
  vtkSmartPointer<vtkIdList> pts = vtkSmartPointer<vtkIdList>::New();

  std::vector<int> types(mesh.n_active_local_elem());
  unsigned active_element_counter = 0;

  vtkSmartPointer<vtkIntArray> elem_id = vtkSmartPointer<vtkIntArray>::New();
  elem_id->SetName("libmesh_elem_id");
  elem_id->SetNumberOfComponents(1);

  vtkSmartPointer<vtkIntArray> subdomain_id = vtkSmartPointer<vtkIntArray>::New();
  subdomain_id->SetName("subdomain_id");
  subdomain_id->SetNumberOfComponents(1);

  MeshBase::const_element_iterator it = mesh.active_local_elements_begin();
  const MeshBase::const_element_iterator end = mesh.active_local_elements_end();
  for (; it != end; ++it, ++active_element_counter)
    {
      Elem *elem = *it;

      pts->SetNumberOfIds(elem->n_nodes());

      // get the connectivity for this element
      std::vector<dof_id_type> conn;
      elem->connectivity(0, VTK, conn);

      for (unsigned int i=0; i<conn.size(); ++i)
        {
          // If the node ID is not found in the _local_node_map, we'll
          // add it to the _vtk_grid.  NOTE[JWP]: none of the examples
          // I have actually enters this section of code...
          if (_local_node_map.find(conn[i]) == _local_node_map.end())
            {
              dof_id_type global_node_id = elem->node(i);

              const Node* the_node = mesh.node_ptr(global_node_id);

              // Error checking...
              if (the_node == NULL)
                libmesh_error_msg("Error getting pointer to node " << global_node_id << "!");

              // InsertNextPoint accepts either a double or float array of length 3.
              Real pt[3] = {0., 0., 0.};
              for (unsigned int d=0; d<LIBMESH_DIM; ++d)
                pt[d] = (*the_node)(d);

              // Insert the point into the _vtk_grid
              vtkIdType local = _vtk_grid->GetPoints()->InsertNextPoint(pt);

              // Update the _local_node_map with the ID returned by VTK
              _local_node_map[global_node_id] =
                cast_int<dof_id_type>(local);
            }

          // Otherwise, the node ID was found in the _local_node_map, so
          // insert it into the vtkIdList.
          pts->InsertId(i, _local_node_map[conn[i]]);
        }

      vtkIdType vtkcellid = cells->InsertNextCell(pts);
      types[active_element_counter] =
        cast_int<int>(this->get_elem_type(elem->type()));
      elem_id->InsertTuple1(vtkcellid, elem->id());
      subdomain_id->InsertTuple1(vtkcellid, elem->subdomain_id());
    } // end loop over active elements

  _vtk_grid->SetCells(&types[0], cells);
  _vtk_grid->GetCellData()->AddArray(elem_id);
  _vtk_grid->GetCellData()->AddArray(subdomain_id);
}
vtkIdType libMesh::VTKIO::get_elem_type ( ElemType  type) [private]

Map libMesh element types to VTK element types

Definition at line 92 of file vtk_io.C.

References libMesh::EDGE2, libMesh::EDGE3, libMesh::EDGE4, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::INFEDGE2, libMesh::INFHEX16, libMesh::INFHEX18, libMesh::INFHEX8, libMesh::INFPRISM12, libMesh::INFPRISM6, libMesh::INFQUAD4, libMesh::INFQUAD6, libMesh::INVALID_ELEM, libMesh::NODEELEM, libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM6, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::TET10, libMesh::TET4, libMesh::TRI3, libMesh::TRI3SUBDIVISION, and libMesh::TRI6.

Referenced by cells_to_vtk().

{
  vtkIdType celltype = VTK_EMPTY_CELL; // initialize to something to avoid compiler warning

  switch(type)
    {
    case EDGE2:
      celltype = VTK_LINE;
      break;
    case EDGE3:
      celltype = VTK_QUADRATIC_EDGE;
      break;// 1
    case TRI3:
    case TRI3SUBDIVISION:
      celltype = VTK_TRIANGLE;
      break;// 3
    case TRI6:
      celltype = VTK_QUADRATIC_TRIANGLE;
      break;// 4
    case QUAD4:
      celltype = VTK_QUAD;
      break;// 5
    case QUAD8:
      celltype = VTK_QUADRATIC_QUAD;
      break;// 6
    case TET4:
      celltype = VTK_TETRA;
      break;// 8
    case TET10:
      celltype = VTK_QUADRATIC_TETRA;
      break;// 9
    case HEX8:
      celltype = VTK_HEXAHEDRON;
      break;// 10
    case HEX20:
      celltype = VTK_QUADRATIC_HEXAHEDRON;
      break;// 12
    case HEX27:
      celltype = VTK_TRIQUADRATIC_HEXAHEDRON;
      break;
    case PRISM6:
      celltype = VTK_WEDGE;
      break;// 13
    case PRISM15:
      celltype = VTK_QUADRATIC_WEDGE;
      break;// 14
    case PRISM18:
      celltype = VTK_BIQUADRATIC_QUADRATIC_WEDGE;
      break;// 15
    case PYRAMID5:
      celltype = VTK_PYRAMID;
      break;// 16
#if VTK_MAJOR_VERSION > 5 || (VTK_MAJOR_VERSION == 5 && VTK_MINOR_VERSION > 0)
    case QUAD9:
      celltype = VTK_BIQUADRATIC_QUAD;
      break;
#else
    case QUAD9:
#endif
    case EDGE4:
    case INFEDGE2:
    case INFQUAD4:
    case INFQUAD6:
    case INFHEX8:
    case INFHEX16:
    case INFHEX18:
    case INFPRISM6:
    case INFPRISM12:
    case NODEELEM:
    case INVALID_ELEM:
    default:
      libmesh_error_msg("Element type " << type << " not implemented.");
    }
  return celltype;
}
vtkUnstructuredGrid * libMesh::VTKIO::get_vtk_grid ( )

Get a pointer to the VTK datastructure

Definition at line 384 of file vtk_io.C.

References _vtk_grid.

{
  return _vtk_grid;
}
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::ExodusII_IO::copy_elemental_solution(), libMesh::ExodusII_IO::copy_nodal_solution(), libMesh::GMVIO::copy_nodal_solution(), libMesh::TetGenIO::element_in(), libMesh::UNVIO::elements_in(), libMesh::UNVIO::elements_out(), libMesh::UNVIO::groups_in(), libMesh::TetGenIO::node_in(), libMesh::UNVIO::nodes_in(), libMesh::UNVIO::nodes_out(), nodes_to_vtk(), libMesh::AbaqusIO::read(), libMesh::NameBasedIO::read(), libMesh::TetGenIO::read(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), libMesh::GMVIO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::read(), read(), libMesh::LegacyXdrIO::read_ascii(), libMesh::CheckpointIO::read_bcs(), libMesh::CheckpointIO::read_connectivity(), libMesh::AbaqusIO::read_elements(), libMesh::UCDIO::read_implementation(), libMesh::UNVIO::read_implementation(), libMesh::GmshIO::read_mesh(), libMesh::LegacyXdrIO::read_mesh(), libMesh::AbaqusIO::read_nodes(), libMesh::CheckpointIO::read_nodes(), libMesh::CheckpointIO::read_nodesets(), libMesh::XdrIO::read_serialized_bcs(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::OFFIO::read_stream(), libMesh::MatlabIO::read_stream(), libMesh::CheckpointIO::read_subdomain_names(), libMesh::NameBasedIO::write(), libMesh::TetGenIO::write(), libMesh::Nemesis_IO::write(), libMesh::ExodusII_IO::write(), libMesh::CheckpointIO::write(), libMesh::XdrIO::write(), libMesh::GMVIO::write_ascii_new_impl(), libMesh::GMVIO::write_ascii_old_impl(), libMesh::CheckpointIO::write_bcs(), libMesh::GMVIO::write_binary(), libMesh::CheckpointIO::write_connectivity(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), libMesh::UCDIO::write_implementation(), libMesh::GmshIO::write_mesh(), libMesh::LegacyXdrIO::write_mesh(), libMesh::UCDIO::write_nodal_data(), write_nodal_data(), libMesh::Nemesis_IO::write_nodal_data(), libMesh::NameBasedIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_common(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::CheckpointIO::write_nodes(), libMesh::CheckpointIO::write_nodesets(), libMesh::XdrIO::write_parallel(), libMesh::GmshIO::write_post(), libMesh::XdrIO::write_serialized_bcs(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::LegacyXdrIO::write_soln(), and libMesh::CheckpointIO::write_subdomain_names().

void libMesh::VTKIO::nodes_to_vtk ( ) [private]

write the nodes from the mesh into a vtkUnstructuredGrid

Definition at line 170 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().

Referenced by write_nodal_data().

{
  const MeshBase& mesh = MeshOutput<MeshBase>::mesh();

  // containers for points and coordinates of points
  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
  vtkSmartPointer<vtkDoubleArray> pcoords = vtkSmartPointer<vtkDoubleArray>::New();
  pcoords->SetNumberOfComponents(LIBMESH_DIM);
  points->SetNumberOfPoints(mesh.n_local_nodes()); // it seems that it needs this to prevent a segfault

  unsigned int local_node_counter = 0;

  MeshBase::const_node_iterator nd = mesh.local_nodes_begin();
  MeshBase::const_node_iterator nd_end = mesh.local_nodes_end();
  for (; nd != nd_end; nd++, ++local_node_counter)
    {
      Node* node = (*nd);

      double pnt[LIBMESH_DIM];
      for (unsigned int i=0; i<LIBMESH_DIM; ++i)
        pnt[i] = (*node)(i);

      // Fill mapping between global and local node numbers
      _local_node_map[node->id()] = local_node_counter;

      // add point
      pcoords->InsertNextTupleValue(pnt);
    }

  // add coordinates to points
  points->SetData(pcoords);

  // add points to grid
  _vtk_grid->SetPoints(points);
}
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 398 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::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 libMesh::VTK.

Referenced by libMesh::NameBasedIO::read().

{
  // This is a serial-only process for now;
  // the Mesh should be read on processor 0 and
  // broadcast later
  libmesh_assert_equal_to (MeshOutput<MeshBase>::mesh().processor_id(), 0);

  // Keep track of what kinds of elements this file contains
  elems_of_dimension.clear();
  elems_of_dimension.resize(4, false);

#ifndef LIBMESH_HAVE_VTK
  libmesh_error_msg("Cannot read VTK file: " << name \
                    << "\nYou must have VTK installed and correctly configured to read VTK meshes.");

#else
  // Use a typedef, because these names are just crazy
  typedef vtkSmartPointer<vtkXMLUnstructuredGridReader> MyReader;
  MyReader reader = MyReader::New();

  // Pass the filename along to the reader
  reader->SetFileName( name.c_str() );

  // Force reading
  reader->Update();

  // read in the grid
  _vtk_grid = reader->GetOutput();
  // _vtk_grid->Update(); // FIXME: Necessary?

  // Get a reference to the mesh
  MeshBase& mesh = MeshInput<MeshBase>::mesh();

  // Clear out any pre-existing data from the Mesh
  mesh.clear();

  // Get the number of points from the _vtk_grid object
  const unsigned int vtk_num_points = static_cast<unsigned int>(_vtk_grid->GetNumberOfPoints());

  // always numbered nicely??, so we can loop like this
  // I'm pretty sure it is numbered nicely
  for (unsigned int i=0; i<vtk_num_points; ++i)
    {
      // add to the id map
      // and add the actual point
      double * pnt = _vtk_grid->GetPoint(static_cast<vtkIdType>(i));
      Point xyz(pnt[0], pnt[1], pnt[2]);
      Node* newnode = mesh.add_point(xyz, i);

      // Add node to the nodes vector &
      // tell the MeshData object the foreign node id.
      if (this->_mesh_data != NULL)
        this->_mesh_data->add_foreign_node_id (newnode, i);
    }

  // Get the number of cells from the _vtk_grid object
  const unsigned int vtk_num_cells = static_cast<unsigned int>(_vtk_grid->GetNumberOfCells());

  for (unsigned int i=0; i<vtk_num_cells; ++i)
    {
      vtkCell* cell = _vtk_grid->GetCell(i);
      Elem* elem = NULL;
      switch (cell->GetCellType())
        {
        case VTK_LINE:
          elem = new Edge2;
          break;
        case VTK_QUADRATIC_EDGE:
          elem = new Edge3;
          break;
        case VTK_TRIANGLE:
          elem = new Tri3();
          break;
        case VTK_QUADRATIC_TRIANGLE:
          elem = new Tri6();
          break;
        case VTK_QUAD:
          elem = new Quad4();
          break;
        case VTK_QUADRATIC_QUAD:
          elem = new Quad8();
          break;
#if VTK_MAJOR_VERSION > 5 || (VTK_MAJOR_VERSION == 5 && VTK_MINOR_VERSION > 0)
        case VTK_BIQUADRATIC_QUAD:
          elem = new Quad9();
          break;
#endif
        case VTK_TETRA:
          elem = new Tet4();
          break;
        case VTK_QUADRATIC_TETRA:
          elem = new Tet10();
          break;
        case VTK_WEDGE:
          elem = new Prism6();
          break;
        case VTK_QUADRATIC_WEDGE:
          elem = new Prism15();
          break;
        case VTK_BIQUADRATIC_QUADRATIC_WEDGE:
          elem = new Prism18();
          break;
        case VTK_HEXAHEDRON:
          elem = new Hex8();
          break;
        case VTK_QUADRATIC_HEXAHEDRON:
          elem = new Hex20();
          break;
        case VTK_TRIQUADRATIC_HEXAHEDRON:
          elem = new Hex27();
          break;
        case VTK_PYRAMID:
          elem = new Pyramid5();
          break;
        default:
          libmesh_error_msg("Element type not implemented in vtkinterface " << cell->GetCellType());
        }

      // get the straightforward numbering from the VTK cells
      for (unsigned int j=0; j<elem->n_nodes(); ++j)
        elem->set_node(j) =
          mesh.node_ptr(cast_int<dof_id_type>(cell->GetPointId(j)));

      // then get the connectivity
      std::vector<dof_id_type> conn;
      elem->connectivity(0, VTK, conn);

      // then reshuffle the nodes according to the connectivity, this
      // two-time-assign would evade the definition of the vtk_mapping
      for (unsigned int j=0; j<conn.size(); ++j)
        elem->set_node(j) = mesh.node_ptr(conn[j]);

      elem->set_id(i);

      elems_of_dimension[elem->dim()] = true;

      mesh.add_elem(elem);
    } // end loop over VTK cells

  // Set the mesh dimension to the largest encountered for an element
  for (unsigned char i=0; i!=4; ++i)
    if (elems_of_dimension[i])
      mesh.set_mesh_dimension(i);

#if LIBMESH_DIM < 3
  if (mesh.mesh_dimension() > LIBMESH_DIM)
    libmesh_error_msg("Cannot open dimension "  \
                      << mesh.mesh_dimension()              \
                      << " mesh file when configured without "  \
                      << mesh.mesh_dimension()                  \
                      << "D support.");
#endif

#endif // LIBMESH_HAVE_VTK
}

Setter for compression flag

Definition at line 391 of file vtk_io.C.

References _compress.

{
  this->_compress = b;
}
void libMesh::MeshInput< MeshBase >::set_n_partitions ( unsigned int  n_parts) [inline, protected, inherited]

Sets the number of partitions in the mesh. Typically this gets done by the partitioner, but some parallel file formats begin "pre-partitioned".

Definition at line 94 of file mesh_input.h.

References libMesh::MeshInput< MT >::mesh().

Referenced by libMesh::Nemesis_IO::read(), and libMesh::XdrIO::read().

{ this->mesh().set_n_partitions() = n_parts; }
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 289 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::sys, libMesh::System::vectors_begin(), and libMesh::System::vectors_end().

{
  if (MeshOutput<MeshBase>::mesh().processor_id() == 0)
    {
      std::map<std::string, std::vector<Number> > vecs;
      for (unsigned int i=0; i<es.n_systems(); ++i)
        {
          const System& sys = es.get_system(i);
          System::const_vectors_iterator v_end = sys.vectors_end();
          System::const_vectors_iterator it = sys.vectors_begin();
          for (; it!= v_end; ++it)
            {
              // for all vectors on this system
              std::vector<Number> values;
              // libMesh::out<<"it "<<it->first<<std::endl;

              it->second->localize_to_one(values, 0);
              // libMesh::out<<"finish localize"<<std::endl;
              vecs[it->first] = values;
            }
        }

      std::map<std::string, std::vector<Number> >::iterator it = vecs.begin();

      for (; it!=vecs.end(); ++it)
        {
          vtkSmartPointer<vtkDoubleArray> data = vtkSmartPointer<vtkDoubleArray>::New();
          data->SetName(it->first.c_str());
          libmesh_assert_equal_to (it->second.size(), es.get_mesh().n_nodes());
          data->SetNumberOfValues(it->second.size());

          for (unsigned int i=0; i<it->second.size(); ++i)
            {
#ifdef LIBMESH_USE_COMPLEX_NUMBERS
              libmesh_do_once (libMesh::err << "Only writing the real part for complex numbers!\n"
                               << "if you need this support contact " << LIBMESH_PACKAGE_BUGREPORT
                               << std::endl);
              data->SetValue(i, it->second[i].real());
#else
              data->SetValue(i, it->second[i]);
#endif

            }
          grid->GetPointData()->AddArray(data);
        }
    }
}
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 672 of file vtk_io.C.

References write_nodal_data().

Referenced by libMesh::NameBasedIO::write().

{
  std::vector<Number> soln;
  std::vector<std::string> names;
  this->write_nodal_data(name, soln, names);
}
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.

Reimplemented in libMesh::NameBasedIO.

Referenced by libMesh::Nemesis_IO::write_timestep(), and libMesh::ExodusII_IO::write_timestep().

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 >.

Definition at line 556 of file vtk_io.C.

References _compress, _local_node_map, _vtk_grid, cells_to_vtk(), data, libMesh::err, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::n_nodes(), libMesh::n_processors(), nodes_to_vtk(), and libMesh::processor_id().

Referenced by write(), and libMesh::NameBasedIO::write_nodal_data().

{
#ifndef LIBMESH_HAVE_VTK

  libmesh_error_msg("Cannot write VTK file: " << fname                  \
                    << "\nYou must have VTK installed and correctly configured to read VTK meshes.");

#else

  const MeshBase & mesh = MeshOutput<MeshBase>::mesh();

  // Is this really important?  If so, it should be more than an assert...
  // libmesh_assert(fname.substr(fname.rfind("."), fname.size()) == ".pvtu");

  // we only use Unstructured grids
  _vtk_grid = vtkUnstructuredGrid::New();
  vtkSmartPointer<vtkXMLPUnstructuredGridWriter> writer = vtkSmartPointer<vtkXMLPUnstructuredGridWriter>::New();

  // add nodes to the grid and update _local_node_map
  _local_node_map.clear();
  this->nodes_to_vtk();

  // add cells to the grid
  this->cells_to_vtk();

  // add nodal solutions to the grid, if solutions are given
  if (names.size() > 0)
    {
      std::size_t num_vars = names.size();
      dof_id_type num_nodes = mesh.n_nodes();

      for (std::size_t variable=0; variable<num_vars; ++variable)
        {
          vtkSmartPointer<vtkDoubleArray> data = vtkSmartPointer<vtkDoubleArray>::New();
          data->SetName(names[variable].c_str());

          // number of local and ghost nodes
          data->SetNumberOfValues(_local_node_map.size());

          // loop over all nodes and get the solution for the current
          // variable, if the node is in the current partition
          for (dof_id_type k=0; k<num_nodes; ++k)
            {
              if (_local_node_map.find(k) == _local_node_map.end())
                continue; // not a local node

              if (!soln.empty())
                {
#ifdef LIBMESH_USE_COMPLEX_NUMBERS
                  libmesh_do_once (libMesh::err << "Only writing the real part for complex numbers!\n"
                                   << "if you need this support contact " << LIBMESH_PACKAGE_BUGREPORT
                                   << std::endl);
                  data->SetValue(_local_node_map[k], soln[k*num_vars + variable].real());
#else
                  data->SetValue(_local_node_map[k], soln[k*num_vars + variable]);
#endif
                }
              else
                {
                  data->SetValue(_local_node_map[k], 0);
                }
            }
          _vtk_grid->GetPointData()->AddArray(data);
        }
    }

  // Tell the writer how many partitions exist and on which processor
  // we are currently
  writer->SetNumberOfPieces(MeshOutput<MeshBase>::mesh().n_processors());
  writer->SetStartPiece(MeshOutput<MeshBase>::mesh().processor_id());
  writer->SetEndPiece(MeshOutput<MeshBase>::mesh().processor_id());

  // partitions overlap by one node
  // FIXME: According to this document
  // http://paraview.org/Wiki/images/5/51/SC07_tut107_ParaView_Handouts.pdf
  // the ghosts are cells rather than nodes.
  writer->SetGhostLevel(1);

  // Not sure exactly when this changed, but SetInput() is not a
  // method on vtkXMLPUnstructuredGridWriter as of VTK 6.1.0
#if VTK_VERSION_LESS_THAN(6,0,0)
  writer->SetInput(_vtk_grid);
#else
  writer->SetInputData(_vtk_grid);
#endif

  writer->SetFileName(fname.c_str());
  writer->SetDataModeToAscii();

  // compress the output, if desired (switches also to binary)
  if (this->_compress)
    {
#if !VTK_VERSION_LESS_THAN(5,6,0)
      writer->SetCompressorTypeToZLib();
#else
      libmesh_do_once(libMesh::err << "Compression not implemented with old VTK libs!" << std::endl;);
#endif
    }

  writer->Write();

  _vtk_grid->Delete();
#endif
}

Member Data Documentation

bool libMesh::VTKIO::_compress [private]

Flag to indicate whether the output should be compressed

Definition at line 153 of file vtk_io.h.

Referenced by set_compression(), and write_nodal_data().

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::FroIO::write(), libMesh::DivaIO::write(), libMesh::PostscriptIO::write(), and libMesh::EnsightIO::write().

maps global node id to node id of partition

Definition at line 158 of file vtk_io.h.

Referenced by cells_to_vtk(), nodes_to_vtk(), and write_nodal_data().

A pointer to the MeshData object you would like to use. with this VTKIO object. Can be NULL.

Definition at line 148 of file vtk_io.h.

Referenced by read().

vtkUnstructuredGrid* libMesh::VTKIO::_vtk_grid [private]

pointer to the VTK grid

Definition at line 142 of file vtk_io.h.

Referenced by cells_to_vtk(), get_vtk_grid(), nodes_to_vtk(), read(), VTKIO(), and write_nodal_data().


The documentation for this class was generated from the following files: