From 406fedff8b792b226ce37f69942f89746dc2a2cc Mon Sep 17 00:00:00 2001 From: =?UTF-8?q?St=C3=A9phane=20Del=20Pino?= <stephane.delpino44@gmail.com> Date: Thu, 11 Feb 2021 16:11:26 +0100 Subject: [PATCH] Add handling of B, N, R^dxd data types in output --- src/output/OutputNamedItemValueSet.hpp | 14 ++++++- src/output/VTKWriter.cpp | 58 ++++++++++++++++++++++++++ src/output/VTKWriter.hpp | 21 ++++++++++ 3 files changed, 91 insertions(+), 2 deletions(-) diff --git a/src/output/OutputNamedItemValueSet.hpp b/src/output/OutputNamedItemValueSet.hpp index c474419ee..c7464135e 100644 --- a/src/output/OutputNamedItemValueSet.hpp +++ b/src/output/OutputNamedItemValueSet.hpp @@ -54,19 +54,29 @@ class NamedItemValue class OutputNamedItemValueSet { public: - using ItemValueVariant = std::variant<NodeValue<const int>, + using ItemValueVariant = std::variant<NodeValue<const bool>, + NodeValue<const int>, NodeValue<const long int>, + NodeValue<const unsigned long int>, NodeValue<const double>, NodeValue<const TinyVector<1, double>>, NodeValue<const TinyVector<2, double>>, NodeValue<const TinyVector<3, double>>, + NodeValue<const TinyMatrix<1, double>>, + NodeValue<const TinyMatrix<2, double>>, + NodeValue<const TinyMatrix<3, double>>, + CellValue<const bool>, CellValue<const int>, CellValue<const long int>, + CellValue<const unsigned long int>, CellValue<const double>, CellValue<const TinyVector<1, double>>, CellValue<const TinyVector<2, double>>, - CellValue<const TinyVector<3, double>>>; + CellValue<const TinyVector<3, double>>, + CellValue<const TinyMatrix<1, double>>, + CellValue<const TinyMatrix<2, double>>, + CellValue<const TinyMatrix<3, double>>>; private: std::map<std::string, ItemValueVariant> m_name_itemvariant_map; diff --git a/src/output/VTKWriter.cpp b/src/output/VTKWriter.cpp index 11a5d4bb5..103807dbc 100644 --- a/src/output/VTKWriter.cpp +++ b/src/output/VTKWriter.cpp @@ -47,6 +47,16 @@ VTKWriter::_write_node_pvtu(std::ofstream& os, << "\"/>\n"; } +template <size_t N, typename DataType> +void +VTKWriter::_write_node_pvtu(std::ofstream& os, + const std::string& name, + const NodeValue<const TinyMatrix<N, DataType>>&) const +{ + os << "<PDataArray type=\"" << VTKType<DataType>::name << "\" Name=\"" << name << "\" NumberOfComponents=\"" << N * N + << "\"/>\n"; +} + template <typename DataType> void VTKWriter::_write_node_pvtu(std::ofstream&, const std::string&, const CellValue<const DataType>&) const @@ -69,6 +79,16 @@ VTKWriter::_write_cell_pvtu(std::ofstream& os, << "\"/>\n"; } +template <size_t N, typename DataType> +void +VTKWriter::_write_cell_pvtu(std::ofstream& os, + const std::string& name, + const CellValue<const TinyMatrix<N, DataType>>&) const +{ + os << "<PDataArray type=\"" << VTKType<DataType>::name << "\" Name=\"" << name << "\" NumberOfComponents=\"" << N * N + << "\"/>\n"; +} + template <typename DataType> void VTKWriter::_write_cell_pvtu(std::ofstream&, const std::string&, const NodeValue<const DataType>&) const @@ -152,6 +172,25 @@ VTKWriter::_write_node_value(std::ofstream& os, os << "\n</DataArray>\n"; } +template <size_t N, typename DataType> +void +VTKWriter::_write_node_value(std::ofstream& os, + const std::string& name, + const NodeValue<const TinyMatrix<N, DataType>>& item_value) const +{ + os << "<DataArray type=\"" << VTKType<DataType>::name << "\" Name=\"" << name << "\" NumberOfComponents=\"" << N * N + << "\">\n"; + for (NodeId i = 0; i < item_value.size(); ++i) { + for (size_t j = 0; j < N; ++j) { + for (size_t k = 0; k < N; ++k) { + // The following '+' enforces integer output for char types + os << +item_value[i](j, k) << ' '; + } + } + } + os << "\n</DataArray>\n"; +} + template <typename DataType> void VTKWriter::_write_node_value(std::ofstream&, const std::string&, const CellValue<const DataType>&) const @@ -188,6 +227,25 @@ VTKWriter::_write_cell_value(std::ofstream& os, os << "\n</DataArray>\n"; } +template <size_t N, typename DataType> +void +VTKWriter::_write_cell_value(std::ofstream& os, + const std::string& name, + const CellValue<const TinyMatrix<N, DataType>>& item_value) const +{ + os << "<DataArray type=\"" << VTKType<DataType>::name << "\" Name=\"" << name << "\" NumberOfComponents=\"" << N * N + << "\">\n"; + for (CellId i = 0; i < item_value.size(); ++i) { + for (size_t j = 0; j < N; ++j) { + for (size_t k = 0; k < N; ++k) { + // The following '+' enforces integer output for char types + os << +item_value[i](j, k) << ' '; + } + } + } + os << "\n</DataArray>\n"; +} + template <typename DataType> void VTKWriter::_write_cell_value(std::ofstream&, const std::string&, const NodeValue<const DataType>&) const diff --git a/src/output/VTKWriter.hpp b/src/output/VTKWriter.hpp index 7661d0a27..fadec4f9d 100644 --- a/src/output/VTKWriter.hpp +++ b/src/output/VTKWriter.hpp @@ -3,6 +3,7 @@ #include <output/IWriter.hpp> +#include <algebra/TinyMatrix.hpp> #include <algebra/TinyVector.hpp> #include <output/OutputNamedItemValueSet.hpp> @@ -30,6 +31,11 @@ class VTKWriter : public IWriter const std::string& name, const NodeValue<const TinyVector<N, DataType>>&) const; + template <size_t N, typename DataType> + void _write_node_pvtu(std::ofstream& os, + const std::string& name, + const NodeValue<const TinyMatrix<N, DataType>>&) const; + template <typename DataType> void _write_node_pvtu(std::ofstream&, const std::string&, const CellValue<const DataType>&) const; @@ -41,6 +47,11 @@ class VTKWriter : public IWriter const std::string& name, const CellValue<const TinyVector<N, DataType>>&) const; + template <size_t N, typename DataType> + void _write_cell_pvtu(std::ofstream& os, + const std::string& name, + const CellValue<const TinyMatrix<N, DataType>>&) const; + template <typename DataType> void _write_cell_pvtu(std::ofstream&, const std::string&, const NodeValue<const DataType>&) const; @@ -61,6 +72,11 @@ class VTKWriter : public IWriter const std::string& name, const NodeValue<const TinyVector<N, DataType>>& item_value) const; + template <size_t N, typename DataType> + void _write_node_value(std::ofstream& os, + const std::string& name, + const NodeValue<const TinyMatrix<N, DataType>>& item_value) const; + template <typename DataType> void _write_node_value(std::ofstream&, const std::string&, const CellValue<const DataType>&) const; @@ -72,6 +88,11 @@ class VTKWriter : public IWriter const std::string& name, const CellValue<const TinyVector<N, DataType>>& item_value) const; + template <size_t N, typename DataType> + void _write_cell_value(std::ofstream& os, + const std::string& name, + const CellValue<const TinyMatrix<N, DataType>>& item_value) const; + template <typename DataType> void _write_cell_value(std::ofstream&, const std::string&, const NodeValue<const DataType>&) const; -- GitLab