PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
vtkWriter.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_VTKWRITER_H
12#define RW_VTKWRITER_H
13
14#include <vtkSmartPointer.h>
15#include <vtkUnstructuredGrid.h>
16#include <vtkXMLUnstructuredGridWriter.h>
17
18#include "util/point.h" // definition of Point
19#include "util/matrix.h" // definition of SymMatrix3
20#include <cstdint> // uint8_t type
21#include <cstring> // string and size_t type
22
23namespace rw {
24
25namespace writer {
26
28class VtkWriter {
29
30public:
40 explicit VtkWriter(const std::string &filename, const std::string &compress_type = "");
41
52 void appendNodes(const std::vector<util::Point> *nodes,
53 const std::vector<util::Point> *u = nullptr);
54
63 void appendMesh(const std::vector<util::Point> *nodes,
64 const size_t &element_type,
65 const std::vector<size_t> *en_con,
66 const std::vector<util::Point> *u = nullptr);
67
80 void appendPointData(const std::string &name,
81 const std::vector<uint8_t> *data);
82
88 void appendPointData(const std::string &name,
89 const std::vector<size_t> *data);
90
96 void appendPointData(const std::string &name, const std::vector<int> *data);
97
103 void appendPointData(const std::string &name, const std::vector<float> *data);
104
110 void appendPointData(const std::string &name,
111 const std::vector<double> *data);
112
118 void appendPointData(const std::string &name,
119 const std::vector<util::Point> *data);
120
127 void appendPointData(const std::string &name,
128 const std::vector<util::SymMatrix3> *data);
129
142 void appendCellData(const std::string &name, const std::vector<float> *data);
143
149 void appendCellData(const std::string &name,
150 const std::vector<util::SymMatrix3> *data);
151
164 void appendFieldData(const std::string &name, const double &data);
165
171 void appendFieldData(const std::string &name, const float &data);
172
177 void addTimeStep(const double &timestep);
178
184 void close();
185
186private:
188 vtkSmartPointer<vtkXMLUnstructuredGridWriter> d_writer_p;
189
191 vtkSmartPointer<vtkUnstructuredGrid> d_grid_p;
192
194 std::string d_compressType;
195};
196
197} // namespace writer
198
199} // namespace rw
200
201#endif // RW_VTKWRITER_H
A vtk writer for simple point data and complex fem mesh data.
Definition vtkWriter.h:28
void appendFieldData(const std::string &name, const double &data)
Writes the scalar field data to the file.
std::string d_compressType
compression_type Specify the compressor (if any)
Definition vtkWriter.h:194
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
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
vtkSmartPointer< vtkUnstructuredGrid > d_grid_p
Unstructured grid.
Definition vtkWriter.h:191
void appendPointData(const std::string &name, const std::vector< uint8_t > *data)
Writes the scalar point data to the file.
Definition vtkWriter.cpp:95
Collection of methods and database related to reading and writing.