LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
vtk_writer.h
1
10#ifndef INCG3e48c7b32a034cb3be3dbca884ff4f6c
11#define INCG3e48c7b32a034cb3be3dbca884ff4f6c
12
13#include <lf/mesh/mesh.h>
14#include <lf/mesh/utils/utils.h>
15
16#include <boost/variant/variant.hpp>
17#include <string>
18#include <utility>
19#include <vector>
20
21namespace lf::io {
22
31class VtkFile {
32 public:
33 using size_type = unsigned int;
34
37 enum class Format { ASCII, BINARY };
38
67
69 public:
70 std::vector<Eigen::Vector3f> points;
71 std::vector<std::vector<size_type>> cells;
72 std::vector<CellType> cell_types;
73 };
74
75 template <class T>
77 public:
78 std::string name;
79 std::vector<T> data;
80
81 FieldDataArray() = default;
82
83 explicit FieldDataArray(std::string name, std::vector<T> data)
84 : name(std::move(name)), data(std::move(data)) {}
85 };
86
88 template <class T>
89 class ScalarData {
90 public:
91 std::string name;
92 std::vector<T> data;
93 std::string lookup_table = "default";
94
95 ScalarData() = default;
96
97 explicit ScalarData(std::string data_name, std::vector<T> data,
98 std::string lookup_table = "default")
99 : name(std::move(data_name)),
100 data(std::move(data)),
101 lookup_table(std::move(lookup_table)) {}
102 };
103
104 template <class T>
106 public:
107 std::string name;
108 std::vector<Eigen::Matrix<T, 3, 1>> data;
109
110 VectorData() = default;
111
112 explicit VectorData(std::string name,
113 std::vector<Eigen::Matrix<T, 3, 1>> data)
114 : name(std::move(name)), data(std::move(data)) {}
115 };
116
117 // WARNING: the type long is interepreted differently depending on the
118 // platorm. I suggest to not use it at all.
119 using Attributes = std::vector<boost::variant<
124
125 using FieldData =
126 std::vector<boost::variant<FieldDataArray<int>, FieldDataArray<float>,
128
129 // Actual members
131
133 std::string header;
136
139
141
142 // Data that is attached to points
144
145 // Data that is attached to the cells of the mesh
147};
148
149void WriteToFile(const VtkFile& vtk_file, const std::string& filename);
150
151// clang-format off
152
271// clang-format on
273 public:
275
276 VtkWriter(const VtkWriter&) = delete;
277 VtkWriter(VtkWriter&&) = delete;
278 VtkWriter& operator=(const VtkWriter&) = delete;
280
295 VtkWriter(std::shared_ptr<const mesh::Mesh> mesh, std::string filename,
296 dim_t codim = 0, unsigned char order = 1);
297
303 void setBinary(bool binary) {
304 if (binary) {
306 } else {
308 }
309 }
310
323 void WritePointData(const std::string& name,
325 unsigned char undefined_value = 0);
326
339 void WritePointData(const std::string& name,
341 char undefined_value = 0);
342
355 void WritePointData(const std::string& name,
357 unsigned int undefined_value = 0);
358
371 void WritePointData(const std::string& name,
373 int undefined_value = 0);
374
387 void WritePointData(const std::string& name,
389 float undefined_value = 0.F);
390
403 void WritePointData(const std::string& name,
405 double undefined_value = 0.);
406
419 void WritePointData(const std::string& name,
421 const Eigen::Vector2d& undefined_value = {0, 0});
422
435 void WritePointData(const std::string& name,
437 const Eigen::Vector2f& undefined_value = {0, 0});
438
451 void WritePointData(const std::string& name,
452 const mesh::utils::MeshDataSet<Eigen::Vector3d>& mds,
453 const Eigen::Vector3d& undefined_value = {0, 0, 0});
454
467 void WritePointData(const std::string& name,
468 const mesh::utils::MeshDataSet<Eigen::Vector3f>& mds,
469 const Eigen::Vector3f& undefined_value = {0, 0, 0});
470
485 void WritePointData(
486 const std::string& name,
487 const mesh::utils::MeshDataSet<Eigen::VectorXd>& mds,
488 const Eigen::VectorXd& undefined_value = Eigen::Vector3d(0, 0, 0));
489
504 void WritePointData(
505 const std::string& name,
506 const mesh::utils::MeshDataSet<Eigen::VectorXf>& mds,
507 const Eigen::VectorXf& undefined_value = Eigen::Vector3f(0, 0, 0));
508
541 template <mesh::utils::MeshFunction MESH_FUNCTION>
542 void WritePointData(const std::string& name,
543 const MESH_FUNCTION& mesh_function);
544
555 void WriteCellData(const std::string& name,
556 const mesh::utils::MeshDataSet<unsigned char>& mds,
557 unsigned char undefined_value = 0);
558
569 void WriteCellData(const std::string& name,
570 const mesh::utils::MeshDataSet<char>& mds,
571 char undefined_value = 0);
572
583 void WriteCellData(const std::string& name,
584 const mesh::utils::MeshDataSet<unsigned int>& mds,
585 unsigned int undefined_value = 0);
586
597 void WriteCellData(const std::string& name,
598 const mesh::utils::MeshDataSet<int>& mds,
599 int undefined_value = 0);
600
611 void WriteCellData(const std::string& name,
612 const mesh::utils::MeshDataSet<float>& mds,
613 float undefined_value = 0);
614
625 void WriteCellData(const std::string& name,
626 const mesh::utils::MeshDataSet<double>& mds,
627 double undefined_value = 0);
628
638 void WriteCellData(const std::string& name,
639 const mesh::utils::MeshDataSet<Eigen::Vector2d>& mds,
640 const Eigen::Vector2d& undefined_value = {0, 0});
641
651 void WriteCellData(const std::string& name,
652 const mesh::utils::MeshDataSet<Eigen::Vector2f>& mds,
653 const Eigen::Vector2f& undefined_value = {0, 0});
654
664 void WriteCellData(const std::string& name,
665 const mesh::utils::MeshDataSet<Eigen::Vector3d>& mds,
666 const Eigen::Vector3d& undefined_value = {0, 0, 0});
667
677 void WriteCellData(const std::string& name,
678 const mesh::utils::MeshDataSet<Eigen::Vector3f>& mds,
679 const Eigen::Vector3f& undefined_value = {0, 0, 0});
680
692 void WriteCellData(
693 const std::string& name,
694 const mesh::utils::MeshDataSet<Eigen::VectorXd>& mds,
695 const Eigen::VectorXd& undefined_value = Eigen::Vector3d(0, 0, 0));
696
708 void WriteCellData(
709 const std::string& name,
710 const mesh::utils::MeshDataSet<Eigen::VectorXf>& mds,
711 const Eigen::VectorXf& undefined_value = Eigen::Vector3f(0, 0, 0));
712
741 template <mesh::utils::MeshFunction MESH_FUNCTION>
742 void WriteCellData(const std::string& name,
743 const MESH_FUNCTION& mesh_function);
744
755 void WriteGlobalData(const std::string& name, std::vector<int> data);
756
767 void WriteGlobalData(const std::string& name, std::vector<float> data);
768
779 void WriteGlobalData(const std::string& name, std::vector<double> data);
780
785
786 private:
787 std::shared_ptr<const mesh::Mesh> mesh_;
789 std::string filename_;
791 unsigned char order_;
792
793 // entry [i] stores the reference coordinates for the auxilliary nodes for
794 // RefEl.Id() == i
795 std::array<Eigen::MatrixXd, 5> aux_nodes_;
796
797 // entry[cd] contains the index of the first auxiliary node for the codim=cd
798 // entity
799 std::array<mesh::utils::CodimMeshDataSet<unsigned int>, 2> aux_node_offset_;
800
801 template <class T>
802 void WriteScalarPointData(const std::string& name,
804 T undefined_value);
805
806 template <int ROWS, class T>
808 const std::string& name,
809 const mesh::utils::MeshDataSet<Eigen::Matrix<T, ROWS, 1>>& mds,
810 const Eigen::Matrix<T, ROWS, 1>& undefined_value);
811
812 template <class T>
813 void WriteScalarCellData(const std::string& name,
815 T undefined_value);
816
817 template <int ROWS, class T>
819 const std::string& name,
820 const mesh::utils::MeshDataSet<Eigen::Matrix<T, ROWS, 1>>& mds,
821 const Eigen::Matrix<T, ROWS, 1>& undefined_value);
822
823 template <class T>
824 void WriteFieldData(const std::string& name, std::vector<T> data);
825
826 template <class DATA>
827 void CheckAttributeSetName(const DATA& data, const std::string& name);
828
829 template <int ROWS, class T>
830 static void PadWithZeros(Eigen::Matrix<T, 3, 1>& out,
831 const Eigen::Matrix<T, ROWS, 1>& in);
832};
833
834template <mesh::utils::MeshFunction MESH_FUNCTION>
835// NOLINTNEXTLINE(readability-function-cognitive-complexity)
836void VtkWriter::WritePointData(const std::string& name,
837 const MESH_FUNCTION& mesh_function) {
838 // we have to evaluate the mesh function at all points
839 // and write this into the vtk file:
842 auto dim_mesh = mesh_->DimMesh();
843
844 const Eigen::Matrix<double, 0, 1> origin;
845
846 if constexpr (std::is_same_v<T, unsigned char> || std::is_same_v<T, char> ||
847 std::is_same_v<T, unsigned> || std::is_same_v<T, int> ||
848 std::is_same_v<T, float> ||
849 std::is_same_v<T, double>) { // NOLINT
850 // MeshFunction is scalar valued:
852 data.data.resize(vtk_file_.unstructured_grid.points.size());
853 data.name = name;
854
855 // evaluate at nodes:
856 for (const auto* n : mesh_->Entities(dim_mesh)) {
857 data.data[mesh_->Index(*n)] = mesh_function(*n, origin)[0];
858 }
859
860 for (int codim = static_cast<int>(dim_mesh - 1);
861 codim >= static_cast<char>(codim_); --codim) {
862 for (const auto* e : mesh_->Entities(codim)) {
863 auto ref_el = e->RefEl();
864 if (order_ == 1 || (order_ == 2 && ref_el == base::RefEl::kTria())) {
865 continue;
866 }
867 auto values = mesh_function(*e, aux_nodes_[ref_el.Id()]);
868 auto offset = aux_node_offset_[codim](*e);
869 for (typename std::vector<T>::size_type i = 0; i < values.size(); ++i) {
870 data.data[offset + i] = values[i];
871 }
872 }
873 }
874 vtk_file_.point_data.push_back(std::move(data));
875 } else if constexpr (base::EigenMatrix<T>) { // NOLINT
876 static_assert(T::ColsAtCompileTime == 1,
877 "The MeshFunction must return row-vectors");
878 static_assert(
879 T::RowsAtCompileTime == Eigen::Dynamic ||
880 (T::RowsAtCompileTime > 1 && T::RowsAtCompileTime < 4),
881 "The Row vectors returned by the MeshFunction must be either dynamic "
882 "or must have 2 or 3 rows (at compile time).");
883 using Scalar = typename T::Scalar;
884 static_assert(
885 std::is_same_v<double, Scalar> || std::is_same_v<float, Scalar>,
886 "The RowVectors returned by the MeshFunction must be either double "
887 "or float valued.");
889 data.data.resize(vtk_file_.unstructured_grid.points.size());
890 data.name = name;
891
892 // evaluate at nodes:
893 for (const auto* n : mesh_->Entities(dim_mesh)) {
894 PadWithZeros<T::RowsAtCompileTime, Scalar>(data.data[mesh_->Index(*n)],
895 mesh_function(*n, origin)[0]);
896 }
897
898 for (int codim = static_cast<int>(dim_mesh - 1);
899 codim >= static_cast<char>(codim_); --codim) {
900 for (const auto* e : mesh_->Entities(codim)) {
901 auto ref_el = e->RefEl();
902 if (order_ < 2 || (order_ < 3 && ref_el == base::RefEl::kTria())) {
903 continue;
904 }
905 auto values = mesh_function(*e, aux_nodes_[ref_el.Id()]);
906 auto offset = aux_node_offset_[codim](*e);
907 for (unsigned long i = 0; i < values.size(); ++i) {
908 PadWithZeros<T::RowsAtCompileTime, Scalar>(data.data[offset + i],
909 values[i]);
910 }
911 }
912 }
913 vtk_file_.point_data.push_back(std::move(data));
914 } else {
915 LF_VERIFY_MSG(false,
916 "MeshFunction values must be one of: unsigned char, char, "
917 "unsigned, int, float, double, Eigen::Vector<double, ...> "
918 "or Eigen::Vector<float, ...>");
919 }
920}
921
922template <mesh::utils::MeshFunction MESH_FUNCTION>
923void VtkWriter::WriteCellData(const std::string& name,
924 const MESH_FUNCTION& mesh_function) {
927
928 // maps from RefEl::Id() -> barycenter of the reference element
929 std::vector<Eigen::VectorXd> barycenters(5);
930 barycenters[base::RefEl::kPoint().Id()] = Eigen::Matrix<double, 0, 1>();
931 barycenters[base::RefEl::kSegment().Id()] =
932 base::RefEl::kSegment().NodeCoords().rowwise().sum() / 2.;
933 barycenters[base::RefEl::kTria().Id()] =
934 base::RefEl::kTria().NodeCoords().rowwise().sum() / 3.;
935 barycenters[base::RefEl::kQuad().Id()] =
936 base::RefEl::kTria().NodeCoords().rowwise().sum() / 4.;
937
938 if constexpr (base::EigenArray<T> || base::EigenMatrix<T>) {
939 static_assert(T::ColsAtCompileTime == 1,
940 "The MeshFunction must return row-vectors");
941 static_assert(
942 T::RowsAtCompileTime == Eigen::Dynamic ||
943 (T::RowsAtCompileTime > 1 && T::RowsAtCompileTime < 4),
944 "The Row vectors returned by the MeshFunction must be either dynamic "
945 "or must have 2 or 3 rows (at compile time).");
946 using Scalar = typename T::Scalar;
947 static_assert(
948 std::is_same_v<double, Scalar> || std::is_same_v<float, Scalar>,
949 "The RowVectors returned by the MeshFunction must be either double "
950 "or float valued.");
952 data.data.resize(mesh_->NumEntities(codim_));
953 data.name = name;
954 constexpr int rows = T::RowsAtCompileTime;
955 for (const auto* p : mesh_->Entities(codim_)) {
956 this->PadWithZeros<rows, Scalar>(
957 data.data[mesh_->Index(*p)],
958 mesh_function(*p, barycenters[p->RefEl().Id()])[0]);
959 }
960 vtk_file_.cell_data.push_back(std::move(data));
961 } else {
963 data.data.resize(mesh_->NumEntities(codim_));
964 data.name = name;
965 for (const auto* e : mesh_->Entities(codim_)) {
966 data.data[mesh_->Index(*e)] =
967 mesh_function(*e, barycenters[e->RefEl().Id()])[0];
968 }
969 vtk_file_.cell_data.push_back(std::move(data));
970 }
971}
972
973template <class DATA>
975 const std::string& name) {
976 if (std::find_if(data.begin(), data.end(), [&](auto& d) {
977 return boost::apply_visitor([&](auto&& d2) { return d2.name; }, d) ==
978 name;
979 }) != data.end()) {
980 throw base::LfException(
981 "There is already another Point/Cell Attribute Set with the "
982 "name " +
983 name);
984 }
985 if (name.find(' ') != std::string::npos) {
986 throw base::LfException(
987 "The name of the attribute set cannot contain spaces!");
988 }
989}
990
991template <int ROWS, class T>
992void VtkWriter::PadWithZeros(Eigen::Matrix<T, 3, 1>& out,
993 const Eigen::Matrix<T, ROWS, 1>& in) {
994 if constexpr (ROWS == 2) { // NOLINT
995 out.template block<2, 1>(0, 0) = in;
996 out(2) = T(0); // NOLINT
997 } else if constexpr (ROWS == 3) { // NOLINT
998 out = in;
999 } else if constexpr (ROWS == Eigen::Dynamic) { // NOLINT
1000 LF_ASSERT_MSG(in.rows() == 2 || in.rows() == 3,
1001 "VtkWriter can only write out 2d or 3d vectors.");
1002 if (in.rows() == 2) {
1003 out(2) = T{0};
1004 out.topRows(2) = in;
1005 } else {
1006 out = in;
1007 }
1008 }
1009}
1010
1011} // namespace lf::io
1012
1013#endif // INCG3e48c7b32a034cb3be3dbca884ff4f6c
A simple, general purpose exception class that is thrown by LehrFEM++ if something is wrong.
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
Definition ref_el.h:153
constexpr unsigned int Id() const
Return a unique id for this reference element.
Definition ref_el.h:494
const Eigen::MatrixXd & NodeCoords() const
Get the coordinates of the nodes of this reference element.
Definition ref_el.h:241
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Definition ref_el.h:144
static constexpr RefEl kTria()
Returns the reference triangle.
Definition ref_el.h:161
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Definition ref_el.h:169
FieldDataArray(std::string name, std::vector< T > data)
Definition vtk_writer.h:83
Represents one set of attribute data (can be attached to points or cells)
Definition vtk_writer.h:89
std::vector< T > data
Definition vtk_writer.h:92
ScalarData(std::string data_name, std::vector< T > data, std::string lookup_table="default")
Definition vtk_writer.h:97
std::vector< Eigen::Vector3f > points
Definition vtk_writer.h:70
std::vector< std::vector< size_type > > cells
Definition vtk_writer.h:71
std::vector< CellType > cell_types
Definition vtk_writer.h:72
VectorData(std::string name, std::vector< Eigen::Matrix< T, 3, 1 > > data)
Definition vtk_writer.h:112
std::vector< Eigen::Matrix< T, 3, 1 > > data
Definition vtk_writer.h:108
Representation of a VTK file (only relevant) features (advanced usage)
Definition vtk_writer.h:31
FieldData field_data
Definition vtk_writer.h:140
std::vector< boost::variant< FieldDataArray< int >, FieldDataArray< float >, FieldDataArray< double > > > FieldData
Definition vtk_writer.h:125
std::vector< boost::variant< ScalarData< char >, ScalarData< unsigned char >, ScalarData< short >, ScalarData< unsigned short >, ScalarData< int >, ScalarData< unsigned int >, ScalarData< long >, ScalarData< unsigned long >, ScalarData< float >, ScalarData< double >, VectorData< float >, VectorData< double > > > Attributes
Definition vtk_writer.h:119
std::string header
The Vtk Header, at most 256 characters, no new lines characters!
Definition vtk_writer.h:133
UnstructuredGrid unstructured_grid
Describes the nodes + cells.
Definition vtk_writer.h:138
Format
Nested types.
Definition vtk_writer.h:37
Attributes cell_data
Definition vtk_writer.h:146
Format format
The format of the file.
Definition vtk_writer.h:135
unsigned int size_type
Definition vtk_writer.h:33
Attributes point_data
Definition vtk_writer.h:143
Write a mesh along with mesh data into a vtk file.
Definition vtk_writer.h:272
void setBinary(bool binary)
Determines whether the Vtk file is written in binary or ASCII mode (default).
Definition vtk_writer.h:303
VtkWriter(VtkWriter &&)=delete
void WriteScalarPointData(const std::string &name, const mesh::utils::MeshDataSet< T > &mds, T undefined_value)
VtkWriter(const VtkWriter &)=delete
base::dim_t dim_t
Definition vtk_writer.h:274
void WriteCellData(const std::string &name, const mesh::utils::MeshDataSet< unsigned char > &mds, unsigned char undefined_value=0)
Add a new unsigned char attribute dataset that attaches data to the cells of the mesh (i....
~VtkWriter()
Destructor, writes everything into the file and closes it.
Definition vtk_writer.h:784
void WriteScalarCellData(const std::string &name, const mesh::utils::MeshDataSet< T > &mds, T undefined_value)
VtkWriter & operator=(VtkWriter &&)=delete
std::array< mesh::utils::CodimMeshDataSet< unsigned int >, 2 > aux_node_offset_
Definition vtk_writer.h:799
static void PadWithZeros(Eigen::Matrix< T, 3, 1 > &out, const Eigen::Matrix< T, ROWS, 1 > &in)
Definition vtk_writer.h:992
void WriteFieldData(const std::string &name, std::vector< T > data)
VtkWriter & operator=(const VtkWriter &)=delete
void WriteGlobalData(const std::string &name, std::vector< int > data)
Write global data into the vtk file that is not related to the mesh at all.
void WriteVectorPointData(const std::string &name, const mesh::utils::MeshDataSet< Eigen::Matrix< T, ROWS, 1 > > &mds, const Eigen::Matrix< T, ROWS, 1 > &undefined_value)
void WriteVectorCellData(const std::string &name, const mesh::utils::MeshDataSet< Eigen::Matrix< T, ROWS, 1 > > &mds, const Eigen::Matrix< T, ROWS, 1 > &undefined_value)
unsigned char order_
Definition vtk_writer.h:791
void CheckAttributeSetName(const DATA &data, const std::string &name)
Definition vtk_writer.h:974
std::string filename_
Definition vtk_writer.h:789
std::shared_ptr< const mesh::Mesh > mesh_
Definition vtk_writer.h:787
std::array< Eigen::MatrixXd, 5 > aux_nodes_
Definition vtk_writer.h:795
void WritePointData(const std::string &name, const mesh::utils::MeshDataSet< unsigned char > &mds, unsigned char undefined_value=0)
Add a new unsigned char attribute dataset that attaches data to the points/nodes of the mesh.
Interface that specifies how data is stored with an entity.
Check if a given type T is a Eigen::Array<...>&
Definition eigen_tools.h:86
Check if a given type T is an Eigen::Matrix.
Definition eigen_tools.h:70
unsigned int dim_t
type for dimensions and co-dimensions and numbers derived from them
Definition types.h:32
Mesh input (from file) and output (in various formats) facilities.
void WriteToFile(const VtkFile &vtk_file, const std::string &filename)
internal::VectorElement_t< internal::MeshFunctionReturnType_t< R > > MeshFunctionReturnType
Determine the type of objects returned by a MeshFunction.