NLMech  0.1.0
vtkWriter.h
1 // Copyright (c) 2019 Prashant K. Jha
3 // Copyright (c) 2019 Patrick Diehl
4 //
5 // Distributed under the Boost Software License, Version 1.0. (See accompanying
6 // file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
8 
9 #ifndef RW_VTKWRITER_H
10 #define RW_VTKWRITER_H
11 
12 #include <vtkSmartPointer.h>
13 #include <vtkUnstructuredGrid.h>
14 #include <vtkXMLUnstructuredGridWriter.h>
15 
16 #include "util/point.h" // definition of Point3
17 #include "util/matrix.h" // definition of SymMatrix3
18 
19 namespace rw {
20 
21 namespace writer {
22 
24 class VtkWriter {
25 
26 public:
36  explicit VtkWriter(const std::string &filename, const std::string &compress_type = "");
37 
48  void appendNodes(const std::vector<util::Point3> *nodes,
49  const std::vector<util::Point3> *u = nullptr);
50 
59  void appendMesh(const std::vector<util::Point3> *nodes,
60  const size_t &element_type,
61  const std::vector<size_t> *en_con,
62  const std::vector<util::Point3> *u = nullptr);
63 
76  void appendPointData(const std::string &name,
77  const std::vector<uint8_t> *data);
78 
84  void appendPointData(const std::string &name,
85  const std::vector<size_t> *data);
86 
92  void appendPointData(const std::string &name, const std::vector<int> *data);
93 
99  void appendPointData(const std::string &name, const std::vector<float> *data);
100 
106  void appendPointData(const std::string &name,
107  const std::vector<double> *data);
108 
114  void appendPointData(const std::string &name,
115  const std::vector<util::Point3> *data);
116 
123  void appendPointData(const std::string &name,
124  const std::vector<util::SymMatrix3> *data);
125 
126 
133  void appendPointData(const std::string &name,
134  const std::vector<blaze::StaticMatrix<double, 3, 3> > *data);
135 
148  void appendCellData(const std::string &name, const std::vector<float> *data);
149 
155  void appendCellData(const std::string &name,
156  const std::vector<util::SymMatrix3> *data);
157 
170  void appendFieldData(const std::string &name, const double &data);
171 
177  void appendFieldData(const std::string &name, const float &data);
178 
183  void addTimeStep(const double &timestep);
184 
190  void close();
191 
192 private:
194  vtkSmartPointer<vtkXMLUnstructuredGridWriter> d_writer_p;
195 
197  vtkSmartPointer<vtkUnstructuredGrid> d_grid_p;
198 
200  std::string d_compressType;
201 };
202 
203 } // namespace writer
204 
205 } // namespace rw
206 
207 #endif // RW_VTKWRITER_H
A vtk writer for simple point data and complex fem mesh data.
Definition: vtkWriter.h:24
void appendFieldData(const std::string &name, const double &data)
Writes the scalar field data to the file.
Definition: vtkWriter.cpp:300
std::string d_compressType
compression_type Specify the compressor (if any)
Definition: vtkWriter.h:200
void close()
Closes the file and store it to the hard disk.
Definition: vtkWriter.cpp:289
void appendMesh(const std::vector< util::Point3 > *nodes, const size_t &element_type, const std::vector< size_t > *en_con, const std::vector< util::Point3 > *u=nullptr)
Writes the mesh data to file.
Definition: vtkWriter.cpp:43
void appendNodes(const std::vector< util::Point3 > *nodes, const std::vector< util::Point3 > *u=nullptr)
Writes the nodes to the file.
Definition: vtkWriter.cpp:29
void appendCellData(const std::string &name, const std::vector< float > *data)
Writes the float data associated to cells to the file.
Definition: vtkWriter.cpp:239
vtkSmartPointer< vtkXMLUnstructuredGridWriter > d_writer_p
XML unstructured grid writer.
Definition: vtkWriter.h:194
void addTimeStep(const double &timestep)
Writes the time step to the file.
Definition: vtkWriter.cpp:281
VtkWriter(const std::string &filename, const std::string &compress_type="")
Constructor.
Definition: vtkWriter.cpp:20
vtkSmartPointer< vtkUnstructuredGrid > d_grid_p
Unstructured grid.
Definition: vtkWriter.h:197
void appendPointData(const std::string &name, const std::vector< uint8_t > *data)
Writes the scalar point data to the file.
Definition: vtkWriter.cpp:88
Data mamanger to share the global simulation data between the classes.
Definition: DataManager.h:55
Collection of methods and database related to reading and writing.
Definition: QuasiStaticModel.h:60