13#include <vtkCellArray.h>
14#include <vtkCellData.h>
15#include <vtkDoubleArray.h>
17#include <vtkIntArray.h>
18#include <vtkPointData.h>
20#include <vtkUnsignedIntArray.h>
23 const std::string &compress_type)
24 : d_compressType(compress_type) {
26 std::string f = filename +
".vtu";
28 d_writer_p = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
29 d_writer_p->SetFileName(
const_cast<char *
>(f.c_str()));
33 const std::vector<util::Point> *u) {
35 auto points = vtkSmartPointer<vtkPoints>::New();
37 for (
size_t i = 0; i < nodes->size(); i++) {
42 points->InsertNextPoint(p.
d_x, p.
d_y, p.
d_z);
45 d_grid_p = vtkSmartPointer<vtkUnstructuredGrid>::New();
46 d_grid_p->SetPoints(points);
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) {
64 this->appendNodes(nodes, u);
68 size_t num_elems = en_con->size() / num_vertex;
74 auto cells = vtkSmartPointer<vtkCellArray>::New();
75 cells->Allocate(num_vertex, num_elems);
78 int cell_types[num_elems];
80 vtkIdType ids[num_vertex];
81 for (
size_t i = 0; i < num_elems; i++) {
84 for (
size_t k = 0; k < num_vertex; k++)
85 ids[k] = (*en_con)[num_vertex*i + k];
87 cells->InsertNextCell(num_vertex, ids);
88 cell_types[i] = element_type;
92 d_grid_p->SetCells(cell_types, cells);
96 const std::vector<uint8_t> *data) {
98 auto array = vtkSmartPointer<vtkDoubleArray>::New();
99 array->SetNumberOfComponents(1);
100 array->SetName(name.c_str());
103 for (
unsigned char i : *data) {
105 array->InsertNextTuple(value);
108 d_grid_p->GetPointData()->AddArray(array);
112 const std::vector<size_t> *data) {
114 auto array = vtkSmartPointer<vtkDoubleArray>::New();
115 array->SetNumberOfComponents(1);
116 array->SetName(name.c_str());
119 for (
unsigned long i : *data) {
121 array->InsertNextTuple(value);
124 d_grid_p->GetPointData()->AddArray(array);
128 const std::vector<int> *data) {
130 auto array = vtkSmartPointer<vtkDoubleArray>::New();
131 array->SetNumberOfComponents(1);
132 array->SetName(name.c_str());
135 for (
int i : *data) {
137 array->InsertNextTuple(value);
140 d_grid_p->GetPointData()->AddArray(array);
144 const std::vector<float> *data) {
146 auto array = vtkSmartPointer<vtkDoubleArray>::New();
147 array->SetNumberOfComponents(1);
148 array->SetName(name.c_str());
151 for (
float i : *data) {
153 array->InsertNextTuple(value);
156 d_grid_p->GetPointData()->AddArray(array);
160 const std::vector<double> *data) {
162 auto array = vtkSmartPointer<vtkDoubleArray>::New();
163 array->SetNumberOfComponents(1);
164 array->SetName(name.c_str());
167 for (
double i : *data) {
169 array->InsertNextTuple(value);
172 d_grid_p->GetPointData()->AddArray(array);
176 const std::string &name,
const std::vector<util::Point> *data) {
178 auto array = vtkSmartPointer<vtkDoubleArray>::New();
179 array->SetNumberOfComponents(3);
180 array->SetName(name.c_str());
182 array->SetComponentName(0,
"x");
183 array->SetComponentName(1,
"y");
184 array->SetComponentName(2,
"z");
187 for (
const auto &i : *data) {
191 array->InsertNextTuple(value);
194 d_grid_p->GetPointData()->AddArray(array);
198 const std::string &name,
const std::vector<util::SymMatrix3> *data) {
200 auto array = vtkSmartPointer<vtkDoubleArray>::New();
201 array->SetNumberOfComponents(6);
202 array->SetName(name.c_str());
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");
212 for (
const auto &i : *data) {
219 array->InsertNextTuple(value);
222 d_grid_p->GetPointData()->AddArray(array);
226 const std::vector<float> *data) {
227 auto array = vtkSmartPointer<vtkDoubleArray>::New();
228 array->SetNumberOfComponents(1);
229 array->SetName(name.c_str());
232 for (
float i : *data) {
234 array->InsertNextTuple(value);
237 d_grid_p->GetCellData()->AddArray(array);
241 const std::string &name,
const std::vector<util::SymMatrix3> *data) {
243 auto array = vtkSmartPointer < vtkDoubleArray > ::New();
244 array->SetNumberOfComponents(6);
245 array->SetName(name.c_str());
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");
255 for (
const auto &i : *data) {
262 array->InsertNextTuple(value);
265 d_grid_p->GetCellData()->AddArray(array);
270 auto t = vtkDoubleArray::New();
272 t->SetNumberOfTuples(1);
273 t->SetTuple1(0, timestep);
274 d_grid_p->GetFieldData()->AddArray(t);
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();
284 d_writer_p->SetCompressor(0);
289 const double &data) {
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);
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);
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.
void addTimeStep(const double ×tep)
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.
VtkWriter(const std::string &filename, const std::string &compress_type="")
Constructor.
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.
void appendPointData(const std::string &name, const std::vector< uint8_t > *data)
Writes the scalar point data to the file.
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.
double d_y
the y coordinate
double d_z
the z coordinate
double d_x
the x coordinate