PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
reader.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_READER_H
12#define RW_READER_H
13
14#include "util/point.h" // definition of Point
15#include "util/matrix.h" // definition of matrices
16#include <cstdint> // uint8_t type
17#include <cstring> // string and size_t type
18#include <vector>
19
26namespace rw {
27
34namespace reader {
35
48void readCsvFile(const std::string &filename, size_t dim,
49 std::vector<util::Point> *nodes,
50 std::vector<double> *volumes);
51
60void readParticleCsvFile(const std::string &filename, size_t dim,
61 std::vector<util::Point> *nodes,
62 std::vector<double> *rads,
63 std::vector<size_t> *zones);
64
73void readParticleCsvFile(const std::string &filename, size_t dim,
74 std::vector<util::Point> *nodes,
75 std::vector<double> *rads,
76 const size_t &zone);
77
87void readParticleWithOrientCsvFile(const std::string &filename, size_t dim,
88 std::vector<util::Point> *nodes,
89 std::vector<double> *rads,
90 std::vector<double> *orients,
91 const size_t &zone);
92
113void readVtuFile(const std::string &filename, size_t dim,
114 std::vector<util::Point> *nodes, size_t &element_type,
115 size_t &num_elem, std::vector<size_t> *enc,
116 std::vector<std::vector<size_t>> *nec,
117 std::vector<double> *volumes, bool is_fd = false);
118
128void readVtuFileNodes(const std::string &filename, size_t dim,
129 std::vector<util::Point> *nodes, bool
130 ref_config = false);
131
141void readVtuFileCells(const std::string &filename, size_t dim,
142 size_t &element_type, size_t &num_elem,
143 std::vector<size_t> *enc,
144 std::vector<std::vector<size_t>> *nec);
145
153void readVtuFileRestart(const std::string &filename,
154 std::vector<util::Point> *u,
155 std::vector<util::Point> *v,
156 const std::vector<util::Point> *X = nullptr);
157
165bool readVtuFilePointData(const std::string &filename,
166 const std::string &tag,
167 std::vector<uint8_t> *data);
168
176bool readVtuFilePointData(const std::string &filename,
177 const std::string &tag,
178 std::vector<size_t> *data);
179
187bool readVtuFilePointData(const std::string &filename,
188 const std::string &tag,
189 std::vector<int> *data);
190
198bool readVtuFilePointData(const std::string &filename,
199 const std::string &tag,
200 std::vector<float> *data);
201
209bool readVtuFilePointData(const std::string &filename,
210 const std::string &tag,
211 std::vector<double> *data);
212
220bool readVtuFilePointData(const std::string &filename,
221 const std::string &tag,
222 std::vector<util::Point> *data);
223
231bool readVtuFilePointData(const std::string &filename,
232 const std::string &tag,
233 std::vector<util::SymMatrix3> *data);
234
242bool readVtuFilePointData(const std::string &filename,
243 const std::string &tag,
244 std::vector<util::Matrix3> *data);
245
253bool readVtuFileCellData(const std::string &filename,
254 const std::string &tag,
255 std::vector<float> *data);
256
264bool readVtuFileCellData(const std::string &filename,
265 const std::string &tag,
266 std::vector<double> *data);
267
275bool readVtuFileCellData(const std::string &filename,
276 const std::string &tag,
277 std::vector<util::Point> *data);
278
286bool readVtuFileCellData(const std::string &filename,
287 const std::string &tag,
288 std::vector<util::SymMatrix3> *data);
289
297bool readVtuFileCellData(const std::string &filename,
298 const std::string &tag,
299 std::vector<util::Matrix3> *data);
300
321void readMshFile(const std::string &filename, size_t dim,
322 std::vector<util::Point> *nodes, size_t &element_type,
323 size_t &num_elem, std::vector<size_t> *enc,
324 std::vector<std::vector<size_t>> *nec,
325 std::vector<double> *volumes, bool is_fd = false);
326
334void readMshFileRestart(const std::string &filename,
335 std::vector<util::Point> *u,
336 std::vector<util::Point> *v,
337 const std::vector<util::Point> *X = nullptr);
338
346bool readMshFilePointData(const std::string &filename,
347 const std::string &tag,
348 std::vector<double> *data);
349
359void readMshFileCells(const std::string &filename, size_t dim,
360 size_t &element_type, size_t &num_elem,
361 std::vector<size_t> *enc,
362 std::vector<std::vector<size_t>> *nec);
363
366} // namespace reader
367
368} // namespace rw
369
370#endif // RW_READER_H
bool readMshFilePointData(const std::string &filename, const std::string &tag, std::vector< double > *data)
Reads data of specified tag from the vtu file.
Definition reader.cpp:390
void readMshFileRestart(const std::string &filename, std::vector< util::Point > *u, std::vector< util::Point > *v, const std::vector< util::Point > *X=nullptr)
Reads mesh data into node file and element file.
Definition reader.cpp:363
void readVtuFileRestart(const std::string &filename, std::vector< util::Point > *u, std::vector< util::Point > *v, const std::vector< util::Point > *X=nullptr)
Reads mesh data into node file and element file.
Definition reader.cpp:181
void readVtuFileCells(const std::string &filename, 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 reader.cpp:168
bool readVtuFileCellData(const std::string &filename, const std::string &tag, std::vector< float > *data)
Reads data of specified tag from the vtu file.
Definition reader.cpp:296
void readParticleWithOrientCsvFile(const std::string &filename, size_t dim, std::vector< util::Point > *nodes, std::vector< double > *rads, std::vector< double > *orients, const size_t &zone)
Reads particles center location, radius, and zone id. In this case, file also provides initial orient...
Definition reader.cpp:101
void readParticleCsvFile(const std::string &filename, size_t dim, std::vector< util::Point > *nodes, std::vector< double > *rads, std::vector< size_t > *zones)
Reads particles center location, radius, and zone id.
Definition reader.cpp:59
void readVtuFile(const std::string &filename, size_t dim, std::vector< util::Point > *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 reader.cpp:125
void readMshFileCells(const std::string &filename, 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 reader.cpp:401
bool readVtuFilePointData(const std::string &filename, const std::string &tag, std::vector< uint8_t > *data)
Reads data of specified tag from the vtu file.
Definition reader.cpp:208
void readVtuFileNodes(const std::string &filename, size_t dim, std::vector< util::Point > *nodes, bool ref_config=false)
Reads nodal coordinates.
Definition reader.cpp:137
void readCsvFile(const std::string &filename, size_t dim, std::vector< util::Point > *nodes, std::vector< double > *volumes)
Reads mesh data into node file and element file.
Definition reader.cpp:16
void readMshFile(const std::string &filename, size_t dim, std::vector< util::Point > *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 reader.cpp:351
Collection of methods and database related to reading and writing.