preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ExportVTK.cpp
Go to the documentation of this file.
1#include "ExportVTK.hpp"
2#include <Eigen/Core>
3#include <filesystem>
4#include <fstream>
5#include <iomanip>
6#include <memory>
7#include "io/Export.hpp"
9#include "mesh/Data.hpp"
10#include "mesh/Edge.hpp"
11#include "mesh/Mesh.hpp"
13#include "mesh/Triangle.hpp"
14#include "mesh/Vertex.hpp"
15#include "utils/IntraComm.hpp"
16#include "utils/assertion.hpp"
17
18namespace precice::io {
19
21 std::string_view participantName,
22 std::string_view location,
23 const mesh::Mesh &mesh,
24 ExportKind kind,
25 int frequency,
26 int rank,
27 int size)
28 : Export(participantName, location, mesh, kind, frequency, rank, size) {};
29
30void ExportVTK::doExport(int index, double time)
31{
34 PRECICE_ASSERT(time >= 0.0);
35 PRECICE_ASSERT(!isParallel(), "ExportVTK only supports serial participants.");
36
37 if (!keepExport(index))
38 return;
39
40 auto filename = fmt::format("{}-{}.{}.vtk", _mesh->getName(), _participantName, formatIndex(index));
41
42 namespace fs = std::filesystem;
43 fs::path outfile(_location);
44 if (not _location.empty())
46 outfile = outfile / filename;
47 std::ofstream outstream(outfile.string(), std::ios::trunc);
48 PRECICE_CHECK(outstream, "VTK export failed to open destination file \"{}\"", outfile.generic_string());
49
50 initializeWriting(outstream);
51 writeHeader(outstream);
52 exportMesh(outstream, *_mesh);
53 exportData(outstream, *_mesh);
54 exportGradient(outstream, *_mesh);
55 outstream.close();
56 recordExport(filename, time);
57}
58
60{
61 if (isParallel())
62 return; // there is no parallel master file
63
64 writeSeriesFile(fmt::format("{}-{}.vtk.series", _mesh->getName(), _participantName));
65}
66
68 std::ofstream &outFile,
69 const mesh::Mesh &mesh)
70{
71 PRECICE_TRACE(mesh.getName());
72
73 // Plot vertices
74 outFile << "POINTS " << mesh.nVertices() << " double \n\n";
75 for (const mesh::Vertex &vertex : mesh.vertices()) {
76 writeVertex(vertex.getCoords(), outFile);
77 }
78 outFile << '\n';
79
80 // Plot edges
81 if (mesh.getDimensions() == 2) {
82 outFile << "CELLS " << mesh.edges().size() << ' ' << mesh.edges().size() * 3 << "\n\n";
83 for (auto const &edge : mesh.edges()) {
84 int internalIndices[2];
85 internalIndices[0] = edge.vertex(0).getID();
86 internalIndices[1] = edge.vertex(1).getID();
87 writeLine(internalIndices, outFile);
88 }
89 outFile << "\nCELL_TYPES " << mesh.edges().size() << "\n\n";
90 for (size_t i = 0; i < mesh.edges().size(); ++i) {
91 outFile << "3\n";
92 }
93 }
94
95 // Plot triangles
96 if (mesh.getDimensions() == 3) {
97 size_t sizeTetrahedra = mesh.tetrahedra().size();
98 size_t sizeTriangles = mesh.triangles().size();
99 size_t sizeEdges = mesh.edges().size();
100 size_t sizeElements = sizeTriangles + sizeEdges + sizeTetrahedra;
101
102 outFile << "CELLS " << sizeElements << ' '
103 << sizeTetrahedra * 5 + sizeTriangles * 4 + sizeEdges * 3 << "\n\n";
104 for (auto const &tetra : mesh.tetrahedra()) {
105 int internalIndices[4];
106 internalIndices[0] = tetra.vertex(0).getID();
107 internalIndices[1] = tetra.vertex(1).getID();
108 internalIndices[2] = tetra.vertex(2).getID();
109 internalIndices[3] = tetra.vertex(3).getID();
110 writeTetrahedron(internalIndices, outFile);
111 }
112 for (auto const &triangle : mesh.triangles()) {
113 int internalIndices[3];
114 internalIndices[0] = triangle.vertex(0).getID();
115 internalIndices[1] = triangle.vertex(1).getID();
116 internalIndices[2] = triangle.vertex(2).getID();
117 writeTriangle(internalIndices, outFile);
118 }
119 for (auto const &edge : mesh.edges()) {
120 int internalIndices[2];
121 internalIndices[0] = edge.vertex(0).getID();
122 internalIndices[1] = edge.vertex(1).getID();
123 writeLine(internalIndices, outFile);
124 }
125
126 outFile << "\nCELL_TYPES " << sizeElements << "\n\n";
127 // See VTK reference for CELL_TYPES
128 for (size_t i = 0; i < sizeTetrahedra; i++) {
129 outFile << "10\n";
130 }
131 for (size_t i = 0; i < sizeTriangles; i++) {
132 outFile << "5\n";
133 }
134 for (size_t i = 0; i < sizeEdges; ++i) {
135 outFile << "3\n";
136 }
137 }
138
139 outFile << '\n';
140}
141
143 std::ofstream &outFile,
144 const mesh::Mesh &mesh)
145{
146 outFile << "POINT_DATA " << mesh.nVertices() << "\n\n";
147
148 outFile << "SCALARS Rank unsigned_int\n";
149 outFile << "LOOKUP_TABLE default\n";
151 outFile << "\n\n";
152
153 for (const mesh::PtrData &data : mesh.data()) { // Plot vertex data
154 if (data->timeStepsStorage().empty()) {
155 continue;
156 }
157 const Eigen::VectorXd &values = data->timeStepsStorage().last().sample.values;
158 if (data->getDimensions() > 1) {
159 Eigen::VectorXd viewTemp(data->getDimensions());
160 outFile << "VECTORS " << data->getName() << " double\n";
161 for (const mesh::Vertex &vertex : mesh.vertices()) {
162 int offset = vertex.getID() * data->getDimensions();
163 for (int i = 0; i < data->getDimensions(); i++) {
164 viewTemp[i] = values(offset + i);
165 }
166 int i = 0;
167 for (; i < data->getDimensions(); i++) {
168 outFile << viewTemp[i] << ' ';
169 }
170 if (i < 3) {
171 outFile << '0';
172 }
173 outFile << '\n';
174 }
175 outFile << '\n';
176 } else if (data->getDimensions() == 1) {
177 outFile << "SCALARS " << data->getName() << " double\n";
178 outFile << "LOOKUP_TABLE default\n";
179 for (const mesh::Vertex &vertex : mesh.vertices()) {
180 outFile << values(vertex.getID()) << '\n';
181 }
182 outFile << '\n';
183 }
184 }
185}
186
188{
189 const int spaceDim = mesh.getDimensions();
190 for (const mesh::PtrData &data : mesh.data()) {
191 if (data->timeStepsStorage().empty() || !data->hasGradient()) {
192 continue;
193 }
194 const auto &gradients = data->timeStepsStorage().last().sample.gradients;
195 if (data->getDimensions() == 1) { // Scalar data, create a vector <dataname>_gradient
196 outFile << "VECTORS " << data->getName() << "_gradient"
197 << " double\n";
198 for (int i = 0; i < gradients.cols(); i++) { // Loop over vertices
199 int j = 0; // Dimension counter
200 for (; j < gradients.rows(); j++) { // Loop over space directions
201 outFile << gradients.coeff(j, i) << " ";
202 }
203 if (j < 3) { // If 2D data add additional zero as third component
204 outFile << '0';
205 }
206 outFile << "\n";
207 }
208 } else { // Vector data, write n vector for n dimension <dataname>_(dx/dy/dz)
209 outFile << "VECTORS " << data->getName() << "_dx"
210 << " double\n";
211 for (int i = 0; i < gradients.cols(); i += spaceDim) { // Loop over vertices
212 int j = 0;
213 for (; j < gradients.rows(); j++) { // Loop over components
214 outFile << gradients.coeff(j, i) << " ";
215 }
216 if (j < 3) { // If 2D data add additional zero as third component
217 outFile << '0';
218 }
219 outFile << "\n";
220 }
221 outFile << "\n";
222
223 outFile << "VECTORS " << data->getName() << "_dy"
224 << " double\n";
225 for (int i = 1; i < gradients.cols(); i += spaceDim) { // Loop over vertices
226 int j = 0;
227 for (; j < gradients.rows(); j++) { // Loop over components
228 outFile << gradients.coeff(j, i) << " ";
229 }
230 if (j < 3) { // If 2D data add additional zero as third component
231 outFile << '0';
232 }
233 outFile << "\n";
234 }
235 outFile << "\n";
236
237 if (spaceDim == 3) { // dz is only for 3D data
238 outFile << "VECTORS " << data->getName() << "_dz"
239 << " double\n";
240 for (int i = 2; i < gradients.cols(); i += spaceDim) { // Loop over vertices
241 int j = 0;
242 for (; j < gradients.rows(); j++) { // Loop over components
243 outFile << gradients.coeff(j, i) << " ";
244 }
245 if (j < 3) { // If 2D data add additional zero as third component
246 outFile << '0';
247 }
248 outFile << "\n";
249 }
250 }
251 }
252 outFile << '\n';
253 }
254}
255
257 std::ofstream &filestream)
258{
259 // size_t pos = fullFilename.rfind(".vtk");
260 // if ((pos == std::string::npos) || (pos != fullFilename.size()-4)){
261 // fullFilename += ".vtk";
262 // }
263 filestream.setf(std::ios::showpoint);
264 filestream.setf(std::ios::scientific);
266}
267
269 std::ostream &outFile)
270{
271 outFile << "# vtk DataFile Version 2.0\n\n"
272 << "ASCII\n\n"
273 << "DATASET UNSTRUCTURED_GRID\n\n";
274}
275
277 const Eigen::VectorXd &position,
278 std::ostream &outFile)
279{
280 if (position.size() == 2) {
281 outFile << position(0) << " " << position(1) << " " << 0.0 << '\n';
282 } else {
283 PRECICE_ASSERT(position.size() == 3);
284 outFile << position(0) << " " << position(1) << " " << position(2) << '\n';
285 }
286}
287
289 int vertexIndices[3],
290 std::ostream &outFile)
291{
292 outFile << 3 << ' ';
293 for (int i = 0; i < 3; i++) {
294 outFile << vertexIndices[i] << ' ';
295 }
296 outFile << '\n';
297}
298
300 int vertexIndices[4],
301 std::ostream &outFile)
302{
303 outFile << 4 << ' ';
304 for (int i = 0; i < 4; i++) {
305 outFile << vertexIndices[i] << ' ';
306 }
307 outFile << '\n';
308}
309
311 int vertexIndices[2],
312 std::ostream &outFile)
313{
314 outFile << 2 << ' ';
315 for (int i = 0; i < 2; i++) {
316 outFile << vertexIndices[i] << ' ';
317 }
318 outFile << '\n';
319}
320
321} // namespace precice::io
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
unsigned int index
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
void exportData(std::ofstream &outFile, const mesh::Mesh &mesh)
ExportVTK(std::string_view participantName, std::string_view location, const mesh::Mesh &mesh, ExportKind kind, int frequency, int rank, int size)
Definition ExportVTK.cpp:20
static void writeVertex(const Eigen::VectorXd &position, std::ostream &outFile)
static void writeHeader(std::ostream &outFile)
static void writeLine(int vertexIndices[2], std::ostream &outFile)
static void initializeWriting(std::ofstream &filestream)
static void writeTetrahedron(int vertexIndices[4], std::ostream &outFile)
void exportSeries() const final override
Definition ExportVTK.cpp:59
void exportMesh(std::ofstream &outFile, const mesh::Mesh &mesh)
Definition ExportVTK.cpp:67
static void writeTriangle(int vertexIndices[3], std::ostream &outFile)
void doExport(int index, double time) final override
Perform writing to VTK file.
Definition ExportVTK.cpp:30
void exportGradient(std::ofstream &outFile, const mesh::Mesh &mesh)
Abstract base class of all classes exporting container data structures.
Definition Export.hpp:14
void writeSeriesFile(std::string_view filename) const
Definition Export.cpp:21
std::string _location
Definition Export.hpp:65
bool isParallel() const
Definition Export.cpp:7
void recordExport(std::string filename, double time)
Definition Export.cpp:44
std::string _participantName
Definition Export.hpp:64
const mesh::Mesh *const _mesh
Definition Export.hpp:66
bool keepExport(int index) const
Definition Export.hpp:59
std::string formatIndex(int index) const
Definition Export.cpp:12
Container and creator for meshes.
Definition Mesh.hpp:38
int getDimensions() const
Definition Mesh.cpp:100
VertexContainer & vertices()
Returns modifieable container holding all vertices.
Definition Mesh.cpp:55
const std::string & getName() const
Returns the name of the mesh, as set in the config file.
Definition Mesh.cpp:220
std::size_t nVertices() const
Returns the number of vertices.
Definition Mesh.cpp:65
TetraContainer & tetrahedra()
Returns modifiable container holding all tetrahedra.
Definition Mesh.cpp:95
const DataContainer & data() const
Allows access to all data.
Definition Mesh.cpp:172
TriangleContainer & triangles()
Returns modifiable container holding all triangles.
Definition Mesh.cpp:80
EdgeContainer & edges()
Returns modifiable container holding all edges.
Definition Mesh.cpp:70
Vertex of a mesh.
Definition Vertex.hpp:16
T close(T... args)
T create_directories(T... args)
T empty(T... args)
T fill_n(T... args)
T generic_string(T... args)
provides Import and Export of the coupling mesh and data.
T setf(T... args)
T setprecision(T... args)
T size(T... args)