PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
mshWriter.h
Go to the documentation of this file.
1// Copyright (c) 2021 Prashant K. Jha
2//
3// Distributed under the GNU GENERAL PUBLIC LICENSE, Version 3.0.
4// (See accompanying file LICENSE.txt)
5
6#ifndef RW_MSHWRITER_H
7#define RW_MSHWRITER_H
8
9#include "util/matrix.h" // definition of SymMatrix3
10#include "util/point.h" // definition of Point
11#include <cstdint> // uint8_t type
12#include <cstring> // string and size_t type
13#include <fstream>
14#include <vector>
15
16namespace rw {
17
18namespace writer {
19
24class MshWriter {
25
26public:
35 explicit MshWriter(const std::string &filename, const std::string &compress_type = "");
36
47 void appendNodes(const std::vector<util::Point> *nodes,
48 const std::vector<util::Point> *u = nullptr);
49
58 void appendMesh(const std::vector<util::Point> *nodes,
59 const size_t &element_type,
60 const std::vector<size_t> *en_con,
61 const std::vector<util::Point> *u = nullptr);
62
75 void appendPointData(const std::string &name,
76 const std::vector<uint8_t> *data);
77
83 void appendPointData(const std::string &name,
84 const std::vector<size_t> *data);
85
91 void appendPointData(const std::string &name, const std::vector<int> *data);
92
98 void appendPointData(const std::string &name, const std::vector<float> *data);
99
105 void appendPointData(const std::string &name,
106 const std::vector<double> *data);
107
113 void appendPointData(const std::string &name,
114 const std::vector<util::Point> *data);
115
122 void appendPointData(const std::string &name,
123 const std::vector<util::SymMatrix3> *data);
124
137 void appendCellData(const std::string &name, const std::vector<float> *data);
138
144 void appendCellData(const std::string &name,
145 const std::vector<util::SymMatrix3> *data);
146
159 void appendFieldData(const std::string &name, const double &data);
160
166 void appendFieldData(const std::string &name, const float &data);
167
172 void addTimeStep(const double &timestep);
173
179 void close();
180
181private:
195 void writeMshDataHeader(const std::string &name, int field_type,
196 size_t num_data, bool is_node_data = true);
197
199 std::string d_filename;
200
202 std::string d_compressType;
203
205 FILE *d_file;
206};
207
208} // namespace writer
209
210} // namespace rw
211
212#endif // RW_MSHWRITER_H
A .msh writer for simple point data and complex fem mesh data.
Definition mshWriter.h:24
std::string d_compressType
compression_type Specify the compressor (if any)
Definition mshWriter.h:202
void addTimeStep(const double &timestep)
Writes the time step to the file.
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.
Definition mshWriter.cpp:87
void writeMshDataHeader(const std::string &name, int field_type, size_t num_data, bool is_node_data=true)
utility function
Definition mshWriter.cpp:22
FILE * d_file
msh file
Definition mshWriter.h:205
void appendFieldData(const std::string &name, const double &data)
Writes the scalar field data to the file.
void close()
Closes the file and store it to the hard disk.
void appendNodes(const std::vector< util::Point > *nodes, const std::vector< util::Point > *u=nullptr)
Writes the nodes to the file.
Definition mshWriter.cpp:54
void appendPointData(const std::string &name, const std::vector< uint8_t > *data)
Writes the scalar point data to the file.
void appendCellData(const std::string &name, const std::vector< float > *data)
Writes the float data associated to cells to the file.
std::string d_filename
filename
Definition mshWriter.h:199
Collection of methods and database related to reading and writing.