PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
vtkWriter.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 "vtkWriter.h"
12#include <util/feElementDefs.h>
13#include <vtkCellArray.h>
14#include <vtkCellData.h>
15#include <vtkDoubleArray.h>
16#include <vtkIdList.h>
17#include <vtkIntArray.h>
18#include <vtkPointData.h>
19#include <vtkPoints.h>
20#include <vtkUnsignedIntArray.h>
21
22rw::writer::VtkWriter::VtkWriter(const std::string &filename,
23 const std::string &compress_type)
24 : d_compressType(compress_type) {
25
26 std::string f = filename + ".vtu";
27
28 d_writer_p = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
29 d_writer_p->SetFileName(const_cast<char *>(f.c_str()));
30}
31
32void rw::writer::VtkWriter::appendNodes(const std::vector<util::Point> *nodes,
33 const std::vector<util::Point> *u) {
34
35 auto points = vtkSmartPointer<vtkPoints>::New();
36
37 for (size_t i = 0; i < nodes->size(); i++) {
38
39 util::Point p = (*nodes)[i];
40 if (u)
41 p = p + (*u)[i];
42 points->InsertNextPoint(p.d_x, p.d_y, p.d_z);
43 }
44
45 d_grid_p = vtkSmartPointer<vtkUnstructuredGrid>::New();
46 d_grid_p->SetPoints(points);
47}
48
50 const std::vector<util::Point> *nodes, const size_t &element_type,
51 const std::vector<size_t> *en_con,
52 const std::vector<util::Point> *u) {
53
54 // we write following things to the file
55 //
56 // Node data
57 // 1. Coordinates of nodes (current)
58 //
59 // Element data
60 // 1. element node connectivity
61 // 2. element type (either triangle or square)
62
63 // add current position of nodes
64 this->appendNodes(nodes, u);
65
66 // get the total number of elements
67 size_t num_vertex = util::vtk_map_element_to_num_nodes[element_type];
68 size_t num_elems = en_con->size() / num_vertex;
69
70 //
71 // process elements data
72 //
73 // element node connectivity
74 auto cells = vtkSmartPointer<vtkCellArray>::New();
75 cells->Allocate(num_vertex, num_elems);
76
77 // element type
78 int cell_types[num_elems];
79
80 vtkIdType ids[num_vertex];
81 for (size_t i = 0; i < num_elems; i++) {
82
83 // get ids of vertex of this element
84 for (size_t k = 0; k < num_vertex; k++)
85 ids[k] = (*en_con)[num_vertex*i + k];
86
87 cells->InsertNextCell(num_vertex, ids);
88 cell_types[i] = element_type;
89 }
90
91 // element node connectivity
92 d_grid_p->SetCells(cell_types, cells);
93}
94
95void rw::writer::VtkWriter::appendPointData(const std::string &name,
96 const std::vector<uint8_t> *data) {
97
98 auto array = vtkSmartPointer<vtkDoubleArray>::New();
99 array->SetNumberOfComponents(1);
100 array->SetName(name.c_str());
101
102 double value[1];
103 for (unsigned char i : *data) {
104 value[0] = i;
105 array->InsertNextTuple(value);
106 }
107
108 d_grid_p->GetPointData()->AddArray(array);
109}
110
111void rw::writer::VtkWriter::appendPointData(const std::string &name,
112 const std::vector<size_t> *data) {
113
114 auto array = vtkSmartPointer<vtkDoubleArray>::New();
115 array->SetNumberOfComponents(1);
116 array->SetName(name.c_str());
117
118 double value[1];
119 for (unsigned long i : *data) {
120 value[0] = i;
121 array->InsertNextTuple(value);
122 }
123
124 d_grid_p->GetPointData()->AddArray(array);
125}
126
127void rw::writer::VtkWriter::appendPointData(const std::string &name,
128 const std::vector<int> *data) {
129
130 auto array = vtkSmartPointer<vtkDoubleArray>::New();
131 array->SetNumberOfComponents(1);
132 array->SetName(name.c_str());
133
134 double value[1];
135 for (int i : *data) {
136 value[0] = i;
137 array->InsertNextTuple(value);
138 }
139
140 d_grid_p->GetPointData()->AddArray(array);
141}
142
143void rw::writer::VtkWriter::appendPointData(const std::string &name,
144 const std::vector<float> *data) {
145
146 auto array = vtkSmartPointer<vtkDoubleArray>::New();
147 array->SetNumberOfComponents(1);
148 array->SetName(name.c_str());
149
150 double value[1];
151 for (float i : *data) {
152 value[0] = i;
153 array->InsertNextTuple(value);
154 }
155
156 d_grid_p->GetPointData()->AddArray(array);
157}
158
159void rw::writer::VtkWriter::appendPointData(const std::string &name,
160 const std::vector<double> *data) {
161
162 auto array = vtkSmartPointer<vtkDoubleArray>::New();
163 array->SetNumberOfComponents(1);
164 array->SetName(name.c_str());
165
166 double value[1];
167 for (double i : *data) {
168 value[0] = i;
169 array->InsertNextTuple(value);
170 }
171
172 d_grid_p->GetPointData()->AddArray(array);
173}
174
176 const std::string &name, const std::vector<util::Point> *data) {
177
178 auto array = vtkSmartPointer<vtkDoubleArray>::New();
179 array->SetNumberOfComponents(3);
180 array->SetName(name.c_str());
181
182 array->SetComponentName(0, "x");
183 array->SetComponentName(1, "y");
184 array->SetComponentName(2, "z");
185
186 double value[3];
187 for (const auto &i : *data) {
188 value[0] = i.d_x;
189 value[1] = i.d_y;
190 value[2] = i.d_z;
191 array->InsertNextTuple(value);
192 }
193
194 d_grid_p->GetPointData()->AddArray(array);
195}
196
198 const std::string &name, const std::vector<util::SymMatrix3> *data) {
199
200 auto array = vtkSmartPointer<vtkDoubleArray>::New();
201 array->SetNumberOfComponents(6);
202 array->SetName(name.c_str());
203
204 array->SetComponentName(0, "xx");
205 array->SetComponentName(1, "yy");
206 array->SetComponentName(2, "zz");
207 array->SetComponentName(3, "yz");
208 array->SetComponentName(4, "xz");
209 array->SetComponentName(5, "xy");
210
211 double value[6];
212 for (const auto &i : *data) {
213 value[0] = i(0,0);
214 value[1] = i(1,1);
215 value[2] = i(2,2);
216 value[3] = i(1,2);
217 value[4] = i(0,2);
218 value[5] = i(0,1);
219 array->InsertNextTuple(value);
220 }
221
222 d_grid_p->GetPointData()->AddArray(array);
223}
224
225void rw::writer::VtkWriter::appendCellData(const std::string &name,
226 const std::vector<float> *data) {
227 auto array = vtkSmartPointer<vtkDoubleArray>::New();
228 array->SetNumberOfComponents(1);
229 array->SetName(name.c_str());
230
231 double value[1];
232 for (float i : *data) {
233 value[0] = i;
234 array->InsertNextTuple(value);
235 }
236
237 d_grid_p->GetCellData()->AddArray(array);
238}
239
241 const std::string &name, const std::vector<util::SymMatrix3> *data) {
242
243 auto array = vtkSmartPointer < vtkDoubleArray > ::New();
244 array->SetNumberOfComponents(6);
245 array->SetName(name.c_str());
246
247 array->SetComponentName(0, "xx");
248 array->SetComponentName(1, "yy");
249 array->SetComponentName(2, "zz");
250 array->SetComponentName(3, "yz");
251 array->SetComponentName(4, "xz");
252 array->SetComponentName(5, "xy");
253
254 double value[6];
255 for (const auto &i : *data) {
256 value[0] = i(0,0);
257 value[1] = i(1,1);
258 value[2] = i(2,2);
259 value[3] = i(1,2);
260 value[4] = i(0,2);
261 value[5] = i(0,1);
262 array->InsertNextTuple(value);
263 }
264
265 d_grid_p->GetCellData()->AddArray(array);
266}
267
268void rw::writer::VtkWriter::addTimeStep(const double &timestep) {
269
270 auto t = vtkDoubleArray::New();
271 t->SetName("TIME");
272 t->SetNumberOfTuples(1);
273 t->SetTuple1(0, timestep);
274 d_grid_p->GetFieldData()->AddArray(t);
275}
276
278 d_writer_p->SetInputData(d_grid_p);
279 d_writer_p->SetDataModeToAppended();
280 d_writer_p->EncodeAppendedDataOn();
281 if (d_compressType == "zlib")
282 d_writer_p->SetCompressorTypeToZLib();
283 else
284 d_writer_p->SetCompressor(0);
285 d_writer_p->Write();
286}
287
288void rw::writer::VtkWriter::appendFieldData(const std::string &name,
289 const double &data) {
290
291 auto t = vtkDoubleArray::New();
292 t->SetName(name.c_str());
293 t->SetNumberOfTuples(1);
294 t->SetTuple1(0, data);
295 d_grid_p->GetFieldData()->AddArray(t);
296}
297
298void rw::writer::VtkWriter::appendFieldData(const std::string &name,
299 const float &data) {
300
301 auto t = vtkDoubleArray::New();
302 t->SetName(name.c_str());
303 t->SetNumberOfTuples(1);
304 t->SetTuple1(0, data);
305 d_grid_p->GetFieldData()->AddArray(t);
306}
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 appendCellData(const std::string &name, const std::vector< float > *data)
Writes the float data associated to cells to the file.
vtkSmartPointer< vtkXMLUnstructuredGridWriter > d_writer_p
XML unstructured grid writer.
Definition vtkWriter.h:188
void addTimeStep(const double &timestep)
Writes the time step to the file.
void appendNodes(const std::vector< util::Point > *nodes, const std::vector< util::Point > *u=nullptr)
Writes the nodes to the file.
Definition vtkWriter.cpp:32
VtkWriter(const std::string &filename, const std::string &compress_type="")
Constructor.
Definition vtkWriter.cpp:22
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 vtkWriter.cpp:49
void appendPointData(const std::string &name, const std::vector< uint8_t > *data)
Writes the scalar point data to the file.
Definition vtkWriter.cpp:95
static int vtk_map_element_to_num_nodes[16]
Map from element type to number of nodes (for vtk)
A structure to represent 3d vectors.
Definition point.h:30
double d_y
the y coordinate
Definition point.h:36
double d_z
the z coordinate
Definition point.h:39
double d_x
the x coordinate
Definition point.h:33