10#ifndef INCG3e48c7b32a034cb3be3dbca884ff4f6c
11#define INCG3e48c7b32a034cb3be3dbca884ff4f6c
13#include <lf/mesh/mesh.h>
14#include <lf/mesh/utils/utils.h>
16#include <boost/variant/variant.hpp>
71 std::vector<std::vector<size_type>>
cells;
99 :
name(std::move(data_name)),
108 std::vector<Eigen::Matrix<T, 3, 1>>
data;
113 std::vector<Eigen::Matrix<T, 3, 1>>
data)
295 VtkWriter(std::shared_ptr<const mesh::Mesh> mesh, std::string filename,
296 dim_t codim = 0,
unsigned char order = 1);
325 unsigned char undefined_value = 0);
341 char undefined_value = 0);
357 unsigned int undefined_value = 0);
373 int undefined_value = 0);
389 float undefined_value = 0.F);
405 double undefined_value = 0.);
421 const Eigen::Vector2d& undefined_value = {0, 0});
437 const Eigen::Vector2f& undefined_value = {0, 0});
452 const mesh::utils::MeshDataSet<Eigen::Vector3d>& mds,
453 const Eigen::Vector3d& undefined_value = {0, 0, 0});
468 const mesh::utils::MeshDataSet<Eigen::Vector3f>& mds,
469 const Eigen::Vector3f& undefined_value = {0, 0, 0});
486 const std::string& name,
487 const mesh::utils::MeshDataSet<Eigen::VectorXd>& mds,
488 const Eigen::VectorXd& undefined_value = Eigen::Vector3d(0, 0, 0));
505 const std::string& name,
506 const mesh::utils::MeshDataSet<Eigen::VectorXf>& mds,
507 const Eigen::VectorXf& undefined_value = Eigen::Vector3f(0, 0, 0));
541 template <mesh::utils::MeshFunction MESH_FUNCTION>
543 const MESH_FUNCTION& mesh_function);
556 const mesh::utils::MeshDataSet<unsigned char>& mds,
557 unsigned char undefined_value = 0);
570 const mesh::utils::MeshDataSet<char>& mds,
571 char undefined_value = 0);
584 const mesh::utils::MeshDataSet<unsigned int>& mds,
585 unsigned int undefined_value = 0);
598 const mesh::utils::MeshDataSet<int>& mds,
599 int undefined_value = 0);
612 const mesh::utils::MeshDataSet<float>& mds,
613 float undefined_value = 0);
626 const mesh::utils::MeshDataSet<double>& mds,
627 double undefined_value = 0);
639 const mesh::utils::MeshDataSet<Eigen::Vector2d>& mds,
640 const Eigen::Vector2d& undefined_value = {0, 0});
652 const mesh::utils::MeshDataSet<Eigen::Vector2f>& mds,
653 const Eigen::Vector2f& undefined_value = {0, 0});
665 const mesh::utils::MeshDataSet<Eigen::Vector3d>& mds,
666 const Eigen::Vector3d& undefined_value = {0, 0, 0});
678 const mesh::utils::MeshDataSet<Eigen::Vector3f>& mds,
679 const Eigen::Vector3f& undefined_value = {0, 0, 0});
693 const std::string& name,
694 const mesh::utils::MeshDataSet<Eigen::VectorXd>& mds,
695 const Eigen::VectorXd& undefined_value = Eigen::Vector3d(0, 0, 0));
709 const std::string& name,
710 const mesh::utils::MeshDataSet<Eigen::VectorXf>& mds,
711 const Eigen::VectorXf& undefined_value = Eigen::Vector3f(0, 0, 0));
741 template <mesh::utils::MeshFunction MESH_FUNCTION>
743 const MESH_FUNCTION& mesh_function);
767 void WriteGlobalData(
const std::string& name, std::vector<float> data);
779 void WriteGlobalData(
const std::string& name, std::vector<double> data);
787 std::shared_ptr<const mesh::Mesh>
mesh_;
806 template <
int ROWS,
class T>
808 const std::string& name,
810 const Eigen::Matrix<T, ROWS, 1>& undefined_value);
817 template <
int ROWS,
class T>
819 const std::string& name,
821 const Eigen::Matrix<T, ROWS, 1>& undefined_value);
824 void WriteFieldData(
const std::string& name, std::vector<T> data);
826 template <
class DATA>
829 template <
int ROWS,
class T>
831 const Eigen::Matrix<T, ROWS, 1>& in);
834template <mesh::utils::MeshFunction MESH_FUNCTION>
837 const MESH_FUNCTION& mesh_function) {
842 auto dim_mesh =
mesh_->DimMesh();
844 const Eigen::Matrix<double, 0, 1> origin;
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>) {
856 for (
const auto* n :
mesh_->Entities(dim_mesh)) {
857 data.data[
mesh_->Index(*n)] = mesh_function(*n, origin)[0];
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();
867 auto values = mesh_function(*e,
aux_nodes_[ref_el.Id()]);
869 for (
typename std::vector<T>::size_type i = 0; i < values.size(); ++i) {
870 data.data[offset + i] = values[i];
876 static_assert(T::ColsAtCompileTime == 1,
877 "The MeshFunction must return row-vectors");
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;
885 std::is_same_v<double, Scalar> || std::is_same_v<float, Scalar>,
886 "The RowVectors returned by the MeshFunction must be either double "
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]);
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();
905 auto values = mesh_function(*e,
aux_nodes_[ref_el.Id()]);
907 for (
unsigned long i = 0; i < values.size(); ++i) {
908 PadWithZeros<T::RowsAtCompileTime, Scalar>(data.data[offset + i],
916 "MeshFunction values must be one of: unsigned char, char, "
917 "unsigned, int, float, double, Eigen::Vector<double, ...> "
918 "or Eigen::Vector<float, ...>");
922template <mesh::utils::MeshFunction MESH_FUNCTION>
924 const MESH_FUNCTION& mesh_function) {
929 std::vector<Eigen::VectorXd> barycenters(5);
939 static_assert(T::ColsAtCompileTime == 1,
940 "The MeshFunction must return row-vectors");
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;
948 std::is_same_v<double, Scalar> || std::is_same_v<float, Scalar>,
949 "The RowVectors returned by the MeshFunction must be either double "
954 constexpr int rows = T::RowsAtCompileTime;
956 this->PadWithZeros<rows, Scalar>(
957 data.data[
mesh_->Index(*p)],
958 mesh_function(*p, barycenters[p->RefEl().Id()])[0]);
966 data.data[
mesh_->Index(*e)] =
967 mesh_function(*e, barycenters[e->RefEl().Id()])[0];
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) ==
981 "There is already another Point/Cell Attribute Set with the "
985 if (name.find(
' ') != std::string::npos) {
987 "The name of the attribute set cannot contain spaces!");
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) {
995 out.template block<2, 1>(0, 0) = in;
997 }
else if constexpr (ROWS == 3) {
999 }
else if constexpr (ROWS == Eigen::Dynamic) {
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) {
1004 out.topRows(2) = in;
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.
constexpr unsigned int Id() const
Return a unique id for this reference element.
const Eigen::MatrixXd & NodeCoords() const
Get the coordinates of the nodes of this reference element.
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
static constexpr RefEl kTria()
Returns the reference triangle.
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
FieldDataArray(std::string name, std::vector< T > data)
Represents one set of attribute data (can be attached to points or cells)
ScalarData(std::string data_name, std::vector< T > data, std::string lookup_table="default")
std::vector< Eigen::Vector3f > points
std::vector< std::vector< size_type > > cells
std::vector< CellType > cell_types
VectorData(std::string name, std::vector< Eigen::Matrix< T, 3, 1 > > data)
std::vector< Eigen::Matrix< T, 3, 1 > > data
Representation of a VTK file (only relevant) features (advanced usage)
std::vector< boost::variant< FieldDataArray< int >, FieldDataArray< float >, FieldDataArray< double > > > FieldData
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
@ VTK_QUADRATIC_HEXAHEDRON
@ VTK_LAGRANGE_HEXAHEDRON
@ VTK_LAGRANGE_TETRAHEDRON
@ VTK_LAGRANGE_QUADRILATERAL
std::string header
The Vtk Header, at most 256 characters, no new lines characters!
UnstructuredGrid unstructured_grid
Describes the nodes + cells.
Format format
The format of the file.
Write a mesh along with mesh data into a vtk file.
void setBinary(bool binary)
Determines whether the Vtk file is written in binary or ASCII mode (default).
VtkWriter(VtkWriter &&)=delete
void WriteScalarPointData(const std::string &name, const mesh::utils::MeshDataSet< T > &mds, T undefined_value)
VtkWriter(const VtkWriter &)=delete
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.
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_
static void PadWithZeros(Eigen::Matrix< T, 3, 1 > &out, const Eigen::Matrix< T, ROWS, 1 > &in)
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)
void CheckAttributeSetName(const DATA &data, const std::string &name)
std::shared_ptr< const mesh::Mesh > mesh_
std::array< Eigen::MatrixXd, 5 > aux_nodes_
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<...>&
Check if a given type T is an Eigen::Matrix.
unsigned int dim_t
type for dimensions and co-dimensions and numbers derived from them
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.