NLMech  0.1.0
vtkReader.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_VTKREADER_H
10 #define RW_VTKREADER_H
11 
12 #include <util/matrixBlaze.h>
13 #include <vtkSmartPointer.h>
14 #include <vtkUnstructuredGrid.h>
15 #include <vtkXMLUnstructuredGridReader.h>
16 
17 #include "util/matrix.h" // definition of matrices
18 #include "util/matrixBlaze.h"
19 #include "util/point.h" // definition of Point3
20 
21 namespace rw {
22 
23 namespace reader {
24 
30 class VtkReader {
31 
32 public:
39  explicit VtkReader(const std::string &filename);
40 
46  bool vtuHasPointData(const std::string &data_tag);
47 
53  bool vtuHasCellData(const std::string &data_tag);
54 
59  std::vector<std::string> readVtuFilePointTags();
60 
65  std::vector<std::string> readVtuFileCellTags();
66 
79  void readMesh(size_t dim, std::vector<util::Point3> *nodes,
80  size_t &element_type, size_t &num_elem,
81  std::vector<size_t> *enc, std::vector<std::vector<size_t>> *nec,
82  std::vector<double> *volumes, bool is_fd = false);
83 
89  void readNodes(std::vector<util::Point3> *nodes);
90 
100  void readCells(size_t dim, size_t &element_type, size_t &num_elem,
101  std::vector<size_t> *enc,
102  std::vector<std::vector<size_t>> *nec);
103 
110  bool readPointData(const std::string &name, std::vector<uint8_t> *data);
111 
118  bool readPointData(const std::string &name, std::vector<size_t> *data);
119 
126  bool readPointData(const std::string &name, std::vector<int> *data);
127 
134  bool readPointData(const std::string &name, std::vector<float> *data);
135 
142  bool readPointData(const std::string &name, std::vector<double> *data);
143 
150  bool readPointData(const std::string &name, std::vector<util::Point3> *data);
151 
158  bool readPointData(const std::string &name,
159  std::vector<util::SymMatrix3> *data);
160 
167  bool readPointData(const std::string &name,
168  std::vector<util::Matrix33> *data);
169 
176  bool readCellData(const std::string &name, std::vector<float> *data);
177 
184  bool readCellData(const std::string &name, std::vector<double> *data);
185 
192  bool readCellData(const std::string &name, std::vector<util::Point3> *data);
193 
200  bool readCellData(const std::string &name,
201  std::vector<util::SymMatrix3> *data);
202 
209  bool readCellData(const std::string &name, std::vector<util::Matrix33> *data);
210 
212  void close();
213 
214 private:
216  static size_t d_count;
217 
219  vtkSmartPointer<vtkXMLUnstructuredGridReader> d_reader_p;
220 
222  vtkSmartPointer<vtkUnstructuredGrid> d_grid_p;
223 };
224 
225 } // namespace reader
226 
227 } // namespace rw
228 
229 #endif // RW_VTKREADER_H
A class to read VTK (.vtu) mesh files.
Definition: vtkReader.h:30
void close()
Close the file.
Definition: vtkReader.cpp:692
static size_t d_count
Counter.
Definition: vtkReader.h:216
void readNodes(std::vector< util::Point3 > *nodes)
Reads nodal position.
Definition: vtkReader.cpp:234
void readCells(size_t dim, size_t &element_type, size_t &num_elem, std::vector< size_t > *enc, std::vector< std::vector< size_t >> *nec)
Reads cell data, i.e. element-node connectivity and node-element connectivity.
Definition: vtkReader.cpp:252
std::vector< std::string > readVtuFileCellTags()
Reads all cell data tags.
Definition: vtkReader.cpp:88
bool vtuHasPointData(const std::string &data_tag)
Checks if file has needed data.
Definition: vtkReader.cpp:48
bool readPointData(const std::string &name, std::vector< uint8_t > *data)
reads point data from .vtu file
Definition: vtkReader.cpp:318
VtkReader(const std::string &filename)
Constructor.
Definition: vtkReader.cpp:29
bool readCellData(const std::string &name, std::vector< float > *data)
reads cell data from .vtu file
Definition: vtkReader.cpp:544
void readMesh(size_t dim, std::vector< util::Point3 > *nodes, size_t &element_type, size_t &num_elem, std::vector< size_t > *enc, std::vector< std::vector< size_t >> *nec, std::vector< double > *volumes, bool is_fd=false)
Reads mesh data into node file and element file.
Definition: vtkReader.cpp:103
vtkSmartPointer< vtkXMLUnstructuredGridReader > d_reader_p
XML unstructured grid writer.
Definition: vtkReader.h:219
std::vector< std::string > readVtuFilePointTags()
Reads all point data tags.
Definition: vtkReader.cpp:73
vtkSmartPointer< vtkUnstructuredGrid > d_grid_p
Unstructured grid.
Definition: vtkReader.h:222
bool vtuHasCellData(const std::string &data_tag)
Checks if file has needed data.
Definition: vtkReader.cpp:61
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