NLMech  0.1.0
legacyVtkWriter.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_LEGACY_VTKWRITER_H
10 #define RW_LEGACY_VTKWRITER_H
11 
12 #include "util/matrix.h" // definition of SymMatrix3
13 #include "util/point.h" // definition of Point3
14 #include <fstream>
15 #include <vector>
16 
17 namespace rw {
18 
19 namespace writer {
20 
23 
24 public:
34  explicit LegacyVtkWriter(const std::string &filename, const std::string &compress_type = "");
35 
45  void appendNodes(const std::vector<util::Point3> *nodes);
46 
52  void appendNodes(const std::vector<util::Point3> *nodes,
53  const std::vector<util::Point3> *u);
54 
63  void appendMesh(const std::vector<util::Point3> *nodes,
64  const size_t &element_type,
65  const std::vector<size_t> *en_con,
66  const std::vector<util::Point3> *u = nullptr);
67 
80  void appendPointData(const std::string &name,
81  const std::vector<uint8_t> *data);
82 
88  void appendPointData(const std::string &name,
89  const std::vector<size_t> *data);
90 
96  void appendPointData(const std::string &name, const std::vector<int> *data);
97 
103  void appendPointData(const std::string &name, const std::vector<float> *data);
104 
110  void appendPointData(const std::string &name,
111  const std::vector<double> *data);
112 
118  void appendPointData(const std::string &name,
119  const std::vector<util::Point3> *data);
120 
127  void appendPointData(const std::string &name,
128  const std::vector<util::SymMatrix3> *data);
129 
142  void appendCellData(const std::string &name, const std::vector<float> *data);
143 
149  void appendCellData(const std::string &name,
150  const std::vector<util::SymMatrix3> *data);
151 
152 
158  void appendPointData(
159  const std::string &name,
160  const std::vector<blaze::StaticMatrix<double, 3, 3> > *data);
161 
174  void appendFieldData(const std::string &name, const double &data);
175 
181  void appendFieldData(const std::string &name, const float &data);
182 
187  void addTimeStep(const double &timestep);
188 
194  void close();
195 
196 private:
198  std::string d_filename;
199 
201  std::string d_compressType;
202 
204  std::ofstream d_file;
205 
207  std::ofstream d_myfile;
208 };
209 
210 } // namespace writer
211 
212 } // namespace rw
213 
214 #endif // RW_LEGACY_VTKWRITER_H
A vtk writer for simple point data and complex fem mesh data.
Definition: legacyVtkWriter.h:22
std::string d_filename
filename
Definition: legacyVtkWriter.h:198
std::ofstream d_file
vtk/vtu file
Definition: legacyVtkWriter.h:204
void appendFieldData(const std::string &name, const double &data)
Writes the scalar field data to the file.
Definition: legacyVtkWriter.cpp:184
void appendPointData(const std::string &name, const std::vector< uint8_t > *data)
Writes the scalar point data to the file.
Definition: legacyVtkWriter.cpp:87
LegacyVtkWriter(const std::string &filename, const std::string &compress_type="")
Constructor.
Definition: legacyVtkWriter.cpp:17
void addTimeStep(const double &timestep)
Writes the time step to the file.
Definition: legacyVtkWriter.cpp:180
void appendNodes(const std::vector< util::Point3 > *nodes)
Writes the nodes to the file.
Definition: legacyVtkWriter.cpp:35
std::ofstream d_myfile
output stream to write the vtu file
Definition: legacyVtkWriter.h:207
std::string d_compressType
compression_type Specify the compressor (if any)
Definition: legacyVtkWriter.h:201
void close()
Closes the file and store it to the hard disk.
Definition: legacyVtkWriter.cpp:182
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: legacyVtkWriter.cpp:79
void appendCellData(const std::string &name, const std::vector< float > *data)
Writes the float data associated to cells to the file.
Definition: legacyVtkWriter.cpp:166
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