11#include <lf/base/base.h>
14#include <boost/fusion/include/adapt_struct.hpp>
15#include <boost/phoenix/phoenix.hpp>
16#include <boost/phoenix/scope/let.hpp>
17#include <boost/spirit/include/karma.hpp>
19#include <unsupported/Eigen/KroneckerProduct>
21#include "eigen_fusion_adapter.h"
23template <
class RESULT_TYPE,
class LAMBDA>
30 return lambda_(std::forward<T>(arg));
39template <
class RESULT_TYPE,
class LAMBDA>
46 (std::vector<Eigen::Vector3f>, points)
47 (std::vector<std::vector<unsigned int>>, cells)
48 (std::vector<lf::io::VtkFile::CellType>, cell_types)
58 (std::string, lookup_table)
92void cell_list_size(
unsigned int& result,
93 const std::vector<std::vector<unsigned int>>& cells) {
94 result = cells.size();
95 for (
const auto& v : cells) {
102BOOST_PHOENIX_ADAPT_FUNCTION(
void, CellListSize, cell_list_size, 2);
120namespace karma = boost::spirit::karma;
121template <
class Iterator,
bool BINARY>
122struct VtkGrammar : karma::grammar<Iterator, VtkFile()> {
129 VtkGrammar() : VtkGrammar::base_type(root) {
130 using boost::phoenix::at_c;
131 using boost::phoenix::size;
141 points_vector %= *(karma::big_bin_float << karma::big_bin_float
142 << karma::big_bin_float)
145 points_vector %= *(karma::float_ <<
' ' << karma::float_ <<
' '
146 << karma::float_ << karma::eol);
149 points = lit(
"POINTS ")
150 << karma::uint_[karma::_1 = boost::phoenix::size(karma::_val)]
151 << lit(
" float") << karma::eol
152 << points_vector[karma::_1 = karma::_val];
156 cells_vector %= karma::big_dword(boost::phoenix::size(karma::_val))
157 << *(karma::big_dword);
162 cells_vector %= lit(size(_val))
163 << lit(
' ') << (karma::uint_ %
" ") << karma::eol;
165 cells = lit(
"CELLS ")
166 << karma::uint_[karma::_1 = boost::phoenix::size(karma::_val)]
167 << lit(
' ') << karma::uint_[CellListSize(karma::_1, karma::_val)]
168 << karma::eol << (*cells_vector)[karma::_1 = karma::_val] << eol;
172 karma::big_dword[karma::_1 =
173 boost::phoenix::static_cast_<int>(karma::_val)];
175 cell_type = karma::int_[karma::_1 = boost::phoenix::static_cast_<int>(
179 cell_types = lit(
"CELL_TYPES ")
180 << karma::uint_[karma::_1 = boost::phoenix::size(karma::_val)]
181 << karma::eol << (*cell_type)[karma::_1 = karma::_val];
183 unstructured_grid %= lit(
"DATASET UNSTRUCTURED_GRID")
184 << karma::eol << points << cells << cell_types;
188 field_data_int %= karma::string << lit(
" 1 ") << lit(size(at_c<1>(_val)))
189 << lit(
" int") << eol << *karma::big_dword
191 field_data_float %= karma::string
192 << lit(
" 1 ") << lit(size(at_c<1>(_val)))
193 << lit(
" float") << eol << *(karma::big_bin_float)
196 field_data_double %= karma::string
197 << lit(
" 1 ") << lit(size(at_c<1>(_val)))
198 << lit(
" double") << eol << *karma::big_bin_double
201 field_data_int %= karma::string << lit(
" 1 ") << lit(size(at_c<1>(_val)))
202 << lit(
" int") << eol
203 << (karma::int_ % lit(
' ')) << eol;
204 field_data_float %= karma::string
205 << lit(
" 1 ") << lit(size(at_c<1>(_val)))
206 << lit(
" float") << eol << (karma::float_ % lit(
' '))
209 field_data_double %= karma::string << lit(
" 1 ")
210 << lit(size(at_c<1>(_val)))
211 << lit(
" double") << eol
212 << (karma::double_ % lit(
' ')) << eol;
214 field_data %= lit(
"FIELD FieldData ")
215 << lit(size(_val)) << eol
216 << *(field_data_int | field_data_float | field_data_double);
221 (lit(
"SCALARS ") << karma::string << lit(
" char 1") << eol
222 << lit(
"LOOKUP_TABLE ") << karma::string << eol
227 (lit(
"SCALARS ") << karma::string << lit(
" char 1") << eol
228 << lit(
"LOOKUP_TABLE ") << karma::string << eol
233 (lit(
"SCALARS ") << karma::string << lit(
" short 1") << eol
234 << lit(
"LOOKUP_TABLE ") << karma::string << eol
238 scalar_data_ushort %=
239 (lit(
"SCALARS ") << karma::string << lit(
" unsigned_short 1") << eol
240 << lit(
"LOOKUP_TABLE ") << karma::string << eol
245 (lit(
"SCALARS ") << karma::string << lit(
" int 1") << eol
246 << lit(
"LOOKUP_TABLE ") << karma::string << eol
247 << *karma::big_dword)
249 scalar_data_unsigned_int %=
250 (lit(
"SCALARS ") << karma::string << lit(
" unsigned_int 1") << eol
251 << lit(
"LOOKUP_TABLE ") << karma::string << eol
252 << *karma::big_dword)
256 if constexpr (
sizeof(long) == 8) {
259 << karma::string << lit(
" long 1") << eol << lit(
"LOOKUP_TABLE ")
260 << karma::string << eol << *karma::big_qword)
262 scalar_data_unsigned_long %=
263 (lit(
"SCALARS ") << karma::string << lit(
" unsigned_long 1") << eol
264 << lit(
"LOOKUP_TABLE ") << karma::string << eol
265 << *karma::big_qword)
270 << karma::string << lit(
" long 1") << eol << lit(
"LOOKUP_TABLE ")
271 << karma::string << eol << *karma::big_dword)
273 scalar_data_unsigned_long %=
274 (lit(
"SCALARS ") << karma::string << lit(
" unsigned_long 1") << eol
275 << lit(
"LOOKUP_TABLE ") << karma::string << eol
276 << *karma::big_dword)
281 (lit(
"SCALARS ") << karma::string << lit(
" float 1") << eol
282 << lit(
"LOOKUP_TABLE ") << karma::string << eol
283 << *karma::big_bin_float)
286 scalar_data_double %=
287 (lit(
"SCALARS ") << karma::string << lit(
" double 1") << eol
288 << lit(
"LOOKUP_TABLE ") << karma::string << eol
289 << *karma::big_bin_double)
292 vector_data_float %= lit(
"VECTORS ")
293 << karma::string << lit(
" float") << eol
294 << *(karma::big_bin_float << karma::big_bin_float
295 << karma::big_bin_float)
298 vector_data_double %= lit(
"VECTORS ")
299 << karma::string << lit(
" double") << eol
300 << *(karma::big_bin_double << karma::big_bin_double
301 << karma::big_bin_double)
305 (lit(
"SCALARS ") << karma::string << lit(
" char 1") << eol
306 << lit(
"LOOKUP_TABLE ") << karma::string << eol
307 << (*(karma::short_ % eol)))
311 (lit(
"SCALARS ") << karma::string << lit(
" char 1") << eol
312 << lit(
"LOOKUP_TABLE ") << karma::string << eol
313 << (*(karma::ushort_ % eol)))
317 (lit(
"SCALARS ") << karma::string << lit(
" short 1") << eol
318 << lit(
"LOOKUP_TABLE ") << karma::string << eol
319 << (*(karma::short_ % eol)))
322 scalar_data_ushort %=
323 (lit(
"SCALARS ") << karma::string << lit(
" unsigned_short 1") << eol
324 << lit(
"LOOKUP_TABLE ") << karma::string << eol
325 << (*(karma::ushort_ % eol)))
329 (lit(
"SCALARS ") << karma::string << lit(
" int 1") << eol
330 << lit(
"LOOKUP_TABLE ") << karma::string << eol
331 << (*(karma::int_ % eol)))
333 scalar_data_unsigned_int %=
334 (lit(
"SCALARS ") << karma::string << lit(
" unsigned_int 1") << eol
335 << lit(
"LOOKUP_TABLE ") << karma::string << eol
336 << (*(karma::uint_ % eol)))
340 (lit(
"SCALARS ") << karma::string << lit(
" long 1") << eol
341 << lit(
"LOOKUP_TABLE ") << karma::string << eol
342 << (*(karma::long_ % eol)))
345 scalar_data_unsigned_long %=
346 (lit(
"SCALARS ") << karma::string << lit(
" unsigned_long 1") << eol
347 << lit(
"LOOKUP_TABLE ") << karma::string << eol
348 << (*(karma::ulong_ % eol)))
352 (lit(
"SCALARS ") << karma::string << lit(
" float 1") << eol
353 << lit(
"LOOKUP_TABLE ") << karma::string << eol
354 << (*(karma::float_ % eol)))
357 scalar_data_double %=
358 (lit(
"SCALARS ") << karma::string << lit(
" double 1") << eol
359 << lit(
"LOOKUP_TABLE ") << karma::string << eol
360 << (*(karma::double_ % eol)))
363 vector_data_float %= lit(
"VECTORS ")
364 << karma::string << lit(
" float") << eol
365 << *(karma::float_ << lit(
' ') << karma::float_
366 << lit(
' ') << karma::float_
369 vector_data_double %= lit(
"VECTORS ")
370 << karma::string << lit(
" double") << eol
371 << *(karma::double_ << lit(
' ') << karma::double_
372 << lit(
' ') << karma::double_
377 (*(scalar_data_char | scalar_data_uchar | scalar_data_short |
378 scalar_data_ushort | scalar_data_int | scalar_data_unsigned_int |
379 scalar_data_long | scalar_data_unsigned_long | scalar_data_float |
380 scalar_data_double | vector_data_float | vector_data_double));
382 root %= lit(
"# vtk DataFile Version 3.0")
383 << karma::eol << karma::string << karma::eol << karma::stream
384 << karma::eol << unstructured_grid << karma::eol << field_data
385 << lit(
"POINT_DATA ") << lit(size(at_c<0>(at_c<2>(_val)))) << eol
386 << attributes << lit(
"CELL_DATA ")
387 << lit(size(at_c<1>(at_c<2>(_val)))) << eol << attributes;
390 karma::rule<Iterator, VtkFile()> root;
391 karma::rule<Iterator, VtkFile()> header;
392 karma::rule<Iterator, VtkFile::UnstructuredGrid()> unstructured_grid;
393 karma::rule<Iterator, std::vector<Eigen::Vector3f>()> points;
394 karma::rule<Iterator, std::vector<Eigen::Vector3f>()> points_vector;
395 karma::rule<Iterator, std::vector<std::vector<unsigned int>>> cells;
396 karma::rule<Iterator, std::vector<unsigned int>> cells_vector;
397 karma::rule<Iterator, std::vector<VtkFile::CellType>> cell_types;
398 karma::rule<Iterator, VtkFile::CellType> cell_type;
401 karma::rule<Iterator, VtkFile::FieldDataArray<int>()> field_data_int;
402 karma::rule<Iterator, VtkFile::FieldDataArray<float>()> field_data_float;
403 karma::rule<Iterator, VtkFile::FieldDataArray<double>()> field_data_double;
407 karma::rule<Iterator, VtkFile::ScalarData<char>()> scalar_data_char;
408 karma::rule<Iterator, VtkFile::ScalarData<unsigned char>()> scalar_data_uchar;
409 karma::rule<Iterator, VtkFile::ScalarData<short>()> scalar_data_short;
410 karma::rule<Iterator, VtkFile::ScalarData<unsigned short>()>
412 karma::rule<Iterator, VtkFile::ScalarData<int>()> scalar_data_int;
413 karma::rule<Iterator, VtkFile::ScalarData<unsigned int>()>
414 scalar_data_unsigned_int;
415 karma::rule<Iterator, VtkFile::ScalarData<long>()> scalar_data_long;
416 karma::rule<Iterator, VtkFile::ScalarData<unsigned long>()>
417 scalar_data_unsigned_long;
418 karma::rule<Iterator, VtkFile::ScalarData<float>()> scalar_data_float;
419 karma::rule<Iterator, VtkFile::ScalarData<double>()> scalar_data_double;
420 karma::rule<Iterator, VtkFile::VectorData<float>()> vector_data_float;
421 karma::rule<Iterator, VtkFile::VectorData<double>()> vector_data_double;
425void ValidateVtkFile(
const VtkFile& vtk_file) {
427 if (vtk_file.header.find(
'\n') != std::string::npos) {
428 throw base::LfException(
429 "Header of vtk file contains a new line character, this is not "
432 if (vtk_file.header.size() > 256) {
433 throw base::LfException(
"Header of vtk file is longer than 256 characters");
435 if (vtk_file.unstructured_grid.cell_types.size() !=
436 vtk_file.unstructured_grid.cells.size()) {
437 throw base::LfException(
438 "Mismatch of size of cell_types and cells in VtkFile.");
440 for (
const auto& d : vtk_file.point_data) {
441 boost::apply_visitor(
443 if (e.data.size() != vtk_file.unstructured_grid.points.size()) {
444 throw base::LfException(
"Length of dataset " + e.name +
445 " does not match the number of points.");
450 for (
const auto& d : vtk_file.cell_data) {
451 boost::apply_visitor(
453 if (e.data.size() != vtk_file.unstructured_grid.cells.size()) {
454 throw base::LfException(
"Length of dataset " + e.name +
455 " does not match the number of cells.");
463unsigned int NumAuxNodes(base::RefEl ref_el,
unsigned char order) {
470 return (2 - 3 * order + order * order) / 2;
472 return (order - 1) * (order - 1);
474 LF_VERIFY_MSG(
false,
"RefElType " << ref_el <<
" not supported.");
479Eigen::MatrixXd AuxNodesSegment(
unsigned char order) {
480 LF_ASSERT_MSG(order > 1,
"order <= 1");
481 return Eigen::VectorXd::LinSpaced(order - 1, 1. / (order),
482 (order - 1.) / order)
488Eigen::MatrixXd AuxNodesTria(
unsigned char order) {
494 result << 1. / 3., 1. / 3.;
495 }
else if (order == 4) {
496 result << 0.25, 0.5, 0.25, 0.25, 0.25, 0.5;
499 result(0, 0) = 1. / order;
500 result(1, 0) = 1. / order;
501 result(0, 1) = (order - 2.) / order;
502 result(1, 1) = 1. / order;
503 result(0, 2) = 1. / order;
504 result(1, 2) = (order - 2.) / order;
507 auto segment_points = AuxNodesSegment(order - 3).eval();
509 result.block(0, 3, 2, order - 4) =
510 Eigen::Vector2d::UnitX() * segment_points * (order - 3.) / order +
511 Eigen::Vector2d(1. / order, 1. / order) *
512 Eigen::MatrixXd::Ones(1, order - 4);
513 result.block(0, order - 1, 2, order - 4) =
514 Eigen::Vector2d(-1, 1) * (order - 3.) / order * segment_points +
515 Eigen::Vector2d((order - 2.) / order, 1. / order) *
516 Eigen::MatrixXd::Ones(1, order - 4);
517 result.block(0, 2 * order - 5, 2, order - 4) =
518 -Eigen::Vector2d::UnitY() * (order - 3.) / order * segment_points +
519 Eigen::Vector2d(1. / order, (order - 2.) / order) *
520 Eigen::MatrixXd::Ones(1, order - 4);
524 auto points = AuxNodesTria(order - 3);
525 result.block(0, 3 * order - 9, 2, points.cols()) =
526 AuxNodesTria(order - 3) * (order - 6.) / order +
527 Eigen::Vector2d(2. / order, 2. / order) *
528 Eigen::MatrixXd::Ones(1, points.cols());
535Eigen::MatrixXd AuxNodesQuad(
unsigned char order) {
538 Eigen::VectorXd::LinSpaced(order - 1, 1. / order, (order - 1.) / order)
540 .replicate(1, order - 1);
542 Eigen::kroneckerProduct(Eigen::VectorXd::LinSpaced(order - 1, 1. / order,
543 (order - 1.) / order),
544 Eigen::VectorXd::Constant(order - 1, 1.))
552 std::ofstream file(filename, std::ios_base::out | std::ios_base::binary |
553 std::ios_base::trunc);
554 ValidateVtkFile(vtk_file);
556 if (!file.is_open()) {
560 karma::ostream_iterator<char> outit(file);
564 const VtkGrammar<
decltype(outit),
true> grammar{};
565 result = karma::generate(outit, grammar, vtk_file);
567 const VtkGrammar<
decltype(outit),
false> grammar{};
568 result = karma::generate(outit, grammar, vtk_file);
579 std::string filename,
dim_t codim,
unsigned char order)
580 : mesh_(std::move(mesh)),
581 filename_(std::move(filename)),
585 auto dim_mesh = mesh_->DimMesh();
586 auto dim_world = mesh_->DimWorld();
587 LF_ASSERT_MSG(dim_world > 0 && dim_world <= 4,
588 "VtkWriter supports only dim_world = 1,2 or 3");
589 LF_ASSERT_MSG(codim >= 0 && codim < dim_mesh,
"codim out of bounds.");
590 LF_ASSERT_MSG(order > 0 && order <= 10,
"order must lie in [1,10]");
600 unsigned int numNodes = mesh_->NumEntities(dim_mesh);
603 if (ref_el.Dimension() > dim_mesh - codim_) {
606 numNodes += NumAuxNodes(ref_el, order_) * mesh_->NumEntities(ref_el);
610 vtk_file_.unstructured_grid.points.resize(numNodes);
611 const Eigen::Matrix<double, 0, 1> zero;
612 for (
const auto* p : mesh_->Entities(dim_mesh)) {
613 auto index = mesh_->Index(*p);
614 Eigen::Vector3f coord;
615 if (dim_world == 1) {
616 coord(0) =
static_cast<float>(p->Geometry()->Global(zero)(0));
619 }
else if (dim_world == 2) {
620 coord.topRows<2>() = p->Geometry()->Global(zero).cast<
float>();
623 coord = p->Geometry()->Global(zero).cast<
float>();
626 vtk_file_.unstructured_grid.points[index] = std::move(coord);
631 auto index_offset = mesh_->NumEntities(dim_mesh);
632 for (
char cd =
static_cast<char>(dim_mesh - 1);
633 cd >=
static_cast<char>(codim); --cd) {
634 for (
const auto* e : mesh_->Entities(cd)) {
635 auto ref_el = e->RefEl();
640 Eigen::MatrixXf coords(3, NumAuxNodes(ref_el, order));
642 if (dim_world == 1) {
644 e->Geometry()->Global(aux_nodes_[ref_el.Id()]).cast<
float>();
645 coords.row(1).setZero();
646 coords.row(2).setZero();
647 }
else if (dim_world == 2) {
649 e->Geometry()->Global(aux_nodes_[ref_el.Id()]).cast<float>();
650 coords.row(2).setZero();
652 coords = e->Geometry()->Global(aux_nodes_[ref_el.Id()]).cast<float>();
654 for (Eigen::Index i = 0; i < coords.cols(); ++i) {
655 vtk_file_.unstructured_grid.points[index_offset + i] = coords.col(i);
657 aux_node_offset_[cd](*e) = index_offset;
658 index_offset += coords.cols();
664 std::array<unsigned int, 5> num_nodes{};
675 vtk_file_.unstructured_grid.cells.resize(mesh_->NumEntities(codim));
676 vtk_file_.unstructured_grid.cell_types.resize(mesh_->NumEntities(codim));
678 for (
const auto* e : mesh_->Entities(codim)) {
679 auto index = mesh_->Index(*e);
680 auto ref_el = e->RefEl();
681 auto& node_indices = vtk_file_.unstructured_grid.cells[index];
682 node_indices.reserve(num_nodes[ref_el.Id()]);
685 for (
const auto* p : e->SubEntities(dim_mesh - codim)) {
686 node_indices.push_back(mesh_->Index(*p));
690 auto addSegmentNodes = [&](
const mesh::Entity& s,
bool invert) {
691 auto start_index = aux_node_offset_[dim_mesh - 1](s);
693 for (
unsigned int i = 0; i < points_per_segment; ++i) {
694 node_indices.push_back(start_index + i);
697 for (
int i =
static_cast<int>(points_per_segment - 1); i >= 0; --i) {
698 node_indices.push_back(start_index + i);
704 addSegmentNodes(*e, false);
707 auto iterator = e->SubEntities(1).begin();
708 auto o_iterator = e->RelativeOrientations().begin();
709 addSegmentNodes(**iterator,
713 addSegmentNodes(**iterator,
717 addSegmentNodes(**iterator,
722 auto iterator = e->SubEntities(1).begin();
723 auto o_iterator = e->RelativeOrientations().begin();
724 addSegmentNodes(**iterator,
728 addSegmentNodes(**iterator,
732 addSegmentNodes(**iterator,
736 addSegmentNodes(**iterator,
741 LF_VERIFY_MSG(
false,
"something is wrong");
745 if (dim_mesh - codim > 1) {
746 auto offset = aux_node_offset_[0](*e);
747 for (
unsigned i = 0; i < NumAuxNodes(e->RefEl(), order); ++i) {
748 node_indices.push_back(offset + i);
754 vtk_file_.unstructured_grid.cell_types[index] =
757 vtk_file_.unstructured_grid.cell_types[index] =
760 vtk_file_.unstructured_grid.cell_types[index] =
766 vtk_file_.unstructured_grid.cell_types[index] =
767 VtkFile::CellType::VTK_LAGRANGE_CURVE;
770 vtk_file_.unstructured_grid.cell_types[index] =
771 VtkFile::CellType::VTK_LAGRANGE_TRIANGLE;
774 vtk_file_.unstructured_grid.cell_types[index] =
775 VtkFile::CellType::VTK_LAGRANGE_QUADRILATERAL;
778 LF_VERIFY_MSG(
false,
"Unsupported RefElType " << ref_el);
786 unsigned char undefined_value) {
792 char undefined_value) {
798 unsigned int undefined_value) {
804 int undefined_value) {
810 float undefined_value) {
816 double undefined_value) {
821 const std::string& name,
823 const Eigen::Vector2d& undefined_value) {
824 WriteVectorPointData<2, double>(name, mds, undefined_value);
828 const std::string& name,
830 const Eigen::Vector2f& undefined_value) {
831 WriteVectorPointData<2, float>(name, mds, undefined_value);
835 const std::string& name,
837 const Eigen::Vector3d& undefined_value) {
838 WriteVectorPointData<3, double>(name, mds, undefined_value);
842 const std::string& name,
844 const Eigen::Vector3f& undefined_value) {
845 WriteVectorPointData<3, float>(name, mds, undefined_value);
849 const std::string& name,
851 const Eigen::VectorXd& undefined_value) {
852 WriteVectorPointData<Eigen::Dynamic, double>(name, mds, undefined_value);
856 const std::string& name,
858 const Eigen::VectorXf& undefined_value) {
859 WriteVectorPointData<Eigen::Dynamic, float>(name, mds, undefined_value);
864 unsigned char undefined_value) {
870 char undefined_value) {
876 unsigned int undefined_value) {
882 int undefined_value) {
888 float undefined_value) {
894 double undefined_value) {
899 const std::string& name,
901 const Eigen::Vector2d& undefined_value) {
902 WriteVectorCellData<2, double>(name, mds, undefined_value);
906 const std::string& name,
908 const Eigen::Vector2f& undefined_value) {
909 WriteVectorCellData<2, float>(name, mds, undefined_value);
913 const std::string& name,
915 const Eigen::Vector3d& undefined_value) {
916 WriteVectorCellData<3, double>(name, mds, undefined_value);
920 const std::string& name,
922 const Eigen::Vector3f& undefined_value) {
923 WriteVectorCellData<3, float>(name, mds, undefined_value);
927 const std::string& name,
929 const Eigen::VectorXd& undefined_value) {
930 WriteVectorCellData<Eigen::Dynamic, double>(name, mds, undefined_value);
934 const std::string& name,
936 const Eigen::VectorXf& undefined_value) {
937 WriteVectorCellData<Eigen::Dynamic, float>(name, mds, undefined_value);
941 std::vector<int> data) {
946 std::vector<float> data) {
951 std::vector<double> data) {
959 LF_ASSERT_MSG(
order_ == 1,
960 "WritePointData accepts MeshDataSets only if order = 1. For "
961 "order > 1 you have to provide MeshFunctions.");
966 for (
const auto* p :
mesh_->Entities(
mesh_->DimMesh())) {
968 data.data[
mesh_->Index(*p)] = mds(*p);
970 data.data[
mesh_->Index(*p)] = undefined_value;
976template <
int ROWS,
class T>
978 const std::string& name,
980 const Eigen::Matrix<T, ROWS, 1>& undefined_value) {
981 LF_ASSERT_MSG(
order_ == 1,
982 "WritePointData accepts MeshDataSets only if order = 1. For "
983 "order > 1 you have to provide MeshFunctions.");
988 Eigen::Matrix<T, 3, 1> undefined_value_padded;
989 PadWithZeros<ROWS, T>(undefined_value_padded, undefined_value);
991 for (
const auto* p :
mesh_->Entities(
mesh_->DimMesh())) {
992 if (mds.DefinedOn(*p)) {
993 PadWithZeros<ROWS, T>(data.data[
mesh_->Index(*p)], mds(*p));
995 data.data[
mesh_->Index(*p)] = undefined_value_padded;
1004 T undefined_value) {
1011 data.data[
mesh_->Index(*e)] = mds(*e);
1013 data.data[
mesh_->Index(*e)] = undefined_value;
1019template <
int ROWS,
class T>
1021 const std::string& name,
1023 const Eigen::Matrix<T, ROWS, 1>& undefined_value) {
1028 Eigen::Matrix<T, 3, 1> undefined_value_padded;
1029 PadWithZeros<ROWS, T>(undefined_value_padded, undefined_value);
1032 if (mds.DefinedOn(*p)) {
1033 PadWithZeros<ROWS, T>(data.data[
mesh_->Index(*p)], mds(*p));
1035 data.data[
mesh_->Index(*p)] = undefined_value_padded;
1045 vtk_data.
name = name;
1046 vtk_data.data = std::move(data);
result_type operator()(T &&arg)
VariantVisitor(LAMBDA lambda)
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.
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.
Represents one set of attribute data (can be attached to points or cells)
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
Format format
The format of the file.
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....
void WriteScalarCellData(const std::string &name, const mesh::utils::MeshDataSet< T > &mds, T undefined_value)
void WriteFieldData(const std::string &name, std::vector< T > data)
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_
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.
virtual bool DefinedOn(const Entity &e) const =0
Does the dataset store information with this entity?
Mesh input (from file) and output (in various formats) facilities.
std::ostream & operator<<(std::ostream &stream, GMshFileV2::ElementType et)
Output the element type onto the console:
void WriteToFile(const VtkFile &vtk_file, const std::string &filename)