Wildmeshing Toolkit
Loading...
Searching...
No Matches
to_points.cpp
Go to the documentation of this file.
1#include "to_points.hpp"
2
3
4#include <bitset>
5
6#include <SimpleBVH/BVH.hpp>
7#include <wmtk/PointMesh.hpp>
11
12namespace wmtk::components {
13
14using namespace simplex;
15
16std::shared_ptr<PointMesh> to_points(
17 const Mesh& mesh,
18 const attribute::MeshAttributeHandle& pts_attr,
19 const ToPtsOptions& options,
20 const std::string& output_pos_attr_name)
21{
22 const auto pts_acc = mesh.create_const_accessor<double>(pts_attr);
23 const auto vs = mesh.get_all(PrimitiveType::Vertex);
24
25 Eigen::AlignedBox<double, Eigen::Dynamic> bbox(pts_acc.dimension());
26
27 int64_t box_size = std::pow(2, pts_acc.dimension());
28
29 Eigen::MatrixXd pts(vs.size() + (options.add_box ? box_size : 0), pts_acc.dimension());
30
31 int64_t index = 0;
32 for (const auto& v : vs) {
33 pts.row(index) << pts_acc.const_vector_attribute(v).transpose();
34 bbox.extend(pts.row(index).transpose());
35 ++index;
36 }
37
38 wmtk::logger().info(
39 "Input bbox max {}, min {}, diag {}, potential target edge length: {}",
40 bbox.max(),
41 bbox.min(),
42 bbox.diagonal().norm(),
43 bbox.diagonal().norm() * options.box_scale);
44
45 if (options.add_box && !options.add_grid) {
46 auto center = bbox.center();
47 auto r = bbox.diagonal() / 2.;
48
49 bbox.min() = center - options.box_scale * r;
50 bbox.max() = center + options.box_scale * r;
51
52 for (int64_t d = 0; d < box_size; ++d) {
53 std::bitset<32> bits = d;
54 for (int64_t i = 0; i < pts.cols(); ++i) {
55 pts(index, i) = bits[i] == 0 ? bbox.min()[i] : bbox.max()[i];
56 }
57 ++index;
58 }
59 }
60
61 if (options.add_grid) {
62 const double bbox_diag = bbox.diagonal().norm();
63 const auto r = Eigen::VectorXd::Ones(pts_acc.dimension()) * bbox_diag;
64 const double grid_spacing = options.grid_spacing;
65
66
67 bbox.min() -= options.box_scale * r;
68 bbox.max() += options.box_scale * r;
69 // Eigen::VectorXi res = (diag / grid_spacing).cast<int>();
70
71 // // TODO: remove the hack
72 // // hack grid spacing as relative
73 const Eigen::VectorXi res = (bbox.diagonal() / (bbox_diag * grid_spacing)).cast<int>() +
74 Eigen::VectorXi::Ones(pts_acc.dimension());
75
76 Eigen::MatrixXd background_V;
77
78
79 if (pts.cols() == 3) {
80 auto v_index = [&res](int i, int j, int k) {
81 return k * (res[0] + 1) * (res[1] + 1) + j * (res[0] + 1) + i;
82 };
83
84 background_V.resize((res[0] + 1) * (res[1] + 1) * (res[2] + 1), 3);
85
86 for (int k = 0; k <= res[2]; ++k) {
87 for (int j = 0; j <= res[1]; ++j) {
88 for (int i = 0; i <= res[0]; ++i) {
89 const Eigen::Vector3d iii(i, j, k);
90 const Eigen::Vector3d ttmp = bbox.diagonal().array() * iii.array();
91 background_V.row(v_index(i, j, k)) =
92 bbox.min().array() + ttmp.array() / res.cast<double>().array();
93 }
94 }
95 }
96 }
97
98 wmtk::logger().info(
99 "Grid bbox max {}, min {}, diag {}, potential target edge length: {}",
100 bbox.max(),
101 bbox.min(),
102 bbox.diagonal().norm(),
103 bbox.diagonal().norm() * options.box_scale);
104 // else 2d and 1d
105
106 if (options.min_dist >= 0) {
107 int64_t count = 0;
108 int64_t findex = 0;
109
110 const std::vector<Tuple>& facest = mesh.get_all(mesh.top_simplex_type());
111 const int64_t dim = int64_t(mesh.top_simplex_type()) + 1;
112
113 Eigen::MatrixXd vertices(dim * facest.size(), pts_acc.dimension());
114 Eigen::MatrixXi faces(facest.size(), dim);
115
116 for (const auto& f : facest) {
118 mesh,
119 simplex::Simplex(mesh, mesh.top_simplex_type(), f),
121
122 assert(tmp.size() == dim);
123 for (int64_t j = 0; j < tmp.size(); ++j) {
124 auto p = pts_acc.const_vector_attribute(tmp[j]);
125 faces(findex, j) = count;
126 vertices.row(dim * findex + j) = p;
127
128 ++count;
129 }
130 ++findex;
131 }
132
133 SimpleBVH::BVH bvh;
134 bvh.init(vertices, faces, 1e-10);
135
136 const double min_dist =
137 options.min_dist > 0 ? (options.min_dist * options.min_dist * bbox_diag * bbox_diag)
138 : (bbox_diag * bbox_diag * grid_spacing * grid_spacing / 4);
139 std::vector<Eigen::VectorXd> good;
140 SimpleBVH::VectorMax3d nearest_point;
141 for (int64_t i = 0; i < background_V.rows(); ++i) {
142 double sq_dist;
143 bvh.nearest_facet(background_V.row(i), nearest_point, sq_dist);
144 if (sq_dist >= min_dist) good.emplace_back(background_V.row(i));
145 }
146 int64_t current_size = pts.rows();
147 pts.conservativeResize(current_size + good.size(), pts.cols());
148 for (const auto& v : good) {
149 pts.row(current_size) = v.transpose();
150 ++current_size;
151 }
152 } else {
153 // TODO
154 }
155 }
156
157 // remove duplicates
158 int64_t old_size = pts.rows();
159 auto remove_duplicated_vertices = [](Eigen::MatrixXd& P) -> Eigen::MatrixXd {
160 std::vector<Eigen::VectorXd> vec;
161 for (int64_t i = 0; i < P.rows(); ++i) vec.push_back(P.row(i));
162
163 std::sort(vec.begin(), vec.end(), [](Eigen::VectorXd const& p1, Eigen::VectorXd const& p2) {
164 return (p1(0) < p2(0)) || (p1(0) == p2(0) && p1(1) < p2(1)) ||
165 (p1(0) == p2(0) && p1(1) == p2(1) && p1(2) < p2(2));
166 });
167
168 auto it = std::unique(vec.begin(), vec.end());
169 vec.resize(std::distance(vec.begin(), it));
170
171 Eigen::MatrixXd new_P(vec.size(), P.cols());
172 for (int64_t i = 0; i < vec.size(); ++i) {
173 new_P.row(i) = vec[i];
174 }
175
176 return new_P;
177 };
178
179 if (options.remove_duplicates) {
180 pts = remove_duplicated_vertices(pts);
181 wmtk::logger().info("removed {} duplicated vertices", pts.rows() - old_size);
182 }
183
184 wmtk::logger().info("generated {} vertices", pts.rows());
185 std::shared_ptr<PointMesh> pts_mesh = std::make_shared<PointMesh>(pts.rows());
186
187 mesh_utils::set_matrix_attribute(pts, output_pos_attr_name, PrimitiveType::Vertex, *pts_mesh);
188
189 return pts_mesh;
190}
191
192} // namespace wmtk::components
const attribute::Accessor< T, Mesh, D > create_const_accessor(const attribute::MeshAttributeHandle &handle) const
std::vector< Tuple > get_all(PrimitiveType type) const
Generate a vector of Tuples from global vertex/edge/triangle/tetrahedron index.
Definition Mesh.cpp:18
PrimitiveType top_simplex_type() const
Definition Mesh.hpp:982
std::shared_ptr< PointMesh > to_points(const Mesh &mesh, const attribute::MeshAttributeHandle &pts_attr, const ToPtsOptions &options, const std::string &output_pos_attr_name)
Definition to_points.cpp:16
attribute::MeshAttributeHandle set_matrix_attribute(const Mat &data, const std::string &name, const PrimitiveType &type, Mesh &mesh)
Definition mesh_utils.hpp:9
std::vector< Tuple > vertices(const Mesh &m, const Simplex &simplex)
std::vector< Tuple > faces_single_dimension_tuples(const Mesh &mesh, const Simplex &simplex, const PrimitiveType face_type)
SimplexCollection faces(const Mesh &mesh, const Simplex &simplex, const bool sort_and_clean)
Returns all faces of a simplex.
Definition faces.cpp:10
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:58