ReUseX  0.0.1
3D Point Cloud Processing for Building Reuse
Loading...
Searching...
No Matches
CellComplexRoomProbabilites.hpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2025 Povl Filip Sonne-Frederiksen
2//
3// SPDX-License-Identifier: GPL-3.0-or-later
4#pragma once
5#include "reusex/core/logging.hpp"
6#include "reusex/core/processing_observer.hpp"
7
8#include <Eigen/Core>
9#include <embree4/rtcore.h>
10#include <fmt/color.h>
11#include <fmt/core.h>
12#include <fmt/ranges.h>
13#include <pcl/filters/uniform_sampling.h>
14
15#include <cmath>
16#include <vector>
17
18// Spherical Fibonacci
19static std::vector<Eigen::Vector3d> sampleSphericalFibonacci(size_t N) {
20 std::vector<Eigen::Vector3d> out;
21 out.reserve(N);
22 if (N == 0)
23 return out;
24 constexpr double golden_angle =
25 M_PI * (3.0 - std::sqrt(5.0)); // ~2.3999632297
26
27 for (size_t i = 0; i < N; ++i) {
28 double t = (double)i;
29 double z = 1.0 - 2.0 * t / (double)(N - 1); // maps to [1, -1] if N>1
30 double phi = golden_angle * t;
31 double r = std::sqrt(std::max(0.0, 1.0 - z * z));
32 double x = std::cos(phi) * r;
33 double y = std::sin(phi) * r;
34 out.emplace_back(x, y, z);
35 }
36 return out;
37}
38
39namespace ReUseX::geometry {
40template <typename PointT, typename PointN, typename PointL>
42 pcl::PointCloud<PointT>::ConstPtr cloud_,
43 pcl::PointCloud<PointN>::ConstPtr normals_,
44 pcl::PointCloud<PointL>::ConstPtr labels_, const double grid_size) -> void {
45
46 using RTCVertex = struct RTCVertex {
47 float x, y, z, r; // x, y, z coordinates and radius
48 };
49
50 using RTCNormal = struct RTCNormal {
51 float x, y, z; // x, y, z components of the normal vector
52 };
53
54 using RTCData = struct RTCData {
55 size_t label_index;
56 };
57
58 ReUseX::core::trace("calling compute_room_probabilities");
59
60 assert(cloud_->size() == labels_->size() &&
61 "Point cloud and label cloud must have the same size");
62 assert(cloud_->size() == normals_->size() &&
63 "Point cloud and normal cloud must have the same size");
64
65 const double radius_ = grid_size * M_SQRT2;
66
68 "Using ray tracing for {} points with grid size {:.3} and radius {:.3}",
69 cloud_->size(), grid_size, radius_);
70
72 "c:room_probabilities")
73 .first;
74
75 // Get unique labels
76 std::vector<unsigned int> labels;
77 labels.reserve(labels_->points.size());
78
79 for (const auto &p : labels_->points)
80 if (p.label != static_cast<unsigned int>(-1))
81 labels.push_back(p.label);
82
83 std::sort(labels.begin(), labels.end());
84 labels.erase(std::unique(labels.begin(), labels.end()), labels.end());
85 ReUseX::core::debug("Room labels: {}", fmt::join(labels, ", "));
86
87 for (auto cit = this->cells_begin(); cit != this->cells_end(); ++cit)
88 c_rp[*cit] = std::vector<double>(labels.size() + 1, 0.0);
89 this->n_rooms = labels.size();
90 ReUseX::core::debug("Number of rooms (labels): {}", labels.size());
91
92 unsigned int old_mxcsr = _mm_getcsr(); // save current flagsq
93 _MM_SET_FLUSH_ZERO_MODE(_MM_FLUSH_ZERO_ON);
94 _MM_SET_DENORMALS_ZERO_MODE(_MM_DENORMALS_ZERO_ON);
95
96 ReUseX::core::trace("Creating device and scene for ray tracing");
97 RTCDevice device_ = rtcNewDevice("verbose=0"); // 0-3
98 RTCScene scene_ = rtcNewScene(device_);
99 rtcSetSceneFlags(scene_, RTC_SCENE_FLAG_ROBUST);
100 rtcSetSceneBuildQuality(scene_, RTC_BUILD_QUALITY_HIGH);
101 assert(device_ != nullptr && "Error creating Embree device");
102 assert(scene_ != nullptr && "Error creating Embree scene");
103 rtcSetDeviceErrorFunction(
104 device_,
105 [](void *, RTCError err, const char *str) {
106 ReUseX::core::error("Embree error {}:", str);
107 switch (err) {
108 case RTC_ERROR_NONE:
109 break;
110 case RTC_ERROR_UNKNOWN:
111 throw std::runtime_error("Embree: An unknown error has occurred.");
112 break;
113 case RTC_ERROR_INVALID_ARGUMENT:
114 throw std::runtime_error(
115 "Embree: An invalid argument was specified.");
116 break;
117 case RTC_ERROR_INVALID_OPERATION:
118 throw std::runtime_error(
119 "Embree: The operation is not allowed for the specified object.");
120 break;
121 case RTC_ERROR_OUT_OF_MEMORY:
122 throw std::runtime_error("Embree: There is not enough memory left to "
123 "complete the operation.");
124 break;
125 case RTC_ERROR_UNSUPPORTED_CPU:
126 throw std::runtime_error(
127 "Embree: The CPU is not supported as it does not support SSE2.");
128 break;
129 case RTC_ERROR_CANCELLED:
130 throw std::runtime_error(
131 "Embree: The operation got cancelled by an Memory Monitor "
132 "Callback or Progress Monitor Callback function.");
133 break;
134 default:
135 throw std::runtime_error("Embree: An invalid error has occurred.");
136 break;
137 }
138 },
139 nullptr);
140
141 // Intersect with ground plane to get line/ INFO: Create Embree scene
142 ReUseX::core::trace("Creating scene geometry for ray tracing");
143
144 pcl::UniformSampling<PointT> us;
145 us.setInputCloud(cloud_);
146 us.setRadiusSearch(grid_size);
147
148 using RTCDataPtr = std::shared_ptr<RTCData>;
149 std::vector<RTCDataPtr> rtc_data_vec{}; // To keep data alive
150 for (size_t i = 0; i < labels.size(); ++i) {
151
152 pcl::IndicesPtr indices(new pcl::Indices);
153 for (size_t j = 0; j < cloud_->size(); ++j)
154 if (labels_->points[j].label == labels[i])
155 indices->push_back(static_cast<int>(j));
156
157 ReUseX::core::trace("Room label {}: {} points", labels[i], indices->size());
158
159 // size_t n_points_before = indices->size();
160 us.setIndices(indices);
161 us.filter(*indices);
162 // ReUseX::core::trace(
163 // "Label {}: Reduced from {} to {} points after uniform sampling",
164 // labels[i], n_points_before, indices->size());
165
166 RTCGeometry geometry_ =
167 rtcNewGeometry(device_, RTC_GEOMETRY_TYPE_ORIENTED_DISC_POINT);
168
169 RTCVertex *vb = (RTCVertex *)rtcSetNewGeometryBuffer(
170 geometry_, RTC_BUFFER_TYPE_VERTEX, 0, RTC_FORMAT_FLOAT4,
171 sizeof(RTCVertex), indices->size());
172
173 RTCNormal *nb = (RTCNormal *)rtcSetNewGeometryBuffer(
174 geometry_, RTC_BUFFER_TYPE_NORMAL, 0, RTC_FORMAT_FLOAT3,
175 sizeof(RTCNormal), indices->size());
176
177 for (size_t i = 0; i < indices->size(); ++i) {
178 const auto &p = cloud_->points[indices->at(i)];
179 const auto &n = normals_->points[indices->at(i)];
180
181 vb[i].x = static_cast<float>(p.x);
182 vb[i].y = static_cast<float>(p.y);
183 vb[i].z = static_cast<float>(p.z);
184 vb[i].r = static_cast<float>(radius_) * 1.05f;
185
186 nb[i].x = static_cast<float>(n.normal_x);
187 nb[i].y = static_cast<float>(n.normal_y);
188 nb[i].z = static_cast<float>(n.normal_z);
189 }
190
191 rtc_data_vec.emplace_back(new RTCData{i});
192 rtcSetGeometryUserData(geometry_, rtc_data_vec.back().get());
193
194 rtcCommitGeometry(geometry_);
195 rtcAttachGeometry(scene_, geometry_);
196 rtcReleaseGeometry(geometry_);
197 }
198
199 ReUseX::core::trace("Committing scene");
200 rtcCommitScene(scene_);
201
202 {
203 const auto dirs = sampleSphericalFibonacci(100);
204
205 auto observer = ReUseX::core::ProgressObserver(
206 "Computing room probabilities", this->num_cells());
207
208#pragma omp parallel for schedule(dynamic)
209 for (size_t cell_idx = 0; cell_idx < this->num_cells(); ++cell_idx) {
210 auto cit = this->cells_begin();
211 std::advance(cit, cell_idx);
212 // const size_t idx = (*this)[*cit].id;
213 // For each cell create n random rays
214 // Count the number of intersections per id nad normalize
215
216 std::vector<double> local_accum(c_rp[*cit].size(), 0.0);
217 double *accum_ptr = local_accum.data();
218 const int N = static_cast<int>(c_rp[*cit].size());
219
220#pragma omp parallel for reduction(+ : accum_ptr[ : N])
221 for (size_t i = 0; i < dirs.size(); ++i) {
222 const auto dir = dirs[i];
223 auto c = (*this)[*cit].pos;
224
225 RTCRayHit rayhit;
226 rayhit.ray.org_x = static_cast<float>(c.x());
227 rayhit.ray.org_y = static_cast<float>(c.y());
228 rayhit.ray.org_z = static_cast<float>(c.z());
229 rayhit.ray.dir_x = static_cast<float>(dir.x());
230 rayhit.ray.dir_y = static_cast<float>(dir.y());
231 rayhit.ray.dir_z = static_cast<float>(dir.z());
232
233 rayhit.ray.tnear = 0.001f; // Start a bit away from the origin
234 rayhit.ray.tfar = std::numeric_limits<float>::infinity();
235 rayhit.hit.geomID = RTC_INVALID_GEOMETRY_ID;
236
237 rayhit.ray.mask = -1;
238 rayhit.ray.flags = 0;
239 rayhit.hit.instID[0] = RTC_INVALID_GEOMETRY_ID;
240
241 rtcIntersect1(scene_, &rayhit);
242
243 if (rayhit.hit.geomID == RTC_INVALID_GEOMETRY_ID) { // No hit
244 accum_ptr[0] += 1;
245 // ReUseX::core::trace(fmt::format(fmt::fg(fmt::color::red), "No
246 // hit"));
247 continue;
248 }
249
250 // Check if backside
251 const auto normal =
252 Eigen::Vector3f(rayhit.hit.Ng_x, rayhit.hit.Ng_y, rayhit.hit.Ng_z)
253 .normalized();
254 const Eigen::Vector3f dir_vec(dir.x(), dir.y(), dir.z());
255 if (normal.dot(dir_vec) > 0) {
256 accum_ptr[0] += 1;
257 // ReUseX::core::trace(
258 // fmt::format(fmt::fg(fmt::color::yellow), "Backside hit"));
259 continue; // Backside
260 }
261 // ReUseX::core::trace(fmt::format(fmt::fg(fmt::color::green),
262 // "Frontside hit"));
263
264 auto geometry = rtcGetGeometry(scene_, rayhit.hit.geomID);
265 rtcGetGeometryUserData(geometry);
266 RTCData *data = (RTCData *)rtcGetGeometryUserData(geometry);
267
268 // const auto label = labels[data->label_index];
269 // ReUseX::core::trace("Cell {:>3} ray {} hit label {} (index {})",
270 // (*this)[*cit].id, i, label, data->label_index);
271
272 // #pragma omp critical
273 accum_ptr[data->label_index + 1] += 1;
274 // ReUseX::core::trace("Cell {:>3} hit label {} (index {})", idx,
275 // labels[data->label_index], data->label_index);
276 }
277
278 // Compute probabilities
279 // double sum = std::accumulate(c_rp[*cit].begin(), c_rp[*cit].end(),
280 // 0.0); ReUseX::core::trace("Cell {:>3} sum = {:.3f}", idx, sum); if (sum
281 // > 0)
282 // for (size_t j = 0; j < c_rp[*cit].size(); ++j)
283 // c_rp[*cit][j] /= sum;
284 //
285 for (auto &p : local_accum)
286 p /= dirs.size();
287 c_rp[*cit] = std::move(local_accum);
288
289 ++observer;
290 }
291 }
292 rtcReleaseScene(scene_);
293 rtcReleaseDevice(device_);
294
295 _mm_setcsr(old_mxcsr); // restore old flags
296
297 // Log sum of probabilities
298 auto sum_results = std::vector<double>(labels.size() + 1, 0.0);
299 for (auto cit = this->cells_begin(); cit != this->cells_end(); ++cit)
300 for (size_t i = 0; i < c_rp[*cit].size(); ++i)
301 sum_results[i] += c_rp[*cit][i];
302 ReUseX::core::trace("Sum probabilities => [{:.3f}]",
303 fmt::join(sum_results, ", "));
304}
305} // namespace ReUseX::geometry
static std::vector< Eigen::Vector3d > sampleSphericalFibonacci(size_t N)
auto compute_room_probabilities(pcl::PointCloud< PointT >::ConstPtr cloud, pcl::PointCloud< PointN >::ConstPtr normals, pcl::PointCloud< PointL >::ConstPtr labels, const double grid_size=0.2) -> void
auto cells_end() const -> CellIterator
auto cells_begin() const -> CellIterator
std::pair< boost::associative_property_map< std::map< Key, T > >, bool > add_property_map(const std::string &name)
Definition Registry.hpp:27