15 #include "mshio/MshSpec.h"
16 #include "predicates.h"
25 extract_vertices<DIM>();
26 extract_simplex_elements<DIM>();
30 for (
int dim = 3; dim > 0; --dim) {
39 const std::filesystem::path& filename,
40 const bool ignore_z_if_zero,
41 const std::vector<std::string>& attrs)
43 m_spec = mshio::load_msh(filename.string());
48 split_attrs.back() = attrs;
55 const std::filesystem::path& filename,
56 const int64_t embedded_dimension,
57 const std::vector<std::vector<std::string>>& extra_attributes)
59 m_spec = mshio::load_msh(filename.string());
67 const std::filesystem::path& filename,
68 const int64_t embedded_dimension)
70 m_spec = mshio::load_msh(filename.string());
80 for (
const auto& block :
m_spec.nodes.entity_blocks) {
81 if (block.entity_dim == DIM) {
90 for (
const auto& block :
m_spec.elements.entity_blocks) {
91 if (block.entity_dim == DIM) {
101 if (block !=
nullptr) {
102 return block->num_nodes_in_block;
111 if (block !=
nullptr) {
112 return block->num_elements_in_block;
122 assert(block !=
nullptr);
124 return block->entity_dim;
131 assert(block !=
nullptr);
133 const size_t num_vertices = block->num_nodes_in_block;
134 const auto embedded_dimension =
136 V.resize(num_vertices, embedded_dimension);
140 assert(num_vertices > 0);
142 assert(embedded_dimension <= 3);
143 assert(embedded_dimension > 0);
146 using ConstMapType =
typename Vector::ConstMapType;
148 const size_t tag_offset = block->tags.front();
149 for (
size_t i = 0; i < num_vertices; i++) {
150 size_t tag = block->tags[i] - tag_offset;
151 ConstMapType vec(block->data.data() + i * 3);
152 V.row(tag) = vec.head(embedded_dimension).transpose();
156 const double max =
V.col(2).array().maxCoeff();
157 const double min =
V.col(2).array().minCoeff();
159 if (min == 0 && max == 0) {
160 V =
V.leftCols(2).eval();
170 if (element_block ==
nullptr)
return;
171 assert(vertex_block !=
nullptr);
173 const size_t num_elements = element_block->num_elements_in_block;
174 if (num_elements == 0)
return;
175 assert(vertex_block->num_nodes_in_block != 0);
177 S.resize(num_elements, DIM + 1);
179 const size_t vert_tag_offset = vertex_block->tags.front();
180 const size_t elem_tag_offset = element_block->data.front();
181 for (
size_t i = 0; i < num_elements; i++) {
182 const size_t tag = element_block->data[i * (DIM + 2)] - elem_tag_offset;
183 assert(tag < num_elements);
184 const size_t* element = element_block->data.data() + i * (DIM + 2) + 1;
187 using ConstMapType =
typename Vector::ConstMapType;
189 ConstMapType element_vec(element);
190 const auto vertex_tag_offset_vec = Vector::Constant(vert_tag_offset);
192 S.row(tag) = (element_vec - vertex_tag_offset_vec).transpose().template cast<int64_t>();
199 auto ret = std::make_shared<MeshType>();
284 const std::string& attr_name,
291 const size_t tag_offset = element_block->data.front();
293 for (
const auto& data :
m_spec.element_data) {
294 if (data.header.string_tags.front() != attr_name) {
297 assert(data.entries.size() == tuples.size());
299 const int64_t attr_dim = data.header.int_tags[1];
304 for (
size_t i = 0; i < tuples.size(); ++i) {
305 const auto& entry = data.entries[i];
306 const size_t tag = entry.tag - tag_offset;
307 assert(tag < element_block->num_elements_in_block);
310 acc.scalar_attribute(tuples[i]) = entry.data[0];
312 acc.vector_attribute(tuples[i]) =
313 Eigen::Map<const Eigen::VectorXd>(entry.data.data(), entry.data.size());
321 const std::optional<std::vector<std::vector<std::string>>>& extra_attributes_opt)
322 -> std::shared_ptr<Mesh>
324 std::shared_ptr<Mesh> res;
325 switch (get_mesh_dimension()) {
326 case 3: res = generateT<3>();
break;
327 case 2: res = generateT<2>();
break;
328 case 1: res = generateT<1>();
break;
330 case 0: res = std::make_shared<PointMesh>();
332 assert(V.size() > 0);
334 if (extra_attributes_opt.has_value()) {
335 const auto& extra_attributes = extra_attributes_opt.value();
337 std::set<std::string> valid_attr_names;
338 for (
const auto& data : m_spec.element_data) {
341 if (data.header.int_tags[2] != get_num_simplex_elements(3)) {
344 valid_attr_names.insert(data.header.string_tags.front());
346 for (
size_t dim = 0; dim < extra_attributes.size(); ++dim) {
347 const auto& attr_names = extra_attributes[dim];
348 for (
const std::string& attr : attr_names) {
349 if (valid_attr_names.count(attr) == 0) {
365 return construct<DIM>();
373 void MshReader::validate<3>()
375 assert(
V.cols() == 3);
380 auto Vs =
V(
S.row(0), Eigen::all);
385 S.col(0).swap(
S.col(1));
387 "Input tet orientation is inverted, swapping col 0 and 1 of TV matirx.");
393 for (int64_t i = 0; i <
S.rows(); i++) {
395 auto Vs =
V(row, Eigen::all);
399 std::swap(row.x(), row.y());
401 }
else if (orient == 0) {
409 for (int64_t i = 0; i <
S.rows(); i++) {
411 auto Vs =
V(row, Eigen::all);
421 }
else if (orient == 0) {
attribute::MeshAttributeHandle register_attribute(const std::string &name, PrimitiveType type, int64_t size, bool replace=false, T default_value=T(0))
std::vector< Tuple > get_all(PrimitiveType type) const
Generate a vector of Tuples from global vertex/edge/triangle/tetrahedron index.
attribute::Accessor< T, Mesh, D > create_accessor(const attribute::MeshAttributeHandle &handle)
PrimitiveType top_simplex_type() const
const mshio::ElementBlock * get_simplex_element_block(int DIM) const
auto construct() -> std::shared_ptr< wmtk::utils::mesh_type_from_dimension_t< DIM >>
size_t get_num_simplex_elements(int DIM) const
int get_mesh_dimension() const
int get_embedded_dimension() const
static const int64_t AUTO_EMBEDDED_DIMENSION
std::shared_ptr< Mesh > read(const std::filesystem::path &filename, const bool ignore_z_if_zero, const std::vector< std::string > &extra_facet_attributes={})
int64_t m_embedded_dimension
size_t get_num_vertices(int DIM) const
const mshio::NodeBlock * get_vertex_block(int DIM) const
auto generateT() -> std::shared_ptr< wmtk::utils::mesh_type_from_dimension_t< DIM >>
void extract_simplex_elements()
std::shared_ptr< Mesh > generate(const std::optional< std::vector< std::vector< std::string >>> &extra_attributes={})
void extract_element_attribute(wmtk::Mesh &m, const std::string &attr_name, PrimitiveType pt)
attribute::MeshAttributeHandle set_matrix_attribute(const Mat &data, const std::string &name, const PrimitiveType &type, Mesh &mesh)
typename mesh_type_from_dimension< DIM >::type mesh_type_from_dimension_t
int wmtk_orient3d(const Eigen::Ref< const Eigen::Vector3< Rational >> &p0, const Eigen::Ref< const Eigen::Vector3< Rational >> &p1, const Eigen::Ref< const Eigen::Vector3< Rational >> &p2, const Eigen::Ref< const Eigen::Vector3< Rational >> &p3)
constexpr PrimitiveType get_primitive_type_from_id(int8_t id)
Get the primitive type corresponding to its unique integer id.
void log_and_throw_error(const std::string &msg)
Eigen::Matrix< T, R, 1 > Vector
constexpr int8_t get_primitive_type_id(PrimitiveType t)
Get a unique integer id corresponding to each primitive type.
spdlog::logger & logger()
Retrieves the current logger.