NLMech  0.1.0
writer.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_WRITER_H
10 #define RW_WRITER_H
11 
12 #include "util/point.h" // definition of Point3
13 #include "util/matrix.h" // definition of SymMatrix3
14 #include <string>
15 #include <vector>
16 
17 // forward declaration
18 namespace rw {
19 namespace writer {
20 class VtkWriter;
21 class LegacyVtkWriter;
22 class MshWriter;
23 }
24 } // namespace rw
25 
26 namespace rw {
27 
34 namespace writer {
35 
41 class Writer {
42 
43 public:
47  Writer();
48 
60  explicit Writer(const std::string &filename, const std::string &format =
61  "vtu", const std::string &compress_type = "");
62 
64  ~Writer();
65 
73  void open(const std::string &filename, const std::string &format = "vtu",
74  const std::string &compress_type = "");
75 
86  void appendNodes(const std::vector<util::Point3> *nodes,
87  const std::vector<util::Point3> *u = nullptr);
88 
97  void appendMesh(const std::vector<util::Point3> *nodes,
98  const size_t &element_type,
99  const std::vector<size_t> *en_con,
100  const std::vector<util::Point3> *u = nullptr);
101 
114  void appendPointData(const std::string &name,
115  const std::vector<uint8_t> *data);
116 
122  void appendPointData(const std::string &name,
123  const std::vector<size_t> *data);
124 
130  void appendPointData(const std::string &name, const std::vector<int> *data);
131 
137  void appendPointData(const std::string &name, const std::vector<float> *data);
138 
144  void appendPointData(const std::string &name,
145  const std::vector<double> *data);
146 
152  void appendPointData(const std::string &name,
153  const std::vector<util::Point3> *data);
154 
161  void appendPointData(const std::string &name,
162  const std::vector<util::SymMatrix3> *data);
163 
170  void appendPointData(const std::string &name,
171  const std::vector<blaze::StaticMatrix<double, 3, 3> > *data);
172 
185  void appendCellData(const std::string &name, const std::vector<float> *data);
186 
192  void appendCellData(const std::string &name,
193  const std::vector<util::SymMatrix3> *data);
194 
207  void appendFieldData(const std::string &name, const double &data);
208 
214  void appendFieldData(const std::string &name, const float &data);
215 
220  void addTimeStep(const double &timestep);
221 
227  void close();
228 
229 private:
232 
235 
238 
240  std::string d_format;
241 
248  inline void checkLength(const size_t length, const std::string &name);
249 
250 }; // class Writer
251 
252 } // namespace writer
253 
254 } // namespace rw
255 
256 #endif // RW_WRITER_H
A vtk writer for simple point data and complex fem mesh data.
Definition: legacyVtkWriter.h:22
A .msh writer for simple point data and complex fem mesh data.
Definition: mshWriter.h:22
A vtk writer for simple point data and complex fem mesh data.
Definition: vtkWriter.h:24
A interface class writing data.
Definition: writer.h:41
void checkLength(const size_t length, const std::string &name)
Check if the field is not zero This check is important for the legacy vtk writer, since empty fields ...
Definition: writer.cpp:219
void appendCellData(const std::string &name, const std::vector< float > *data)
Writes the float data associated to cells to the file.
Definition: writer.cpp:159
rw::writer::MshWriter * d_mshWriter_p
Pointer to the vtk writer class.
Definition: writer.h:237
void appendFieldData(const std::string &name, const double &data)
Writes the scalar field data to the file.
Definition: writer.cpp:190
void addTimeStep(const double &timestep)
Writes the time step to the file.
Definition: writer.cpp:181
void appendPointData(const std::string &name, const std::vector< uint8_t > *data)
Writes the scalar point data to the file.
Definition: writer.cpp:67
std::string d_format
Format of output file.
Definition: writer.h:240
~Writer()
Destructor.
Definition: writer.cpp:43
void open(const std::string &filename, const std::string &format="vtu", const std::string &compress_type="")
Open a .vtu file.
Definition: writer.cpp:30
Writer()
Constructor.
Definition: writer.cpp:14
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: writer.cpp:55
void close()
Closes the file and store it to the hard disk.
Definition: writer.cpp:210
rw::writer::VtkWriter * d_vtkWriter_p
Pointer to the vtk writer class.
Definition: writer.h:231
void appendNodes(const std::vector< util::Point3 > *nodes, const std::vector< util::Point3 > *u=nullptr)
Writes the nodes to the file.
Definition: writer.cpp:45
rw::writer::LegacyVtkWriter * d_legacyVtkWriter_p
Pointer to the vtk writer class.
Definition: writer.h:234
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