PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
legacyVtkWriter.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_LEGACY_VTKWRITER_H
12#define RW_LEGACY_VTKWRITER_H
13
14#include "util/matrix.h" // definition of SymMatrix3
15#include "util/point.h" // definition of Point
16#include <cstdint> // uint8_t type
17#include <cstring> // string and size_t type
18#include <fstream>
19#include <vector>
20
21namespace rw {
22
23namespace writer {
24
27
28public:
38 explicit LegacyVtkWriter(const std::string &filename, const std::string &compress_type = "");
39
49 void appendNodes(const std::vector<util::Point> *nodes);
50
56 void appendNodes(const std::vector<util::Point> *nodes,
57 const std::vector<util::Point> *u);
58
67 void appendMesh(const std::vector<util::Point> *nodes,
68 const size_t &element_type,
69 const std::vector<size_t> *en_con,
70 const std::vector<util::Point> *u = nullptr);
71
84 void appendPointData(const std::string &name,
85 const std::vector<uint8_t> *data);
86
92 void appendPointData(const std::string &name,
93 const std::vector<size_t> *data);
94
100 void appendPointData(const std::string &name, const std::vector<int> *data);
101
107 void appendPointData(const std::string &name, const std::vector<float> *data);
108
114 void appendPointData(const std::string &name,
115 const std::vector<double> *data);
116
122 void appendPointData(const std::string &name,
123 const std::vector<util::Point> *data);
124
131 void appendPointData(const std::string &name,
132 const std::vector<util::SymMatrix3> *data);
133
146 void appendCellData(const std::string &name, const std::vector<float> *data);
147
153 void appendCellData(const std::string &name,
154 const std::vector<util::SymMatrix3> *data);
155
168 void appendFieldData(const std::string &name, const double &data);
169
175 void appendFieldData(const std::string &name, const float &data);
176
181 void addTimeStep(const double &timestep);
182
188 void close();
189
190private:
192 std::string d_filename;
193
195 std::string d_compressType;
196
198 std::ofstream d_file;
199};
200
201} // namespace writer
202
203} // namespace rw
204
205#endif // RW_LEGACY_VTKWRITER_H
A vtk writer for simple point data and complex fem mesh data.
void appendMesh(const std::vector< util::Point > *nodes, const size_t &element_type, const std::vector< size_t > *en_con, const std::vector< util::Point > *u=nullptr)
Writes the mesh data to file.
std::string d_filename
filename
std::ofstream d_file
vtk/vtu file
void appendFieldData(const std::string &name, const double &data)
Writes the scalar field data to the file.
void appendPointData(const std::string &name, const std::vector< uint8_t > *data)
Writes the scalar point data to the file.
void addTimeStep(const double &timestep)
Writes the time step to the file.
void appendNodes(const std::vector< util::Point > *nodes)
Writes the nodes to the file.
std::string d_compressType
compression_type Specify the compressor (if any)
void close()
Closes the file and store it to the hard disk.
void appendCellData(const std::string &name, const std::vector< float > *data)
Writes the float data associated to cells to the file.
Collection of methods and database related to reading and writing.