PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
reader.cpp
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#include "reader.h"
12#include "mshReader.h"
13#include "vtkReader.h"
14#include <csv/csv.h>
15
16void rw::reader::readCsvFile(const std::string &filename, size_t dim,
17 std::vector<util::Point> *nodes,
18 std::vector<double> *volumes) {
19 nodes->clear();
20 volumes->clear();
21 if (dim == 1) {
22 io::CSVReader<3> in(filename);
23 in.read_header(io::ignore_extra_column, "id", "x", "volume");
24
25 double x;
26 double volume;
27 int id;
28 while (in.read_row(id, x, volume)) {
29 volumes->emplace_back(volume);
30 nodes->emplace_back(util::Point(x, 0., 0.));
31 }
32 }
33
34 if (dim == 2) {
35 io::CSVReader<4> in(filename);
36 in.read_header(io::ignore_extra_column, "id", "x", "y", "volume");
37
38 double x, y, volume;
39 int id;
40 while (in.read_row(id, x, y, volume)) {
41 volumes->emplace_back(volume);
42 nodes->emplace_back(util::Point(x, y, 0.));
43 }
44 }
45
46 if (dim == 3) {
47 io::CSVReader<5> in(filename);
48 in.read_header(io::ignore_extra_column, "id", "x", "y", "z", "volume");
49
50 double x, y, z, volume;
51 int id;
52 while (in.read_row(id, x, y, z, volume)) {
53 volumes->emplace_back(volume);
54 nodes->emplace_back(util::Point(x, y, z));
55 }
56 }
57}
58
59void rw::reader::readParticleCsvFile(const std::string &filename, size_t dim,
60 std::vector<util::Point> *nodes,
61 std::vector<double> *rads,
62 std::vector<size_t> *zones) {
63
64 nodes->clear();
65 zones->clear();
66 rads->clear();
67
68 io::CSVReader<5> in(filename);
69 in.read_header(io::ignore_extra_column, "i", "x", "y", "z", "r");
70
71 double x, y, z, r;
72 int id;
73 while (in.read_row(id, x, y, z, r)) {
74 rads->emplace_back(r);
75 nodes->emplace_back(util::Point(x, y, z));
76 zones->emplace_back(id);
77 }
78}
79
80void rw::reader::readParticleCsvFile(const std::string &filename, size_t dim,
81 std::vector<util::Point> *nodes,
82 std::vector<double> *rads,
83 const size_t &zone) {
84
85 nodes->clear();
86 rads->clear();
87
88 io::CSVReader<5> in(filename);
89 in.read_header(io::ignore_extra_column, "i", "x", "y", "z", "r");
90
91 double x, y, z, r;
92 int id;
93 while (in.read_row(id, x, y, z, r)) {
94 if (id == zone) {
95 rads->emplace_back(r);
96 nodes->emplace_back(util::Point(x, y, z));
97 }
98 }
99}
100
101void rw::reader::readParticleWithOrientCsvFile(const std::string &filename, size_t dim,
102 std::vector<util::Point> *nodes,
103 std::vector<double> *rads,
104 std::vector<double> *orients,
105 const size_t &zone) {
106
107 nodes->clear();
108 rads->clear();
109 orients->clear();
110
111 io::CSVReader<6> in(filename);
112 in.read_header(io::ignore_extra_column, "i", "x", "y", "z", "r", "o");
113
114 double x, y, z, r, o;
115 int id;
116 while (in.read_row(id, x, y, z, r, o)) {
117 if (id == zone) {
118 rads->emplace_back(r);
119 nodes->emplace_back(util::Point(x, y, z));
120 orients->emplace_back(o);
121 }
122 }
123}
124
125void rw::reader::readVtuFile(const std::string &filename, size_t dim,
126 std::vector<util::Point> *nodes,
127 size_t &element_type, size_t &num_elem,
128 std::vector<size_t> *enc,
129 std::vector<std::vector<size_t>> *nec,
130 std::vector<double> *volumes, bool is_fd) {
131 // call vtk reader
132 auto rdr = rw::reader::VtkReader(filename);
133 rdr.readMesh(dim, nodes, element_type, num_elem, enc, nec, volumes, is_fd);
134 rdr.close();
135}
136
137void rw::reader::readVtuFileNodes(const std::string &filename, size_t dim,
138 std::vector<util::Point> *nodes,
139 bool ref_config) {
140 // call vtk reader
141 auto rdr = rw::reader::VtkReader(filename);
142
143 // below will read the current position of nodes
144 rdr.readNodes(nodes);
145
146 // need to subtract the displacement to get reference configuration of nodes
147 if (ref_config) {
148 std::vector<util::Point> u;
149 if (rdr.readPointData("Displacement", &u)) {
150 std::cerr << "Error: Did not find displacement in the vtu file."
151 << std::endl;
152 exit(1);
153 }
154
155 if (u.size() != nodes->size()) {
156 std::cerr << "Error: Displacement data and node data size do not match."
157 << std::endl;
158 exit(1);
159 }
160
161 for (size_t i = 0; i < u.size(); i++)
162 (*nodes)[i] -= u[i];
163 }
164
165 rdr.close();
166}
167
168void rw::reader::readVtuFileCells(const std::string &filename, size_t dim,
169 size_t &element_type, size_t &num_elem,
170 std::vector<size_t> *enc,
171 std::vector<std::vector<size_t>> *nec) {
172 // call vtk reader
173 auto rdr = rw::reader::VtkReader(filename);
174
175 // below will read the current position of nodes
176 rdr.readCells(dim, element_type, num_elem, enc, nec);
177
178 rdr.close();
179}
180
181void rw::reader::readVtuFileRestart(const std::string &filename,
182 std::vector<util::Point> *u,
183 std::vector<util::Point> *v,
184 const std::vector<util::Point> *X) {
185 // call vtk reader
186 auto rdr = rw::reader::VtkReader(filename);
187 // if displacement is not in input file, use reference coordinate to get
188 // displacement
189 if (!rdr.readPointData("Displacement", u)) {
190 std::vector<util::Point> y;
191 rdr.readNodes(&y);
192 if (y.size() != X->size()) {
193 std::cerr << "Error: Number of nodes in input file = " << filename
194 << " and number nodes in data X are not same.\n";
195 exit(1);
196 }
197
198 u->resize(y.size());
199 for (size_t i = 0; i < y.size(); i++)
200 (*u)[i] = y[i] - (*X)[i];
201 }
202
203 // get velocity
204 rdr.readPointData("Velocity", v);
205 rdr.close();
206}
207
208bool rw::reader::readVtuFilePointData(const std::string &filename,
209 const std::string &tag,
210 std::vector<uint8_t> *data) {
211 // call vtk reader
212 auto rdr = rw::reader::VtkReader(filename);
213 // read data
214 auto st = rdr.readPointData(tag, data);
215 rdr.close();
216 return st;
217}
218
219bool rw::reader::readVtuFilePointData(const std::string &filename,
220 const std::string &tag,
221 std::vector<size_t> *data) {
222 // call vtk reader
223 auto rdr = rw::reader::VtkReader(filename);
224 // read data
225 auto st = rdr.readPointData(tag, data);
226 rdr.close();
227 return st;
228}
229
230bool rw::reader::readVtuFilePointData(const std::string &filename,
231 const std::string &tag,
232 std::vector<int> *data) {
233 // call vtk reader
234 auto rdr = rw::reader::VtkReader(filename);
235 // read data
236 auto st = rdr.readPointData(tag, data);
237 rdr.close();
238 return st;
239}
240
241bool rw::reader::readVtuFilePointData(const std::string &filename,
242 const std::string &tag,
243 std::vector<float> *data) {
244 // call vtk reader
245 auto rdr = rw::reader::VtkReader(filename);
246 // read data
247 auto st = rdr.readPointData(tag, data);
248 rdr.close();
249 return st;
250}
251
252bool rw::reader::readVtuFilePointData(const std::string &filename,
253 const std::string &tag,
254 std::vector<double> *data) {
255 // call vtk reader
256 auto rdr = rw::reader::VtkReader(filename);
257 // read data
258 auto st = rdr.readPointData(tag, data);
259 rdr.close();
260 return st;
261}
262
263bool rw::reader::readVtuFilePointData(const std::string &filename,
264 const std::string &tag,
265 std::vector<util::Point> *data) {
266 // call vtk reader
267 auto rdr = rw::reader::VtkReader(filename);
268 // read data
269 auto st = rdr.readPointData(tag, data);
270 rdr.close();
271 return st;
272}
273
274bool rw::reader::readVtuFilePointData(const std::string &filename,
275 const std::string &tag,
276 std::vector<util::SymMatrix3> *data) {
277 // call vtk reader
278 auto rdr = rw::reader::VtkReader(filename);
279 // read data
280 auto st = rdr.readPointData(tag, data);
281 rdr.close();
282 return st;
283}
284
285bool rw::reader::readVtuFilePointData(const std::string &filename,
286 const std::string &tag,
287 std::vector<util::Matrix3> *data) {
288 // call vtk reader
289 auto rdr = rw::reader::VtkReader(filename);
290 // read data
291 auto st = rdr.readPointData(tag, data);
292 rdr.close();
293 return st;
294}
295
296bool rw::reader::readVtuFileCellData(const std::string &filename,
297 const std::string &tag,
298 std::vector<float> *data) {
299 // call vtk reader
300 auto rdr = rw::reader::VtkReader(filename);
301 // read data
302 auto st = rdr.readCellData(tag, data);
303 rdr.close();
304 return st;
305}
306
307bool rw::reader::readVtuFileCellData(const std::string &filename,
308 const std::string &tag,
309 std::vector<double> *data) {
310 // call vtk reader
311 auto rdr = rw::reader::VtkReader(filename);
312 // read data
313 auto st = rdr.readCellData(tag, data);
314 rdr.close();
315 return st;
316}
317
318bool rw::reader::readVtuFileCellData(const std::string &filename,
319 const std::string &tag,
320 std::vector<util::Point> *data) {
321 // call vtk reader
322 auto rdr = rw::reader::VtkReader(filename);
323 // read data
324 auto st = rdr.readCellData(tag, data);
325 rdr.close();
326 return st;
327}
328
329bool rw::reader::readVtuFileCellData(const std::string &filename,
330 const std::string &tag,
331 std::vector<util::SymMatrix3> *data) {
332 // call vtk reader
333 auto rdr = rw::reader::VtkReader(filename);
334 // read data
335 auto st = rdr.readCellData(tag, data);
336 rdr.close();
337 return st;
338}
339
340bool rw::reader::readVtuFileCellData(const std::string &filename,
341 const std::string &tag,
342 std::vector<util::Matrix3> *data) {
343 // call vtk reader
344 auto rdr = rw::reader::VtkReader(filename);
345 // read data
346 auto st = rdr.readCellData(tag, data);
347 rdr.close();
348 return st;
349}
350
351void rw::reader::readMshFile(const std::string &filename, size_t dim,
352 std::vector<util::Point> *nodes,
353 size_t &element_type, size_t &num_elem,
354 std::vector<size_t> *enc,
355 std::vector<std::vector<size_t>> *nec,
356 std::vector<double> *volumes, bool is_fd) {
357 // call msh reader
358 auto rdr = rw::reader::MshReader(filename);
359 rdr.readMesh(dim, nodes, element_type, num_elem, enc, nec, volumes, is_fd);
360 rdr.close();
361}
362
363void rw::reader::readMshFileRestart(const std::string &filename,
364 std::vector<util::Point> *u,
365 std::vector<util::Point> *v,
366 const std::vector<util::Point> *X) {
367 // call msh reader
368 auto rdr = rw::reader::MshReader(filename);
369 // if displacement is not in input file, use reference coordinate to get
370 // displacement
371 if (!rdr.readPointData("Displacement", u)) {
372 std::vector<util::Point> y;
373 rdr.readNodes(&y);
374 if (y.size() != X->size()) {
375 std::cerr << "Error: Number of nodes in input file = " << filename
376 << " and number nodes in data X are not same.\n";
377 exit(1);
378 }
379
380 u->resize(y.size());
381 for (size_t i = 0; i < y.size(); i++)
382 (*u)[i] = y[i] - (*X)[i];
383 }
384
385 // get velocity
386 rdr.readPointData("Velocity", v);
387 rdr.close();
388}
389
390bool rw::reader::readMshFilePointData(const std::string &filename,
391 const std::string &tag,
392 std::vector<double> *data) {
393 // call msh reader
394 auto rdr = rw::reader::MshReader(filename);
395 // get velocity
396 auto st = rdr.readPointData(tag, data);
397 rdr.close();
398 return st;
399}
400
401void rw::reader::readMshFileCells(const std::string &filename, size_t dim,
402 size_t &element_type, size_t &num_elem,
403 std::vector<size_t> *enc,
404 std::vector<std::vector<size_t>> *nec) {
405 // call msh reader
406 auto rdr = rw::reader::MshReader(filename);
407 rdr.readCells(dim, element_type, num_elem, enc, nec);
408 rdr.close();
409}
A class to read Gmsh (msh) mesh files.
Definition mshReader.h:28
A class to read VTK (.vtu) mesh files.
Definition vtkReader.h:32
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
A structure to represent 3d vectors.
Definition point.h:30