ACloudViewer  3.9.4
A Modern Library for 3D Data Processing
raycasting_scene.cpp
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 // - CloudViewer: www.cloudViewer.org -
3 // ----------------------------------------------------------------------------
4 // Copyright (c) 2018-2024 www.cloudViewer.org
5 // SPDX-License-Identifier: MIT
6 // ----------------------------------------------------------------------------
7 
11 
12 namespace cloudViewer {
13 namespace t {
14 namespace geometry {
15 
16 void pybind_raycasting_scene(py::module& m) {
17  py::class_<RaycastingScene> raycasting_scene(m, "RaycastingScene", R"doc(
18 A scene class with basic ray casting and closest point queries.
19 
20 The RaycastingScene allows to compute ray intersections with triangle meshes
21 or compute the closest point on the surface of a mesh with respect to one
22 or more query points.
23 It builds an internal acceleration structure to speed up those queries.
24 
25 This class supports only the CPU device.
26 
27 The following shows how to create a scene and compute ray intersections::
28 
29  import cloudViewer as cv3d
30  import matplotlib.pyplot as plt
31 
32  cube = cv3d.t.geometry.TriangleMesh.from_legacy(
33  cv3d.geometry.ccMesh.create_box())
34 
35  # Create scene and add the cube mesh
36  scene = cv3d.t.geometry.RaycastingScene()
37  scene.add_triangles(cube)
38 
39  # Rays are 6D vectors with origin and ray direction.
40  # Here we use a helper function to create rays for a pinhole camera.
41  rays = scene.create_rays_pinhole(fov_deg=60,
42  center=[0.5,0.5,0.5],
43  eye=[-1,-1,-1],
44  up=[0,0,1],
45  width_px=320,
46  height_px=240)
47 
48  # Compute the ray intersections.
49  ans = scene.cast_rays(rays)
50 
51  # Visualize the hit distance (depth)
52  plt.imshow(ans['t_hit'].numpy())
53 
54 )doc");
55 
56  // Constructors.
57  raycasting_scene.def(py::init<int64_t, core::Device>(), "nthreads"_a = 0,
58  "device"_a = core::Device("CPU:0"), R"doc(
59 Create a RaycastingScene.
60 
61 Args:
62  nthreads (int): The number of threads to use for building the scene. Set to 0 for automatic.
63  device (cloudViewer.core.Device): The device to use. Currently CPU and SYCL devices are supported.
64 )doc");
65 
66  raycasting_scene.def(
67  "add_triangles",
68  py::overload_cast<const core::Tensor&, const core::Tensor&>(
70  "vertex_positions"_a, "triangle_indices"_a, R"doc(
71 Add a triangle mesh to the scene.
72 
73 Args:
74  vertex_positions (cloudViewer.core.Tensor): Vertices as Tensor of dim {N,3} and dtype
75  Float32.
76  triangle_indices (cloudViewer.core.Tensor): Triangles as Tensor of dim {M,3} and dtype
77  UInt32.
78 
79 Returns:
80  The geometry ID of the added mesh.
81 )doc");
82 
83  raycasting_scene.def("add_triangles",
84  py::overload_cast<const TriangleMesh&>(
86  "mesh"_a, R"doc(
87 Add a triangle mesh to the scene.
88 
89 Args:
90  mesh (cloudViewer.t.geometry.TriangleMesh): A triangle mesh.
91 
92 Returns:
93  The geometry ID of the added mesh.
94 )doc");
95 
96  raycasting_scene.def("cast_rays", &RaycastingScene::CastRays, "rays"_a,
97  "nthreads"_a = 0,
98  R"doc(
99 Computes the first intersection of the rays with the scene.
100 
101 Args:
102  rays (cloudViewer.core.Tensor): A tensor with >=2 dims, shape {.., 6}, and Dtype
103  Float32 describing the rays.
104  {..} can be any number of dimensions, e.g., to organize rays for
105  creating an image the shape can be {height, width, 6}. The last
106  dimension must be 6 and has the format [ox, oy, oz, dx, dy, dz]
107  with [ox,oy,oz] as the origin and [dx,dy,dz] as the direction. It is
108  not necessary to normalize the direction but the returned hit distance
109  uses the length of the direction vector as unit.
110 
111  nthreads (int): The number of threads to use. Set to 0 for automatic.
112 
113 Returns:
114  A dictionary which contains the following keys
115 
116  t_hit
117  A tensor with the distance to the first hit. The shape is {..}. If there
118  is no intersection the hit distance is *inf*.
119 
120  geometry_ids
121  A tensor with the geometry IDs. The shape is {..}. If there
122  is no intersection the ID is *INVALID_ID*.
123 
124  primitive_ids
125  A tensor with the primitive IDs, which corresponds to the triangle
126  index. The shape is {..}. If there is no intersection the ID is
127  *INVALID_ID*.
128 
129  primitive_uvs
130  A tensor with the barycentric coordinates of the hit points within the
131  hit triangles. The shape is {.., 2}.
132 
133  primitive_normals
134  A tensor with the normals of the hit triangles. The shape is {.., 3}.
135 )doc");
136 
137  raycasting_scene.def("test_occlusions", &RaycastingScene::TestOcclusions,
138  "rays"_a, "tnear"_a = 0.f,
140  "nthreads"_a = 0,
141  R"doc(
142 Checks if the rays have any intersection with the scene.
143 
144 Args:
145  rays (cloudViewer.core.Tensor): A tensor with >=2 dims, shape {.., 6}, and Dtype
146  Float32 describing the rays.
147  {..} can be any number of dimensions, e.g., to organize rays for
148  creating an image the shape can be {height, width, 6}.
149  The last dimension must be 6 and has the format [ox, oy, oz, dx, dy, dz]
150  with [ox,oy,oz] as the origin and [dx,dy,dz] as the direction. It is not
151  necessary to normalize the direction.
152 
153  tnear (float): The tnear offset for the rays. The default is 0.
154 
155  tfar (float): The tfar value for the ray. The default is infinity.
156 
157  nthreads (int): The number of threads to use. Set to 0 for automatic.
158 
159 Returns:
160  A boolean tensor which indicates if the ray is occluded by the scene (true)
161  or not (false).
162 )doc");
163 
164  raycasting_scene.def("count_intersections",
166  "nthreads"_a = 0, R"doc(
167 Computes the number of intersection of the rays with the scene.
168 
169 Args:
170  rays (cloudViewer.core.Tensor): A tensor with >=2 dims, shape {.., 6}, and Dtype
171  Float32 describing the rays.
172  {..} can be any number of dimensions, e.g., to organize rays for
173  creating an image the shape can be {height, width, 6}.
174  The last dimension must be 6 and has the format [ox, oy, oz, dx, dy, dz]
175  with [ox,oy,oz] as the origin and [dx,dy,dz] as the direction. It is not
176  necessary to normalize the direction.
177 
178  nthreads (int): The number of threads to use. Set to 0 for automatic.
179 
180 Returns:
181  A tensor with the number of intersections. The shape is {..}.
182 )doc");
183 
184  raycasting_scene.def("list_intersections",
186  "nthreads"_a = 0, R"doc(
187 Lists the intersections of the rays with the scene::
188 
189  import cloudViewer as cv3d
190  import numpy as np
191 
192  # Create scene and add the monkey model.
193  scene = cv3d.t.geometry.RaycastingScene()
194  d = cv3d.data.MonkeyModel()
195  mesh = cv3d.t.io.read_triangle_mesh(d.path)
196  mesh_id = scene.add_triangles(mesh)
197 
198  # Create a grid of rays covering the bounding box
199  bb_min = mesh.vertex['positions'].min(dim=0).numpy()
200  bb_max = mesh.vertex['positions'].max(dim=0).numpy()
201  x,y = np.linspace(bb_min, bb_max, num=10)[:,:2].T
202  xv, yv = np.meshgrid(x,y)
203  orig = np.stack([xv, yv, np.full_like(xv, bb_min[2]-1)], axis=-1).reshape(-1,3)
204  dest = orig + np.full(orig.shape, (0,0,2+bb_max[2]-bb_min[2]),dtype=np.float32)
205  rays = np.concatenate([orig, dest-orig], axis=-1).astype(np.float32)
206 
207  # Compute the ray intersections.
208  lx = scene.list_intersections(rays)
209  lx = {k:v.numpy() for k,v in lx.items()}
210 
211  # Calculate intersection coordinates using the primitive uvs and the mesh
212  v = mesh.vertex['positions'].numpy()
213  t = mesh.triangle['indices'].numpy()
214  tidx = lx['primitive_ids']
215  uv = lx['primitive_uvs']
216  w = 1 - np.sum(uv, axis=1)
217  c = \
218  v[t[tidx, 1].flatten(), :] * uv[:, 0][:, None] + \
219  v[t[tidx, 2].flatten(), :] * uv[:, 1][:, None] + \
220  v[t[tidx, 0].flatten(), :] * w[:, None]
221 
222  # Calculate intersection coordinates using ray_ids
223  c = rays[lx['ray_ids']][:,:3] + rays[lx['ray_ids']][:,3:]*lx['t_hit'][...,None]
224 
225  # Visualize the rays and intersections.
226  lines = cv3d.t.geometry.LineSet()
227  lines.point.positions = np.hstack([orig,dest]).reshape(-1,3)
228  lines.line.indices = np.arange(lines.point.positions.shape[0]).reshape(-1,2)
229  lines.line.colors = np.full((lines.line.indices.shape[0],3), (1,0,0))
230  x = cv3d.t.geometry.PointCloud(positions=c)
231  cv3d.visualization.draw([mesh, lines, x], point_size=8)
232 
233 
234 Args:
235  rays (cloudViewer.core.Tensor): A tensor with >=2 dims, shape {.., 6}, and Dtype
236  Float32 describing the rays; {..} can be any number of dimensions.
237  The last dimension must be 6 and has the format [ox, oy, oz, dx, dy, dz]
238  with [ox,oy,oz] as the origin and [dx,dy,dz] as the direction. It is not
239  necessary to normalize the direction although it should be normalised if
240  t_hit is to be calculated in coordinate units.
241 
242  nthreads (int): The number of threads to use. Set to 0 for automatic.
243 
244 Returns:
245  The returned dictionary contains
246 
247  ray_splits
248  A tensor with ray intersection splits. Can be used to iterate over all intersections for each ray. The shape is {num_rays + 1}.
249 
250  ray_ids
251  A tensor with ray IDs. The shape is {num_intersections}.
252 
253  t_hit
254  A tensor with the distance to the hit. The shape is {num_intersections}.
255 
256  geometry_ids
257  A tensor with the geometry IDs. The shape is {num_intersections}.
258 
259  primitive_ids
260  A tensor with the primitive IDs, which corresponds to the triangle
261  index. The shape is {num_intersections}.
262 
263  primitive_uvs
264  A tensor with the barycentric coordinates of the intersection points within
265  the triangles. The shape is {num_intersections, 2}.
266 
267 
268 An example of using ray_splits::
269 
270  ray_splits: [0, 2, 3, 6, 6, 8] # note that the length of this is num_rays+1
271  t_hit: [t1, t2, t3, t4, t5, t6, t7, t8]
272 
273  for ray_id, (start, end) in enumerate(zip(ray_splits[:-1], ray_splits[1:])):
274  for i,t in enumerate(t_hit[start:end]):
275  print(f'ray {ray_id}, intersection {i} at {t}')
276 
277 
278 )doc");
279 
280  raycasting_scene.def("compute_closest_points",
282  "query_points"_a, "nthreads"_a = 0, R"doc(
283 Computes the closest points on the surfaces of the scene.
284 
285 Args:
286  query_points (cloudViewer.core.Tensor): A tensor with >=2 dims, shape {.., 3},
287  and Dtype Float32 describing the query points.
288  {..} can be any number of dimensions, e.g., to organize the query_point
289  to create a 3D grid the shape can be {depth, height, width, 3}.
290  The last dimension must be 3 and has the format [x, y, z].
291 
292  nthreads (int): The number of threads to use. Set to 0 for automatic.
293 
294 Returns:
295  The returned dictionary contains
296 
297  points
298  A tensor with the closest surface points. The shape is {..}.
299 
300  geometry_ids
301  A tensor with the geometry IDs. The shape is {..}.
302 
303  primitive_ids
304  A tensor with the primitive IDs, which corresponds to the triangle
305  index. The shape is {..}.
306 
307  primitive_uvs
308  A tensor with the barycentric coordinates of the closest points within
309  the triangles. The shape is {.., 2}.
310 
311  primitive_normals
312  A tensor with the normals of the closest triangle . The shape is
313  {.., 3}.
314 
315 )doc");
316 
317  raycasting_scene.def("compute_distance", &RaycastingScene::ComputeDistance,
318  "query_points"_a, "nthreads"_a = 0, R"doc(
319 Computes the distance to the surface of the scene.
320 
321 Args:
322  query_points (cloudViewer.core.Tensor): A tensor with >=2 dims, shape {.., 3},
323  and Dtype Float32 describing the query points.
324  {..} can be any number of dimensions, e.g., to organize the
325  query points to create a 3D grid the shape can be
326  {depth, height, width, 3}.
327  The last dimension must be 3 and has the format [x, y, z].
328 
329  nthreads (int): The number of threads to use. Set to 0 for automatic.
330 
331 Returns:
332  A tensor with the distances to the surface. The shape is {..}.
333 )doc");
334 
335  raycasting_scene.def(
336  "compute_signed_distance", &RaycastingScene::ComputeSignedDistance,
337  "query_points"_a, "nthreads"_a = 0, "nsamples"_a = 1, R"doc(
338 Computes the signed distance to the surface of the scene.
339 
340 This function computes the signed distance to the meshes in the scene.
341 The function assumes that all meshes are watertight and that there are
342 no intersections between meshes, i.e., inside and outside must be well
343 defined. The function determines the sign of the distance by counting
344 the intersections of a rays starting at the query points.
345 
346 Args:
347  query_points (cloudViewer.core.Tensor): A tensor with >=2 dims, shape {.., 3},
348  and Dtype Float32 describing the query_points.
349  {..} can be any number of dimensions, e.g., to organize the
350  query points to create a 3D grid the shape can be
351  {depth, height, width, 3}.
352  The last dimension must be 3 and has the format [x, y, z].
353 
354  nthreads (int): The number of threads to use. Set to 0 for automatic.
355 
356  nsamples (int): The number of rays used for determining the inside.
357  This must be an odd number. The default is 1. Use a higher value if you
358  notice sign flipping, which can occur when rays hit exactly an edge or
359  vertex in the scene.
360 
361 Returns:
362  A tensor with the signed distances to the surface. The shape is {..}.
363  Negative distances mean a point is inside a closed surface.
364 )doc");
365 
366  raycasting_scene.def("compute_occupancy",
367  &RaycastingScene::ComputeOccupancy, "query_points"_a,
368  "nthreads"_a = 0, "nsamples"_a = 1,
369  R"doc(
370 Computes the occupancy at the query point positions.
371 
372 This function computes whether the query points are inside or outside.
373 The function assumes that all meshes are watertight and that there are
374 no intersections between meshes, i.e., inside and outside must be well
375 defined. The function determines if a point is inside by counting the
376 intersections of a rays starting at the query points.
377 
378 Args:
379  query_points (cloudViewer.core.Tensor): A tensor with >=2 dims, shape {.., 3},
380  and Dtype Float32 describing the query points.
381  {..} can be any number of dimensions, e.g., to organize the
382  query points to create a 3D grid the shape can be
383  {depth, height, width, 3}.
384  The last dimension must be 3 and has the format [x, y, z].
385 
386  nthreads (int): The number of threads to use. Set to 0 for automatic.
387 
388  nsamples (int): The number of rays used for determining the inside.
389  This must be an odd number. The default is 1. Use a higher value if you
390  notice errors in the occupancy values. Errors can occur when rays hit
391  exactly an edge or vertex in the scene.
392 
393 Returns:
394  A tensor with the occupancy values. The shape is {..}. Values are either 0
395  or 1. A point is occupied or inside if the value is 1.
396 )doc");
397 
398  raycasting_scene.def_static(
399  "create_rays_pinhole",
400  py::overload_cast<const core::Tensor&, const core::Tensor&, int,
402  "intrinsic_matrix"_a, "extrinsic_matrix"_a, "width_px"_a,
403  "height_px"_a, R"doc(
404 Creates rays for the given camera parameters.
405 
406 Args:
407  intrinsic_matrix (cloudViewer.core.Tensor): The upper triangular intrinsic matrix
408  with shape {3,3}.
409  extrinsic_matrix (cloudViewer.core.Tensor): The 4x4 world to camera SE(3)
410  transformation matrix.
411  width_px (int): The width of the image in pixels.
412  height_px (int): The height of the image in pixels.
413 
414 Returns:
415  A tensor of shape {height_px, width_px, 6} with the rays.
416 )doc");
417 
418  raycasting_scene.def_static(
419  "create_rays_pinhole",
420  py::overload_cast<double, const core::Tensor&, const core::Tensor&,
421  const core::Tensor&, int, int>(
423  "fov_deg"_a, "center"_a, "eye"_a, "up"_a, "width_px"_a,
424  "height_px"_a, R"doc(
425 Creates rays for the given camera parameters.
426 
427 Args:
428  fov_deg (float): The horizontal field of view in degree.
429  center (cloudViewer.core.Tensor): The point the camera is looking at with shape
430  {3}.
431  eye (cloudViewer.core.Tensor): The position of the camera with shape {3}.
432  up (cloudViewer.core.Tensor): The up-vector with shape {3}.
433  width_px (int): The width of the image in pixels.
434  height_px (int): The height of the image in pixels.
435 
436 Returns:
437  A tensor of shape {height_px, width_px, 6} with the rays.
438 )doc");
439 
440  raycasting_scene.def_property_readonly_static(
441  "INVALID_ID",
442  [](py::object /* self */) -> uint32_t {
444  },
445  R"doc(
446 The value for invalid IDs
447 )doc");
448 }
449 } // namespace geometry
450 } // namespace t
451 } // namespace cloudViewer
static core::Tensor CreateRaysPinhole(const core::Tensor &intrinsic_matrix, const core::Tensor &extrinsic_matrix, int width_px, int height_px)
Creates rays for the given camera parameters.
std::unordered_map< std::string, core::Tensor > ListIntersections(const core::Tensor &rays, const int nthreads=0)
Lists the intersections of the rays with the scene.
uint32_t AddTriangles(const core::Tensor &vertex_positions, const core::Tensor &triangle_indices)
Add a triangle mesh to the scene.
std::unordered_map< std::string, core::Tensor > CastRays(const core::Tensor &rays, const int nthreads=0) const
Computes the first intersection of the rays with the scene.
core::Tensor CountIntersections(const core::Tensor &rays, const int nthreads=0)
Computes the number of intersection of the rays with the scene.
core::Tensor ComputeOccupancy(const core::Tensor &query_points, const int nthreads=0, const int nsamples=1)
Computes the occupancy at the query point positions.
core::Tensor ComputeDistance(const core::Tensor &query_points, const int nthreads=0)
Computes the distance to the surface of the scene.
static uint32_t INVALID_ID()
The value for invalid IDs.
core::Tensor ComputeSignedDistance(const core::Tensor &query_points, const int nthreads=0, const int nsamples=1)
Computes the signed distance to the surface of the scene.
core::Tensor TestOcclusions(const core::Tensor &rays, const float tnear=0.f, const float tfar=std::numeric_limits< float >::infinity(), const int nthreads=0)
Checks if the rays have any intersection with the scene.
std::unordered_map< std::string, core::Tensor > ComputeClosestPoints(const core::Tensor &query_points, const int nthreads=0)
Computes the closest points on the surfaces of the scene.
__device__ __forceinline__ float infinity()
Definition: result_set.h:36
void pybind_raycasting_scene(py::module &m)
Generic file read and write utility for python interface.