PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
|
A vtk writer for simple point data and complex fem mesh data. More...
#include <vtkParticleWriter.h>
Public Member Functions | |
VtkParticleWriter (const std::string &filename, const std::string &compress_type="") | |
Constructor. | |
void | addTimeStep (const double ×tep) |
Writes the time step to the file. | |
void | close () |
Closes the file and store it to the hard disk. | |
Mesh data | |
void | appendNodes (const model::ModelData *model, const std::vector< std::string > &tags) |
Writes the nodes to the file. | |
void | appendMesh (const model::ModelData *model, const std::vector< std::string > &tags) |
Writes the nodes to the file. | |
void | appendContactData (const model::ModelData *model, const std::vector< size_t > *processed_nodes, const std::vector< std::pair< size_t, size_t > > *processed_elems) |
Prepares contact data that is set of nodes in contact and line-element connecting two contacting nodes. | |
void | appendStrainStress (const model::ModelData *model) |
Writes strain/stress. | |
Private Attributes | |
vtkSmartPointer< vtkXMLUnstructuredGridWriter > | d_writer_p |
XML unstructured grid writer. | |
vtkSmartPointer< vtkUnstructuredGrid > | d_grid_p |
Unstructured grid. | |
std::string | d_compressType |
compression_type Specify the compressor (if any) | |
A vtk writer for simple point data and complex fem mesh data.
Definition at line 28 of file vtkParticleWriter.h.
|
explicit |
Constructor.
Creates and opens .vtu file of name given by filename. The file remains open till the close() function is invoked.
filename | Name of file which will be created |
compress_type | Compression method (optional) |
Definition at line 29 of file vtkParticleWriter.cpp.
References d_writer_p.
void rw::writer::VtkParticleWriter::addTimeStep | ( | const double & | timestep | ) |
Writes the time step to the file.
timestep | Current time step of the simulation |
Definition at line 341 of file vtkParticleWriter.cpp.
void rw::writer::VtkParticleWriter::appendContactData | ( | const model::ModelData * | model, |
const std::vector< size_t > * | processed_nodes, | ||
const std::vector< std::pair< size_t, size_t > > * | processed_elems | ||
) |
Prepares contact data that is set of nodes in contact and line-element connecting two contacting nodes.
model | ModelData class object |
processed_nodes | Nodes in the list for which we are writing contact data |
processed_elems | Number of line elements (two nodes) |
Definition at line 361 of file vtkParticleWriter.cpp.
void rw::writer::VtkParticleWriter::appendMesh | ( | const model::ModelData * | model, |
const std::vector< std::string > & | tags | ||
) |
Writes the nodes to the file.
model | ModelData class object |
tags | Vector of tags (name of data, e.g., 'Displacement', 'Velocity') to append to the file |
Definition at line 275 of file vtkParticleWriter.cpp.
References util::vtk_map_element_to_num_nodes.
void rw::writer::VtkParticleWriter::appendNodes | ( | const model::ModelData * | model, |
const std::vector< std::string > & | tags | ||
) |
Writes the nodes to the file.
model | ModelData class object |
tags | Vector of tags (name of data, e.g., 'Displacement', 'Velocity') to append to the file |
Definition at line 39 of file vtkParticleWriter.cpp.
References util::methods::isTagInList().
void rw::writer::VtkParticleWriter::appendStrainStress | ( | const model::ModelData * | model | ) |
Writes strain/stress.
model | ModelData class object |
Definition at line 444 of file vtkParticleWriter.cpp.
void rw::writer::VtkParticleWriter::close | ( | ) |
Closes the file and store it to the hard disk.
Definition at line 350 of file vtkParticleWriter.cpp.
|
private |
compression_type Specify the compressor (if any)
Definition at line 103 of file vtkParticleWriter.h.
|
private |
Unstructured grid.
Definition at line 100 of file vtkParticleWriter.h.
|
private |
XML unstructured grid writer.
Definition at line 97 of file vtkParticleWriter.h.
Referenced by VtkParticleWriter().