MercuryDPM  Beta
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
VTKUnstructuredGrid< T > Class Template Reference

#include <VTKData.h>

Public Member Functions

 VTKUnstructuredGrid (std::string filename, const VTKPointDescriptor< T > *descr)
 
 operator bool () const
 
template<typename C >
void write (const C &container)
 

Private Attributes

const VTKPointDescriptor< T > * descriptor_
 
std::ofstream outFile_
 

Detailed Description

template<typename T>
class VTKUnstructuredGrid< T >

This function writes single output frames. this requires an filename and a descriptor. To actually write it, we need a collection.

Definition at line 294 of file VTKData.h.

Constructor & Destructor Documentation

template<typename T >
VTKUnstructuredGrid< T >::VTKUnstructuredGrid ( std::string  filename,
const VTKPointDescriptor< T > *  descr 
)
inline

Create a new VTK Unstructured grid file.

Parameters
filenameThe name of the output file
descrThe typedescriptor for T

Definition at line 303 of file VTKData.h.

304  : descriptor_(descr), outFile_(filename)
305  {
306  }
const VTKPointDescriptor< T > * descriptor_
Definition: VTKData.h:296
std::ofstream outFile_
Definition: VTKData.h:297

Member Function Documentation

template<typename T >
VTKUnstructuredGrid< T >::operator bool ( ) const
inline

Returns true if this VTK output file is valid and has no errorbits set

Definition at line 309 of file VTKData.h.

References VTKUnstructuredGrid< T >::outFile_.

310  {
311  return outFile_.good();
312  }
std::ofstream outFile_
Definition: VTKData.h:297
template<typename T >
template<typename C >
void VTKUnstructuredGrid< T >::write ( const C &  container)
inline

Writes out the container C containing T's, which all get written using the information in the descriptor. C should contain the forward iterable traits, as well as a method size() to give the number of elements stored in this container. All STL containers have this property. Yours should too.

Definition at line 320 of file VTKData.h.

References VTKUnstructuredGrid< T >::descriptor_, Detail::VTKPointDescriptorEntry< T >::emit(), Detail::VTKPointDescriptorEntry< T >::getName(), Detail::VTKPointDescriptorEntry< T >::getNumberOfComponents(), Detail::VTKPointDescriptorEntry< T >::getTypeName(), and VTKUnstructuredGrid< T >::outFile_.

321  {
322  outFile_ <<
323  "<?xml version=\"1.0\"?>\n"
324  "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n"
325  " <UnstructuredGrid>\n"
326  " <Piece NumberOfPoints=\"" << container.size() << "\" NumberOfCells=\"0\">\n"
327  " <Cells>\n"
328  " <DataArray type=\"Int32\" name=\"connectivity\" format=\"ascii\">\n"
329  " 0\n"
330  " </DataArray>\n"
331  " <DataArray type=\"Float32\" name=\"offset\" format=\"ascii\">\n"
332  " 0\n"
333  " </DataArray>\n"
334  " <DataArray type=\"UInt8\" name=\"types\" format=\"ascii\">\n"
335  " 1\n"
336  " </DataArray>\n"
337  " </Cells>\n"
338  " <Points>\n";
339  Detail::VTKPointDescriptorEntry<T>* pDescr = descriptor_->positionEntry_;
340  outFile_ << "<DataArray type=\"" << pDescr->getTypeName() << "\" "
341  "NumberOfComponents=\"" << pDescr->getNumberOfComponents() << "\" "
342  "format=\"ascii\">\n";
343  for (const T& mem : container)
344  {
345  pDescr->emit(outFile_, mem);
346  }
347  outFile_ << "\n</DataArray>\n";
348 
349  outFile_ <<
350  " </Points>\n"
351  " <PointData>\n";
352 
353  for (Detail::VTKPointDescriptorEntry<T>* descr : descriptor_->entries_)
354  {
355  outFile_ << "<DataArray type=\"" << descr->getTypeName() << "\" "
356  "Name=\"" << descr->getName() << "\" "
357  "NumberOfComponents=\"" << descr->getNumberOfComponents() << "\" "
358  "format=\"ascii\">\n";
359  for (const T& mem : container)
360  {
361  descr->emit(outFile_, mem);
362  }
363  outFile_ << "\n</DataArray>\n";
364  }
365  outFile_ <<
366  " </PointData>\n"
367  " <CellData/>\n"
368  " </Piece>\n"
369  " </UnstructuredGrid>\n"
370  "</VTKFile>\n";
371  }
const VTKPointDescriptor< T > * descriptor_
Definition: VTKData.h:296
virtual void emit(std::ostream &out, const T &t) const =0
writes this VTKData to the given output stream for a single T
virtual std::size_t getNumberOfComponents() const =0
Returns the number of components in this type.
std::string getName() const
Returns the name associated with this field.
Definition: VTKData.h:75
std::ofstream outFile_
Definition: VTKData.h:297
virtual std::string getTypeName() const =0
Gives the VTKDataType for VTK.

Member Data Documentation

template<typename T >
const VTKPointDescriptor<T>* VTKUnstructuredGrid< T >::descriptor_
private

Definition at line 296 of file VTKData.h.

Referenced by VTKUnstructuredGrid< T >::write().

template<typename T >
std::ofstream VTKUnstructuredGrid< T >::outFile_
private

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