PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
vtkReader.h
Go to the documentation of this file.
1/*
2 * -------------------------------------------
3 * Copyright (c) 2021 - 2024 Prashant K. Jha
4 * -------------------------------------------
5 * PeriDEM https://github.com/prashjha/PeriDEM
6 *
7 * Distributed under the Boost Software License, Version 1.0. (See accompanying
8 * file LICENSE)
9 */
10
11#ifndef RW_VTKREADER_H
12#define RW_VTKREADER_H
13
14#include <vtkSmartPointer.h>
15#include <vtkUnstructuredGrid.h>
16#include <vtkXMLUnstructuredGridReader.h>
17
18#include "util/point.h" // definition of Point
19#include "util/matrix.h" // definition of matrices
20#include <cstdint> // uint8_t type
21#include <cstring> // string and size_t type
22
23namespace rw {
24
25namespace reader {
26
32class VtkReader {
33
34public:
39 explicit VtkReader(const std::string &filename);
40
53 void readMesh(size_t dim, std::vector<util::Point> *nodes,
54 size_t &element_type, size_t &num_elems,
55 std::vector<size_t> *enc, std::vector<std::vector<size_t>> *nec,
56 std::vector<double> *volumes, bool is_fd = false);
57
63 void readNodes(std::vector<util::Point> *nodes);
64
73 void readCells(size_t dim, size_t &element_type,
74 size_t &num_elems, std::vector<size_t> *enc,
75 std::vector<std::vector<size_t>> *nec);
76
83 bool readPointData(const std::string &name, std::vector<uint8_t> *data);
84
86 bool readPointData(const std::string &name, std::vector<size_t> *data);
87
89 bool readPointData(const std::string &name, std::vector<int> *data);
90
92 bool readPointData(const std::string &name, std::vector<float> *data);
93
95 bool readPointData(const std::string &name, std::vector<double> *data);
96
98 bool readPointData(const std::string &name, std::vector<util::Point> *data);
99
101 bool readPointData(const std::string &name, std::vector<util::SymMatrix3> *data);
102
104 bool readPointData(const std::string &name,
105 std::vector<util::Matrix3> *data);
106
113 bool readCellData(const std::string &name, std::vector<float> *data);
114
116 bool readCellData(const std::string &name, std::vector<double> *data);
117
119 bool readCellData(const std::string &name, std::vector<util::Point> *data);
120
122 bool readCellData(const std::string &name, std::vector<util::SymMatrix3> *data);
123
125 bool readCellData(const std::string &name, std::vector<util::Matrix3> *data);
126
128 void close();
129
130private:
131
133 static size_t d_count;
134
136 vtkSmartPointer<vtkXMLUnstructuredGridReader> d_reader_p;
137
139 vtkSmartPointer<vtkUnstructuredGrid> d_grid_p;
140};
141
142} // namespace reader
143
144} // namespace rw
145
146#endif // RW_VTKREADER_H
A class to read VTK (.vtu) mesh files.
Definition vtkReader.h:32
void close()
Close the file.
static size_t d_count
Counter.
Definition vtkReader.h:133
bool readPointData(const std::string &name, std::vector< uint8_t > *data)
reads point data from .vtu file
void readNodes(std::vector< util::Point > *nodes)
Reads nodal position.
bool readCellData(const std::string &name, std::vector< float > *data)
reads cell data from .vtu file
void readMesh(size_t dim, std::vector< util::Point > *nodes, size_t &element_type, size_t &num_elems, 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:41
void readCells(size_t dim, size_t &element_type, size_t &num_elems, std::vector< size_t > *enc, std::vector< std::vector< size_t > > *nec)
Reads cell data, i.e. element-node connectivity and node-element connectivity.
vtkSmartPointer< vtkXMLUnstructuredGridReader > d_reader_p
XML unstructured grid writer.
Definition vtkReader.h:136
vtkSmartPointer< vtkUnstructuredGrid > d_grid_p
Unstructured grid.
Definition vtkReader.h:139
Collection of methods and database related to reading and writing.