0.1.22
Loading...
Searching...
No Matches
LiDARCloud.py
Go to the documentation of this file.
1"""
2LiDARCloud - High-level interface for LiDAR simulation and point cloud processing
3
4Provides Python interface to Helios LiDAR plugin for:
5- Synthetic LiDAR scanning
6- Point cloud management and filtering
7- Triangulation and mesh generation
8- Leaf area density calculations
9"""
10
11from enum import IntEnum
12from typing import List, Tuple, Optional, Union
13from .wrappers import ULiDARWrapper as lidar_wrapper
14from .Context import Context
15from .plugins.registry import get_plugin_registry
16from .exceptions import HeliosError
17from .wrappers.DataTypes import vec3, RGBcolor, SphericalCoord
18from .validation.datatypes import validate_vec3
19from .validation.core import validate_positive_value
20
21
23 """Exception raised for LiDAR-specific errors"""
24 pass
25
26
27class ScanPattern(IntEnum):
28 """Scan pattern returned by :meth:`LiDARCloud.getScanPattern`.
29
30 RASTER is the uniform-angular-grid pattern produced by :meth:`LiDARCloud.addScan`;
31 SPINNING_MULTIBEAM is the rotating multi-channel pattern produced by
32 :meth:`LiDARCloud.addScanMultibeam`.
33 """
34 RASTER = 0
35 SPINNING_MULTIBEAM = 1
36
37
38class LiDARCloud:
39 """
40 High-level interface for LiDAR point cloud operations.
41
42 Supports synthetic scanning, point cloud filtering, triangulation,
43 and leaf area density calculations.
44
45 Example:
46 >>> from pyhelios import LiDARCloud
47 >>> from pyhelios.types import vec3
48 >>>
49 >>> with LiDARCloud() as lidar:
50 ... # Add a scan
51 ... scan_id = lidar.addScan(
52 ... origin=vec3(0, 0, 1),
53 ... Ntheta=100, theta_range=(0, 1.57),
54 ... Nphi=100, phi_range=(-3.14, 3.14),
55 ... exit_diameter=0.01, beam_divergence=0.001
56 ... )
57 ...
58 ... # Add hit points
59 ... lidar.addHitPoint(scan_id, vec3(1, 0, 0), vec3(1, 0, 0))
60 ...
61 ... # Export point cloud
62 ... lidar.exportPointCloud("output.xyz")
63 """
64
65 def __init__(self):
66 """
67 Initialize LiDARCloud.
68
69 Raises:
70 LiDARError: If plugin not available in current build
71 RuntimeError: If cloud initialization fails
72 """
73 # Check plugin availability
74 registry = get_plugin_registry()
75 if not registry.is_plugin_available('lidar'):
76 raise LiDARError(
77 "LiDAR plugin not available. Rebuild PyHelios with LiDAR:\n"
78 " build_scripts/build_helios --plugins lidar\n"
79 "\n"
80 "System requirements:\n"
81 " - Platforms: Windows, Linux, macOS\n"
82 " - GPU: Optional (enables GPU acceleration)"
83 )
84
85 self._cloud_ptr = lidar_wrapper.createLiDARcloud()
86 if not self._cloud_ptr:
87 raise LiDARError("Failed to create LiDAR cloud")
88
89 def __enter__(self):
90 """Context manager entry"""
91 return self
92
93 def __exit__(self, exc_type, exc_val, exc_tb):
94 """Context manager exit - cleanup resources"""
95 if hasattr(self, '_cloud_ptr') and self._cloud_ptr:
96 lidar_wrapper.destroyLiDARcloud(self._cloud_ptr)
97 self._cloud_ptr = None
98
99 def __del__(self):
100 """Fallback destructor for cleanup without context manager"""
101 if hasattr(self, '_cloud_ptr') and self._cloud_ptr is not None:
102 try:
103 lidar_wrapper.destroyLiDARcloud(self._cloud_ptr)
104 self._cloud_ptr = None
105 except Exception as e:
106 import warnings
107 warnings.warn(f"Error in LiDARCloud.__del__: {e}")
108
109 def addScan(self, origin: Union[vec3, List[float], Tuple[float, float, float]],
110 Ntheta: int, theta_range: Tuple[float, float],
111 Nphi: int, phi_range: Tuple[float, float],
112 exit_diameter: float, beam_divergence: float,
113 column_format: Optional[List[str]] = None,
114 range_noise_stddev: float = 0.0, angle_noise_stddev: float = 0.0,
115 scan_tilt_roll: float = 0.0, scan_tilt_pitch: float = 0.0) -> int:
116 """
117 Add a LiDAR scan to the point cloud.
118
119 Args:
120 origin: Scanner position (vec3 or 3-element list/tuple)
121 Ntheta: Number of scan points in zenith direction
122 theta_range: Zenith angle range (min, max) in radians
123 Nphi: Number of scan points in azimuthal direction
124 phi_range: Azimuthal angle range (min, max) in radians
125 exit_diameter: Laser beam exit diameter (meters)
126 beam_divergence: Beam divergence angle (radians)
127 column_format: Optional list of column-format labels. Non-standard labels
128 (anything other than geometry/standard tokens like x/y/z/r/g/b/raydir)
129 cause syntheticScan to sample that named primitive data from the struck
130 primitive onto each hit's data map, retrievable via getHitData(). Defaults
131 to None (empty format).
132
133 One label is special: "reflectivity_lidar" modulates each hit's "intensity"
134 (intensity *= reflectivity) rather than being stored as its own hit-data
135 key, so getHitData(i, "reflectivity_lidar") will NOT return it.
136 range_noise_stddev: Standard deviation of Gaussian range (along-beam) measurement
137 noise in meters. Only affects synthetic-scan generation. Defaults to 0.0
138 (noise disabled).
139 angle_noise_stddev: Standard deviation of Gaussian angular (beam-pointing) jitter
140 in radians. Only affects synthetic-scan generation. Defaults to 0.0 (jitter
141 disabled).
142 scan_tilt_roll: Global scanner tilt roll angle in radians, modeling residual tilt of
143 the scanner spin axis away from plumb (right-hand rotation about the body lateral
144 axis). Only affects synthetic-scan generation. Defaults to 0.0 (level).
145 scan_tilt_pitch: Global scanner tilt pitch angle in radians (right-hand rotation about
146 the body forward/azimuth-zero axis). Only affects synthetic-scan generation.
147 Defaults to 0.0 (level).
148
149 Returns:
150 Scan ID for referencing this scan
151
152 Example:
153 >>> scan_id = lidar.addScan(
154 ... origin=vec3(0, 0, 1),
155 ... Ntheta=100, theta_range=(0, 1.57),
156 ... Nphi=100, phi_range=(-3.14, 3.14),
157 ... exit_diameter=0.01, beam_divergence=0.001,
158 ... column_format=["my_scalar"]
159 ... )
160 """
161 # Convert origin to vec3 if needed
162 if isinstance(origin, (list, tuple)):
163 if len(origin) != 3:
164 raise ValueError("Origin must have 3 elements [x, y, z]")
165 origin = vec3(*origin)
166 elif not hasattr(origin, 'x'):
167 raise ValueError("Origin must be vec3 or 3-element list/tuple")
168
169 origin_list = [origin.x, origin.y, origin.z]
170
171 # Validate scan parameters
172 validate_positive_value(Ntheta, 'Ntheta', 'addScan')
173 validate_positive_value(Nphi, 'Nphi', 'addScan')
174
175 if not isinstance(theta_range, (list, tuple)) or len(theta_range) != 2:
176 raise ValueError("theta_range must be a tuple (min, max)")
177 if not isinstance(phi_range, (list, tuple)) or len(phi_range) != 2:
178 raise ValueError("phi_range must be a tuple (min, max)")
179
180 if column_format is not None:
181 if not isinstance(column_format, (list, tuple)) or \
182 not all(isinstance(c, str) for c in column_format):
183 raise ValueError("column_format must be a list of strings")
184 column_format = list(column_format)
185
186 if range_noise_stddev < 0:
187 raise ValueError("range_noise_stddev must be non-negative")
188 if angle_noise_stddev < 0:
189 raise ValueError("angle_noise_stddev must be non-negative")
190
191 return lidar_wrapper.addLiDARScan(
192 self._cloud_ptr, origin_list, Ntheta, theta_range,
193 Nphi, phi_range, exit_diameter, beam_divergence, column_format,
194 range_noise_stddev, angle_noise_stddev,
195 scan_tilt_roll, scan_tilt_pitch
196 )
197
198 def addScanMultibeam(self, origin: Union[vec3, List[float], Tuple[float, float, float]],
199 beam_zenith_angles: List[float],
200 Nphi: int, phi_range: Tuple[float, float],
201 exit_diameter: float, beam_divergence: float,
202 column_format: Optional[List[str]] = None,
203 range_noise_stddev: float = 0.0, angle_noise_stddev: float = 0.0,
204 scan_tilt_roll: float = 0.0, scan_tilt_pitch: float = 0.0) -> int:
205 """
206 Add a spinning multibeam LiDAR scan (rotating multi-channel sensor, e.g. Velodyne/Ouster/Hesai).
207
208 Each laser channel is fired at a fixed zenith angle (taken from ``beam_zenith_angles``) as the
209 sensor head rotates through ``Nphi`` uniform azimuth steps. The scan is stored as an
210 (len(beam_zenith_angles) x Nphi) table, so all downstream processing is shared with raster scans.
211
212 Args:
213 origin: Scanner position (vec3 or 3-element list/tuple)
214 beam_zenith_angles: Per-channel zenith angles in radians (0 = upward, pi/2 = horizontal,
215 pi = downward). Manufacturer spec sheets typically list channel angles as elevation
216 above the horizon; zenith = pi/2 - elevation. Its length sets Ntheta (number of channels).
217 Nphi: Number of azimuth steps (columns) per rotation
218 phi_range: Azimuthal angle range (min, max) in radians
219 exit_diameter: Laser beam exit diameter (meters)
220 beam_divergence: Beam divergence angle (radians)
221 column_format: Optional list of column-format labels (see addScan)
222 range_noise_stddev: Std. dev. of Gaussian range noise in meters (default 0)
223 angle_noise_stddev: Std. dev. of Gaussian angular jitter in radians (default 0)
224 scan_tilt_roll: Global scanner tilt roll angle in radians (default 0, level)
225 scan_tilt_pitch: Global scanner tilt pitch angle in radians (default 0, level)
226
227 Returns:
228 Scan ID for referencing this scan
229 """
230 if isinstance(origin, (list, tuple)):
231 if len(origin) != 3:
232 raise ValueError("Origin must have 3 elements [x, y, z]")
233 origin = vec3(*origin)
234 elif not hasattr(origin, 'x'):
235 raise ValueError("Origin must be vec3 or 3-element list/tuple")
236
237 origin_list = [origin.x, origin.y, origin.z]
238
239 if not isinstance(beam_zenith_angles, (list, tuple)) or len(beam_zenith_angles) == 0:
240 raise ValueError("beam_zenith_angles must be a non-empty list of per-channel angles")
241 if not all(isinstance(a, (int, float)) for a in beam_zenith_angles):
242 raise ValueError("beam_zenith_angles must be a list of numbers (radians)")
243
244 validate_positive_value(Nphi, 'Nphi', 'addScanMultibeam')
245
246 if not isinstance(phi_range, (list, tuple)) or len(phi_range) != 2:
247 raise ValueError("phi_range must be a tuple (min, max)")
248
249 if column_format is not None:
250 if not isinstance(column_format, (list, tuple)) or \
251 not all(isinstance(c, str) for c in column_format):
252 raise ValueError("column_format must be a list of strings")
253 column_format = list(column_format)
254
255 if range_noise_stddev < 0:
256 raise ValueError("range_noise_stddev must be non-negative")
257 if angle_noise_stddev < 0:
258 raise ValueError("angle_noise_stddev must be non-negative")
259
260 return lidar_wrapper.addLiDARScanMultibeam(
261 self._cloud_ptr, origin_list, list(beam_zenith_angles),
262 Nphi, phi_range, exit_diameter, beam_divergence, column_format,
263 range_noise_stddev, angle_noise_stddev,
264 scan_tilt_roll, scan_tilt_pitch
265 )
266
267 def getScanCount(self) -> int:
268 """Get total number of scans in the cloud"""
269 return lidar_wrapper.getLiDARScanCount(self._cloud_ptr)
270
271 def getScanOrigin(self, scanID: int) -> vec3:
272 """Get origin of a specific scan"""
273 if scanID < 0:
274 raise ValueError("Scan ID must be non-negative")
275 origin_list = lidar_wrapper.getLiDARScanOrigin(self._cloud_ptr, scanID)
276 return vec3(*origin_list)
277
278 def getScanSizeTheta(self, scanID: int) -> int:
279 """Get number of zenith scan points for a scan"""
280 if scanID < 0:
281 raise ValueError("Scan ID must be non-negative")
282 return lidar_wrapper.getLiDARScanSizeTheta(self._cloud_ptr, scanID)
283
284 def getScanSizePhi(self, scanID: int) -> int:
285 """Get number of azimuthal scan points for a scan"""
286 if scanID < 0:
287 raise ValueError("Scan ID must be non-negative")
288 return lidar_wrapper.getLiDARScanSizePhi(self._cloud_ptr, scanID)
289
290 def getScanRangeNoiseStdDev(self, scanID: int) -> float:
291 """Get the range (along-beam) measurement noise standard deviation for a scan (meters).
292
293 Returns the value supplied to addScan() as ``range_noise_stddev`` (0.0 if disabled).
294 """
295 if scanID < 0:
296 raise ValueError("Scan ID must be non-negative")
297 return lidar_wrapper.getLiDARScanRangeNoiseStdDev(self._cloud_ptr, scanID)
298
299 def getScanAngleNoiseStdDev(self, scanID: int) -> float:
300 """Get the angular (beam-pointing) jitter standard deviation for a scan (radians).
301
302 Returns the value supplied to addScan() as ``angle_noise_stddev`` (0.0 if disabled).
303 """
304 if scanID < 0:
305 raise ValueError("Scan ID must be non-negative")
306 return lidar_wrapper.getLiDARScanAngleNoiseStdDev(self._cloud_ptr, scanID)
307
308 def getScanTiltRoll(self, scanID: int) -> float:
309 """Get the global scanner tilt roll angle for a scan (radians; 0.0 if level)."""
310 if scanID < 0:
311 raise ValueError("Scan ID must be non-negative")
312 return lidar_wrapper.getLiDARScanTiltRoll(self._cloud_ptr, scanID)
313
314 def getScanTiltPitch(self, scanID: int) -> float:
315 """Get the global scanner tilt pitch angle for a scan (radians; 0.0 if level)."""
316 if scanID < 0:
317 raise ValueError("Scan ID must be non-negative")
318 return lidar_wrapper.getLiDARScanTiltPitch(self._cloud_ptr, scanID)
319
320 def getScanPattern(self, scanID: int) -> int:
321 """Get the scan pattern for a scan.
322
323 Returns an integer: 0 = raster (uniform angular grid), 1 = spinning multibeam
324 (rotating multi-channel sensor). Compare against ``ScanPattern.RASTER`` /
325 ``ScanPattern.SPINNING_MULTIBEAM``.
326 """
327 if scanID < 0:
328 raise ValueError("Scan ID must be non-negative")
329 return lidar_wrapper.getLiDARScanPattern(self._cloud_ptr, scanID)
330
331 def getScanBeamZenithAngles(self, scanID: int) -> List[float]:
332 """Get the per-channel beam zenith angles (radians) for a multibeam scan.
333
334 Returns an empty list for a raster scan.
335 """
336 if scanID < 0:
337 raise ValueError("Scan ID must be non-negative")
338 return lidar_wrapper.getLiDARScanBeamZenithAngles(self._cloud_ptr, scanID)
339
340 def addHitPoint(self, scanID: int,
341 xyz: Union[vec3, List[float], Tuple[float, float, float]],
342 direction: Union[vec3, SphericalCoord, List[float], Tuple[float, float]],
343 color: Optional[Union[RGBcolor, List[float], Tuple[float, float, float]]] = None):
344 """
345 Add a hit point to the point cloud.
346
347 Args:
348 scanID: Scan ID this hit belongs to
349 xyz: Hit point coordinates (vec3 or 3-element list)
350 direction: Ray direction (vec3/SphericalCoord or 2-3 element list)
351 color: Optional RGB color (RGBcolor or 3-element list)
352 """
353 # Convert xyz to list
354 if isinstance(xyz, (list, tuple)):
355 if len(xyz) != 3:
356 raise ValueError("XYZ must have 3 elements")
357 xyz_list = list(xyz)
358 elif hasattr(xyz, 'x'):
359 xyz_list = [xyz.x, xyz.y, xyz.z]
360 else:
361 raise ValueError("XYZ must be vec3 or 3-element list/tuple")
362
363 # Convert direction to list
364 if isinstance(direction, (list, tuple)):
365 if len(direction) < 2:
366 raise ValueError("Direction must have at least 2 elements [radius, elevation]")
367 direction_list = list(direction)
368 elif hasattr(direction, 'radius'): # SphericalCoord
369 direction_list = [direction.radius, direction.elevation, direction.azimuth]
370 elif hasattr(direction, 'x'): # vec3
371 direction_list = [direction.x, direction.y, direction.z]
372 else:
373 raise ValueError("Direction must be vec3/SphericalCoord or 2-3 element list")
374
375 # Add with or without color
376 if color is not None:
377 if isinstance(color, (list, tuple)):
378 if len(color) != 3:
379 raise ValueError("Color must have 3 elements [r, g, b]")
380 color_list = list(color)
381 elif hasattr(color, 'r'):
382 color_list = [color.r, color.g, color.b]
383 else:
384 raise ValueError("Color must be RGBcolor or 3-element list")
385
386 lidar_wrapper.addLiDARHitPointRGB(self._cloud_ptr, scanID, xyz_list, direction_list, color_list)
387 else:
388 lidar_wrapper.addLiDARHitPoint(self._cloud_ptr, scanID, xyz_list, direction_list)
389
390 def addHitPoints(self, scanID: int, xyz_array, direction_array, color_array=None):
391 """
392 Add many hit points to the point cloud in a single bulk call.
393
394 This skips the per-point Python loop by passing contiguous buffers
395 straight to the native library in one FFI call.
396
397 Args:
398 scanID: Scan ID these hits belong to
399 xyz_array: Hit point coordinates, shape (N, 3) [x, y, z]
400 direction_array: Ray directions, shape (N, 3) [radius, elevation, azimuth]
401 (azimuth is currently ignored, matching addHitPoint)
402 color_array: Optional RGB colors, shape (N, 3) [r, g, b]
403 """
404 import numpy as np
405
406 xyz_array = np.ascontiguousarray(xyz_array, dtype=np.float32)
407 direction_array = np.ascontiguousarray(direction_array, dtype=np.float32)
408
409 if xyz_array.ndim != 2 or xyz_array.shape[1] != 3:
410 raise ValueError("xyz_array must have shape (N, 3)")
411 if direction_array.ndim != 2 or direction_array.shape[1] != 3:
412 raise ValueError("direction_array must have shape (N, 3)")
413
414 count = xyz_array.shape[0]
415 if direction_array.shape[0] != count:
416 raise ValueError("xyz_array and direction_array must have the same number of rows")
417
418 if color_array is not None:
419 color_array = np.ascontiguousarray(color_array, dtype=np.float32)
420 if color_array.ndim != 2 or color_array.shape[1] != 3:
421 raise ValueError("color_array must have shape (N, 3)")
422 if color_array.shape[0] != count:
423 raise ValueError("color_array must have the same number of rows as xyz_array")
424
425 lidar_wrapper.addLiDARHitPoints(self._cloud_ptr, scanID,
426 xyz_array, direction_array, count, color_array)
427
428 def addHitPointsWithData(self, scanID: int, xyz_array, direction_array,
429 data_labels=None, data_values=None, color_array=None):
430 """
431 Add many hit points carrying a per-hit data map in a single bulk call.
432
433 Like addHitPoints, but also populates each hit's named-scalar data map —
434 the in-memory equivalent of what the ASCII loader does for non-standard
435 columns. This is the path multi-return LAD needs (timestamp/target_index/
436 target_count land in the map so gapfillMisses() can group beams by pulse).
437
438 Args:
439 scanID: Scan ID these hits belong to (the scan must already exist)
440 xyz_array: Hit point coordinates, shape (N, 3) [x, y, z]
441 direction_array: Ray directions, shape (N, 3) [radius, elevation, azimuth].
442 Pass cart2sphere(xyz - origin) to match loadASCIIFile;
443 the full SphericalCoord (incl. radius) is used.
444 data_labels: Optional list of data-map key names (length k)
445 data_values: Optional (N, k) values for those keys (float64)
446 color_array: Optional RGB colors, shape (N, 3) [r, g, b]
447 """
448 import numpy as np
449
450 xyz_array = np.ascontiguousarray(xyz_array, dtype=np.float32)
451 direction_array = np.ascontiguousarray(direction_array, dtype=np.float32)
452
453 if xyz_array.ndim != 2 or xyz_array.shape[1] != 3:
454 raise ValueError("xyz_array must have shape (N, 3)")
455 if direction_array.ndim != 2 or direction_array.shape[1] != 3:
456 raise ValueError("direction_array must have shape (N, 3)")
457
458 count = xyz_array.shape[0]
459 if direction_array.shape[0] != count:
460 raise ValueError("xyz_array and direction_array must have the same number of rows")
461
462 labels = list(data_labels or [])
463 if labels:
464 data_values = np.ascontiguousarray(data_values, dtype=np.float64)
465 if data_values.ndim != 2 or data_values.shape != (count, len(labels)):
466 raise ValueError("data_values must have shape (N, len(data_labels))")
467 else:
468 data_values = None
469
470 if color_array is not None:
471 color_array = np.ascontiguousarray(color_array, dtype=np.float32)
472 if color_array.ndim != 2 or color_array.shape[1] != 3:
473 raise ValueError("color_array must have shape (N, 3)")
474 if color_array.shape[0] != count:
475 raise ValueError("color_array must have the same number of rows as xyz_array")
476
477 lidar_wrapper.addLiDARHitPointsWithData(
478 self._cloud_ptr, scanID, xyz_array, direction_array, count,
479 color_array, labels, data_values)
480
481 def getHitCount(self) -> int:
482 """Get total number of hit points in cloud"""
483 return lidar_wrapper.getLiDARHitCount(self._cloud_ptr)
484
485 def getHitXYZ(self, index: int) -> vec3:
486 """Get coordinates of a hit point"""
487 if index < 0:
488 raise ValueError("Index must be non-negative")
489 xyz_list = lidar_wrapper.getLiDARHitXYZ(self._cloud_ptr, index)
490 return vec3(*xyz_list)
491
492 def getHitRaydir(self, index: int) -> SphericalCoord:
493 """Get ray direction of a hit point"""
494 if index < 0:
495 raise ValueError("Index must be non-negative")
496 direction_list = lidar_wrapper.getLiDARHitRaydir(self._cloud_ptr, index)
497 # direction_list is [radius, elevation, azimuth]; preserve azimuth (was previously dropped).
498 return SphericalCoord(direction_list[0], direction_list[1], direction_list[2])
499
500 def getHitColor(self, index: int) -> RGBcolor:
501 """Get color of a hit point"""
502 if index < 0:
503 raise ValueError("Index must be non-negative")
504 color_list = lidar_wrapper.getLiDARHitColor(self._cloud_ptr, index)
505 return RGBcolor(*color_list)
506
507 def getHitScanID(self, index: int) -> int:
508 """Get the scan ID a hit point belongs to"""
509 if index < 0:
510 raise ValueError("Index must be non-negative")
511 return lidar_wrapper.getLiDARHitScanID(self._cloud_ptr, index)
512
513 def doesHitDataExist(self, index: int, label: str) -> bool:
514 """Check whether a named scalar data value exists for a hit point.
515
516 Per-hit data computed by syntheticScan includes 'intensity', 'distance',
517 'timestamp', 'target_index', 'target_count', 'deviation', 'nRaysHit', plus any
518 primitive-data labels listed in the scan's column_format.
519 """
520 if index < 0:
521 raise ValueError("Index must be non-negative")
522 return lidar_wrapper.doesLiDARHitDataExist(self._cloud_ptr, index, label)
523
524 def getHitData(self, index: int, label: str) -> float:
525 """Get a named scalar data value for a hit point.
526
527 Raises HeliosError if the label does not exist for this hit; guard with
528 doesHitDataExist() when unsure.
529 """
530 if index < 0:
531 raise ValueError("Index must be non-negative")
532 return lidar_wrapper.getLiDARHitData(self._cloud_ptr, index, label)
533
534 def getHitDataAll(self, label: str) -> List[float]:
535 """Bulk-export a named scalar data value for all hits in a single FFI call.
536
537 Returns a list of length getHitCount(); entries are NaN where the label is
538 absent for that hit. Much faster than looping getHitData() for large clouds.
539
540 Note: values are returned at float32 precision (vs. getHitData(), which returns
541 full float64). Use getHitData() per-hit if full precision is required.
542 """
543 n = self.getHitCount()
544 if n == 0:
545 return []
546 return lidar_wrapper.getLiDARHitData_all(self._cloud_ptr, label, n)
547
548 def getHitsXYZRGB(self) -> Tuple[List[vec3], List[RGBcolor]]:
549 """Bulk-export coordinates and colors for all hits in a single FFI call.
550
551 Returns (positions, colors) where positions is a list of vec3 and colors a list
552 of RGBcolor, each of length getHitCount(). Much faster than looping
553 getHitXYZ()/getHitColor() for large clouds.
554 """
555 n = self.getHitCount()
556 if n == 0:
557 return [], []
558 xyz_flat, rgb_flat = lidar_wrapper.getLiDARHitsXYZRGB_all(self._cloud_ptr, n)
559 positions = [vec3(xyz_flat[3 * i], xyz_flat[3 * i + 1], xyz_flat[3 * i + 2]) for i in range(n)]
560 colors = [RGBcolor(rgb_flat[3 * i], rgb_flat[3 * i + 1], rgb_flat[3 * i + 2]) for i in range(n)]
561 return positions, colors
562
563 def deleteHitPoint(self, index: int):
564 """Delete a hit point from the cloud"""
565 if index < 0:
566 raise ValueError("Index must be non-negative")
567 lidar_wrapper.deleteLiDARHitPoint(self._cloud_ptr, index)
568
569 def isHitMiss(self, index: int) -> bool:
570 """Return True if a hit is a "miss" (a fired pulse that returned nothing).
571
572 Misses are the transmitted beams that form the denominator of the per-voxel
573 transmission probability used by :meth:`calculateLeafArea`. They are produced by
574 ``syntheticScan(..., record_misses=True)`` and by :meth:`gapfillMisses`.
575 """
576 if index < 0:
577 raise ValueError("Index must be non-negative")
578 return lidar_wrapper.isLiDARHitMiss(self._cloud_ptr, index)
579
580 def hasMisses(self) -> bool:
581 """Return True if the cloud contains at least one miss.
582
583 :meth:`calculateLeafArea` requires misses and fails fast without them.
584 """
585 return lidar_wrapper.lidarHasMisses(self._cloud_ptr)
586
587 @staticmethod
588 def getMissDistance() -> float:
589 """Return the LIDAR_MISS_DISTANCE constant (meters): the distance at which a
590 miss point is placed along its beam."""
591 return lidar_wrapper.getLiDARMissDistance()
592
593 def coordinateShift(self, shift: Union[vec3, List[float], Tuple[float, float, float]]):
594 """
595 Translate all hit points by a shift vector.
596
597 Args:
598 shift: Translation vector (vec3 or 3-element list)
599 """
600 if isinstance(shift, (list, tuple)):
601 if len(shift) != 3:
602 raise ValueError("Shift must have 3 elements [x, y, z]")
603 shift_list = list(shift)
604 elif hasattr(shift, 'x'):
605 shift_list = [shift.x, shift.y, shift.z]
606 else:
607 raise ValueError("Shift must be vec3 or 3-element list/tuple")
608
609 lidar_wrapper.lidarCoordinateShift(self._cloud_ptr, shift_list)
610
611 def coordinateRotation(self, rotation: Union[SphericalCoord, List[float], Tuple[float, float]]):
612 """
613 Rotate all hit points by spherical rotation angles.
614
615 Args:
616 rotation: Rotation angles (SphericalCoord or 2-3 element list)
617 """
618 if isinstance(rotation, (list, tuple)):
619 if len(rotation) < 2:
620 raise ValueError("Rotation must have at least 2 elements [radius, elevation]")
621 rotation_list = list(rotation)
622 elif hasattr(rotation, 'radius'):
623 rotation_list = [rotation.radius, rotation.elevation, rotation.azimuth]
624 else:
625 raise ValueError("Rotation must be SphericalCoord or 2-3 element list")
626
627 lidar_wrapper.lidarCoordinateRotation(self._cloud_ptr, rotation_list)
628
629 def triangulateHitPoints(self, Lmax: float, max_aspect_ratio: float = 4.0):
630 """
631 Generate triangle mesh from hit points using Delaunay triangulation.
632
633 Args:
634 Lmax: Maximum triangle edge length
635 max_aspect_ratio: Maximum triangle aspect ratio (default 4.0)
636 """
637 validate_positive_value(Lmax, 'Lmax', 'triangulateHitPoints')
638 validate_positive_value(max_aspect_ratio, 'max_aspect_ratio', 'triangulateHitPoints')
639 lidar_wrapper.lidarTriangulateHitPoints(self._cloud_ptr, Lmax, max_aspect_ratio)
640
641 def getTriangleCount(self) -> int:
642 """Get number of triangles in the mesh"""
643 return lidar_wrapper.getLiDARTriangleCount(self._cloud_ptr)
644
645 def getTriangulationStats(self) -> dict:
646 """Filter diagnostics from the most recent triangulateHitPoints() call.
647
648 Returns a dict::
649
650 {"candidates", "dropped_lmax", "dropped_aspect", "dropped_degenerate"}
651
652 Each dropped triangle is attributed to one primary reason (Lmax, then
653 aspect, then degenerate), so ``candidates == getTriangleCount() +
654 dropped_lmax + dropped_aspect + dropped_degenerate``. All zero if
655 triangulation has not been run. Use this to tell whether an empty or
656 sparse mesh is data-limited (few candidates) or filter-limited (many
657 candidates dropped by Lmax/aspect).
658 """
659 return lidar_wrapper.getLiDARTriangulationStats(self._cloud_ptr)
660
661 def getTriangleVerticesAll(self):
662 """Bulk-export every triangle's vertices and source scan in one call.
663
664 Returns (xyz_flat, scan_ids): xyz_flat is a (T*9,) float32 array laid out
665 [v0x,v0y,v0z, v1x,v1y,v1z, v2x,v2y,v2z] per triangle, scan_ids is a (T,)
666 int32 array. Avoids the Context round-trip and the per-triangle
667 getPrimitiveVertices loop.
668 """
669 return lidar_wrapper.getLiDARTriangleVertices_all(
670 self._cloud_ptr, self.getTriangleCount())
671
672 def distanceFilter(self, maxdistance: float):
673 """Filter hit points by maximum distance from scanner"""
674 validate_positive_value(maxdistance, 'maxdistance', 'distanceFilter')
675 lidar_wrapper.lidarDistanceFilter(self._cloud_ptr, maxdistance)
676
677 def reflectanceFilter(self, minreflectance: float):
678 """Filter hit points by minimum reflectance value"""
679 lidar_wrapper.lidarReflectanceFilter(self._cloud_ptr, minreflectance)
680
681 def firstHitFilter(self):
682 """Keep only first return hit points"""
683 lidar_wrapper.lidarFirstHitFilter(self._cloud_ptr)
684
685 def lastHitFilter(self):
686 """Keep only last return hit points"""
687 lidar_wrapper.lidarLastHitFilter(self._cloud_ptr)
688
689 def exportPointCloud(self, filename: str, write_header: bool = True):
690 """Export point cloud to ASCII file.
691
692 Args:
693 filename: Output file path.
694 write_header: If True (default), prepend a ``#``-prefixed comment line listing the
695 column field names (CloudCompare convention). The loader skips ``#``-prefixed
696 lines, so headered files round-trip through ``loadXML()``. Set False for a
697 bare data file.
698 """
699 if not filename:
700 raise ValueError("Filename cannot be empty")
701 lidar_wrapper.exportLiDARPointCloud(self._cloud_ptr, filename, write_header)
702
703 def exportLeafAreaUncertainty(self, filename: str):
704 """Export per-voxel leaf-area sampling uncertainty to a self-describing ASCII file.
705
706 The file has a ``#``-prefixed header and one row per grid cell:
707 ``cell_index leaf_area beam_count I_rdi LAD_std_error ci_valid``. Requires that
708 :meth:`calculateLeafArea` has been run with an ``element_width`` (the uncertainty
709 overload).
710 """
711 if not filename:
712 raise ValueError("Filename cannot be empty")
713 lidar_wrapper.exportLiDARLeafAreaUncertainty(self._cloud_ptr, filename)
714
715 def exportScans(self, filename: str):
716 """Export all scans to an XML metadata file plus one ASCII data file per scan.
717
718 Args:
719 filename: Path of the XML metadata file to write (e.g. "output/scans.xml").
720 One ASCII data file is auto-generated per scan, named by stripping the XML
721 extension and appending "_<scanID>.xyz" (e.g. "output/scans_0.xyz"). The
722 resulting XML can be re-loaded with loadXML() from the same working directory.
723 """
724 if not filename:
725 raise ValueError("Filename cannot be empty")
726 lidar_wrapper.exportLiDARScans(self._cloud_ptr, filename)
727
728 def loadXML(self, filename: str):
729 """Load scan metadata from XML file"""
730 if not filename:
731 raise ValueError("Filename cannot be empty")
732 lidar_wrapper.loadLiDARXML(self._cloud_ptr, filename)
733
734 def disableMessages(self):
735 """Disable console output messages"""
736 lidar_wrapper.lidarDisableMessages(self._cloud_ptr)
737
738 def enableMessages(self):
739 """Enable console output messages"""
740 lidar_wrapper.lidarEnableMessages(self._cloud_ptr)
741
742 def addGrid(self, center: Union[vec3, List[float], Tuple[float, float, float]],
743 size: Union[vec3, List[float], Tuple[float, float, float]],
744 ndiv: Union[List[int], Tuple[int, int, int]],
745 rotation: float = 0.0):
746 """
747 Add a rectangular grid of voxel cells.
748
749 Args:
750 center: Grid center position (vec3 or 3-element list)
751 size: Grid dimensions [x, y, z] (vec3 or 3-element list)
752 ndiv: Number of divisions [nx, ny, nz] (3-element list)
753 rotation: Azimuthal rotation angle (radians, default 0.0)
754
755 Example:
756 >>> lidar.addGrid(
757 ... center=vec3(0, 0, 0.5),
758 ... size=vec3(10, 10, 1),
759 ... ndiv=[10, 10, 5],
760 ... rotation=0.0
761 ... )
762 """
763 # Convert center to list
764 if isinstance(center, (list, tuple)):
765 if len(center) != 3:
766 raise ValueError("Center must have 3 elements [x, y, z]")
767 center_list = list(center)
768 elif hasattr(center, 'x'):
769 center_list = [center.x, center.y, center.z]
770 else:
771 raise ValueError("Center must be vec3 or 3-element list/tuple")
772
773 # Convert size to list
774 if isinstance(size, (list, tuple)):
775 if len(size) != 3:
776 raise ValueError("Size must have 3 elements [x, y, z]")
777 size_list = list(size)
778 elif hasattr(size, 'x'):
779 size_list = [size.x, size.y, size.z]
780 else:
781 raise ValueError("Size must be vec3 or 3-element list/tuple")
782
783 # Validate ndiv
784 if not isinstance(ndiv, (list, tuple)) or len(ndiv) != 3:
785 raise ValueError("Ndiv must be a 3-element list [nx, ny, nz]")
786
787 lidar_wrapper.addLiDARGrid(self._cloud_ptr, center_list, size_list, list(ndiv), rotation)
788
789 def addGridCell(self, center: Union[vec3, List[float], Tuple[float, float, float]],
790 size: Union[vec3, List[float], Tuple[float, float, float]],
791 rotation: float = 0.0):
792 """
793 Add a single grid cell.
794
795 Args:
796 center: Cell center position (vec3 or 3-element list)
797 size: Cell dimensions [x, y, z] (vec3 or 3-element list)
798 rotation: Azimuthal rotation angle (radians, default 0.0)
799 """
800 # Convert center to list
801 if isinstance(center, (list, tuple)):
802 if len(center) != 3:
803 raise ValueError("Center must have 3 elements [x, y, z]")
804 center_list = list(center)
805 elif hasattr(center, 'x'):
806 center_list = [center.x, center.y, center.z]
807 else:
808 raise ValueError("Center must be vec3 or 3-element list/tuple")
809
810 # Convert size to list
811 if isinstance(size, (list, tuple)):
812 if len(size) != 3:
813 raise ValueError("Size must have 3 elements [x, y, z]")
814 size_list = list(size)
815 elif hasattr(size, 'x'):
816 size_list = [size.x, size.y, size.z]
817 else:
818 raise ValueError("Size must be vec3 or 3-element list/tuple")
819
820 lidar_wrapper.addLiDARGridCell(self._cloud_ptr, center_list, size_list, rotation)
821
822 def getGridCellCount(self) -> int:
823 """Get total number of grid cells"""
824 return lidar_wrapper.getLiDARGridCellCount(self._cloud_ptr)
825
826 def getCellCenter(self, index: int) -> vec3:
827 """Get center position of a grid cell"""
828 if index < 0:
829 raise ValueError("Index must be non-negative")
830 center_list = lidar_wrapper.getLiDARCellCenter(self._cloud_ptr, index)
831 return vec3(*center_list)
832
833 def getCellSize(self, index: int) -> vec3:
834 """Get size of a grid cell"""
835 if index < 0:
836 raise ValueError("Index must be non-negative")
837 size_list = lidar_wrapper.getLiDARCellSize(self._cloud_ptr, index)
838 return vec3(*size_list)
839
840 def getCellLeafArea(self, index: int) -> float:
841 """Get leaf area of a grid cell (m²)"""
842 if index < 0:
843 raise ValueError("Index must be non-negative")
844 return lidar_wrapper.getLiDARCellLeafArea(self._cloud_ptr, index)
845
846 def getCellLeafAreaDensity(self, index: int) -> float:
847 """Get leaf area density of a grid cell (m²/m³)"""
848 if index < 0:
849 raise ValueError("Index must be non-negative")
850 return lidar_wrapper.getLiDARCellLeafAreaDensity(self._cloud_ptr, index)
851
852 def getCellBeamCount(self, index: int) -> int:
853 """Get the beam count N that entered a grid cell during the leaf-area inversion.
854
855 Returns -1 if :meth:`calculateLeafArea` has not been run for this cell.
856 """
857 if index < 0:
858 raise ValueError("Index must be non-negative")
859 return lidar_wrapper.getLiDARCellBeamCount(self._cloud_ptr, index)
860
861 def getCellRelativeDensityIndex(self, index: int) -> float:
862 """Get the relative density index (I_rdi) for a grid cell."""
863 if index < 0:
864 raise ValueError("Index must be non-negative")
865 return lidar_wrapper.getLiDARCellRelativeDensityIndex(self._cloud_ptr, index)
866
867 def getCellMeanPathLength(self, index: int) -> float:
868 """Get the mean beam path length (m) through a grid cell."""
869 if index < 0:
870 raise ValueError("Index must be non-negative")
871 return lidar_wrapper.getLiDARCellMeanPathLength(self._cloud_ptr, index)
872
873 def getCellLADVariance(self, index: int) -> float:
874 """Get the per-voxel LAD sampling variance for a grid cell.
875
876 Returns -1 if uncertainty has not been computed (call :meth:`calculateLeafArea`
877 with an ``element_width``).
878 """
879 if index < 0:
880 raise ValueError("Index must be non-negative")
881 return lidar_wrapper.getLiDARCellLADVariance(self._cloud_ptr, index)
882
883 def getCellLeafAreaConfidenceInterval(self, index: int, confidence_level: float = 0.95):
884 """Get the leaf-area confidence interval for a single grid cell.
885
886 Returns a ``(valid, lower, upper)`` tuple. ``valid`` is False when the interval is
887 gated out by the Pimont validity envelope (single-voxel intervals are often
888 untrustworthy; prefer :meth:`getGroupLADConfidenceInterval`). Requires
889 :meth:`calculateLeafArea` to have been run with an ``element_width``.
890 """
891 if index < 0:
892 raise ValueError("Index must be non-negative")
893 return lidar_wrapper.getLiDARCellLeafAreaConfidenceInterval(
894 self._cloud_ptr, index, confidence_level)
895
896 def getGroupLADConfidenceInterval(self, indices: List[int], confidence_level: float = 0.95):
897 """Get the group-scale LAD confidence interval over a set of grid cells (recommended).
898
899 Returns a ``(valid, mean_lad, lower, upper)`` tuple (Pimont et al. 2018, Eq. 39,
900 assuming voxel independence). Requires :meth:`calculateLeafArea` to have been run
901 with an ``element_width``.
902 """
903 if not indices:
904 raise ValueError("indices must contain at least one cell index")
905 if any(i < 0 for i in indices):
906 raise ValueError("Cell indices must be non-negative")
907 return lidar_wrapper.getLiDARGroupLADConfidenceInterval(
908 self._cloud_ptr, indices, confidence_level)
909
910 def getCellGtheta(self, index: int) -> float:
911 """Get G(theta) value for a grid cell"""
912 if index < 0:
913 raise ValueError("Index must be non-negative")
914 return lidar_wrapper.getLiDARCellGtheta(self._cloud_ptr, index)
915
916 def setCellGtheta(self, Gtheta: float, index: int):
917 """Set G(theta) value for a grid cell"""
918 if index < 0:
919 raise ValueError("Index must be non-negative")
920 lidar_wrapper.setLiDARCellGtheta(self._cloud_ptr, Gtheta, index)
921
922 def calculateHitGridCell(self):
923 """Calculate hit point grid cell assignments"""
924 lidar_wrapper.calculateLiDARHitGridCell(self._cloud_ptr)
925
926 def gapfillMisses(self):
927 """
928 Gapfill sky/miss points where rays didn't hit geometry.
929
930 Important for accurate leaf area calculations with real LiDAR data.
931 Should be called before triangulation when processing real data.
932 """
933 lidar_wrapper.gapfillLiDARMisses(self._cloud_ptr)
934
935 def syntheticScan(self, context: Context,
936 rays_per_pulse: Optional[int] = None,
937 pulse_distance_threshold: Optional[float] = None,
938 scan_grid_only: bool = False,
939 record_misses: bool = True,
940 append: bool = False):
941 """
942 Perform synthetic LiDAR scan of geometry in Context.
943
944 Requires scan metadata to be defined first via addScan() or loadXML().
945 Uses ray tracing to simulate LiDAR instrument measurements.
946
947 Args:
948 context: Helios Context containing geometry to scan
949 rays_per_pulse: Number of rays per pulse (None=discrete-return, typical: 100)
950 pulse_distance_threshold: Distance threshold for aggregating hits (meters, required for waveform)
951 scan_grid_only: If True, only scan within defined grid cells
952 record_misses: If True, record miss/sky points where rays don't hit geometry
953 append: If True, append to existing hits; if False, clear existing hits
954
955 Example (Discrete-return):
956 >>> from pyhelios import Context, LiDARCloud
957 >>> from pyhelios.types import vec3
958 >>> with Context() as context:
959 ... # Add geometry
960 ... context.addPatch(center=vec3(0, 0, 0.5), size=vec2(1, 1))
961 ...
962 ... with LiDARCloud() as lidar:
963 ... # Define scan parameters
964 ... scan_id = lidar.addScan(
965 ... origin=vec3(0, 0, 2),
966 ... Ntheta=100, theta_range=(0, 1.57),
967 ... Nphi=100, phi_range=(0, 6.28),
968 ... exit_diameter=0, beam_divergence=0
969 ... )
970 ...
971 ... # Perform discrete-return scan
972 ... lidar.syntheticScan(context)
973
974 Example (Full-waveform):
975 >>> lidar.syntheticScan(
976 ... context,
977 ... rays_per_pulse=100,
978 ... pulse_distance_threshold=0.02,
979 ... record_misses=True
980 ... )
981 """
982 if not isinstance(context, Context):
983 raise TypeError("context must be a Context instance")
984
985 context_ptr = context.getNativePtr()
986
987 # Discrete-return mode (single ray per pulse)
988 if rays_per_pulse is None:
989 # Honor scan_grid_only and record_misses for discrete scans too. record_misses
990 # defaults to True so the cloud carries the transmitted beams that
991 # calculateLeafArea() requires.
992 lidar_wrapper.syntheticLiDARScanDiscrete(
993 self._cloud_ptr, context_ptr, scan_grid_only, record_misses, append)
994 else:
995 # Full-waveform mode (multiple rays per pulse)
996 if pulse_distance_threshold is None:
997 raise ValueError("pulse_distance_threshold required for full-waveform scanning")
998
999 validate_positive_value(rays_per_pulse, 'rays_per_pulse', 'syntheticScan')
1000 validate_positive_value(pulse_distance_threshold, 'pulse_distance_threshold', 'syntheticScan')
1001
1002 lidar_wrapper.syntheticLiDARScanFull(
1003 self._cloud_ptr, context_ptr,
1004 rays_per_pulse, pulse_distance_threshold,
1005 scan_grid_only, record_misses, append
1006 )
1007
1008 def calculateLeafArea(self, context: Context, min_voxel_hits: Optional[int] = None,
1009 element_width: Optional[float] = None):
1010 """
1011 Calculate leaf area for each grid cell.
1012
1013 Requires triangulation to have been performed first.
1014
1015 .. note::
1016 The cloud must contain misses (transmitted beams that returned nothing) — the
1017 inversion fails fast without them. Misses are produced by
1018 ``syntheticScan(..., record_misses=True)`` (the default) or by
1019 :meth:`gapfillMisses`. Use :meth:`hasMisses` to check.
1020
1021 Args:
1022 context: Helios Context instance
1023 min_voxel_hits: Optional minimum number of hits required per voxel
1024 element_width: Optional characteristic vegetation element width (meters). When
1025 provided, per-voxel sampling uncertainty (Pimont et al. 2018) is computed
1026 alongside the leaf-area estimate and becomes available via
1027 :meth:`getCellLADVariance`, :meth:`getCellLeafAreaConfidenceInterval`, and
1028 :meth:`getGroupLADConfidenceInterval`. ``element_width <= 0`` yields a
1029 sampling-only variance. Requires ``min_voxel_hits`` to also be specified.
1030
1031 Example:
1032 >>> from pyhelios import Context, LiDARCloud
1033 >>> with Context() as context:
1034 ... with LiDARCloud() as lidar:
1035 ... # ... load data, add grid, triangulate ...
1036 ... lidar.calculateLeafArea(context)
1037 """
1038 if not isinstance(context, Context):
1039 raise TypeError("context must be a Context instance")
1040
1041 context_ptr = context.getNativePtr()
1042 if element_width is not None:
1043 if min_voxel_hits is None:
1044 raise ValueError(
1045 "element_width requires min_voxel_hits to also be specified "
1046 "(the uncertainty overload takes both)")
1047 lidar_wrapper.calculateLiDARLeafAreaUncertainty(
1048 self._cloud_ptr, context_ptr, min_voxel_hits, element_width)
1049 elif min_voxel_hits is None:
1050 lidar_wrapper.calculateLiDARLeafArea(self._cloud_ptr, context_ptr)
1051 else:
1052 lidar_wrapper.calculateLiDARLeafAreaMinHits(self._cloud_ptr, context_ptr, min_voxel_hits)
1053
1054 def calculateSyntheticLeafArea(self, context: Context):
1055 """
1056 Calculate synthetic leaf area (for validation of synthetic scans).
1057
1058 Uses exact primitive geometry to calculate leaf area, useful for
1059 validating synthetic scan accuracy.
1060
1061 Args:
1062 context: Helios Context instance containing primitive geometry
1063 """
1064 if not isinstance(context, Context):
1065 raise TypeError("context must be a Context instance")
1066 context_ptr = context.getNativePtr()
1067 lidar_wrapper.calculateSyntheticLiDARLeafArea(self._cloud_ptr, context_ptr)
1068
1069 def calculateSyntheticGtheta(self, context: Context):
1070 """
1071 Calculate synthetic G(theta) (for validation of synthetic scans).
1072
1073 Uses exact primitive geometry to calculate G(theta), useful for
1074 validating synthetic scan accuracy.
1075
1076 Args:
1077 context: Helios Context instance containing primitive geometry
1078 """
1079 if not isinstance(context, Context):
1080 raise TypeError("context must be a Context instance")
1081 context_ptr = context.getNativePtr()
1082 lidar_wrapper.calculateSyntheticLiDARGtheta(self._cloud_ptr, context_ptr)
1083
1084 def exportTriangleNormals(self, filename: str):
1085 """Export triangle normal vectors to file"""
1086 if not filename:
1087 raise ValueError("Filename cannot be empty")
1088 lidar_wrapper.exportLiDARTriangleNormals(self._cloud_ptr, filename)
1089
1090 def exportTriangleAreas(self, filename: str):
1091 """Export triangle areas to file"""
1092 if not filename:
1093 raise ValueError("Filename cannot be empty")
1094 lidar_wrapper.exportLiDARTriangleAreas(self._cloud_ptr, filename)
1095
1096 def exportLeafAreas(self, filename: str):
1097 """Export leaf areas for each grid cell to file"""
1098 if not filename:
1099 raise ValueError("Filename cannot be empty")
1100 lidar_wrapper.exportLiDARLeafAreas(self._cloud_ptr, filename)
1101
1102 def exportLeafAreaDensities(self, filename: str):
1103 """Export leaf area densities for each grid cell to file"""
1104 if not filename:
1105 raise ValueError("Filename cannot be empty")
1106 lidar_wrapper.exportLiDARLeafAreaDensities(self._cloud_ptr, filename)
1107
1108 def exportGtheta(self, filename: str):
1109 """Export G(theta) values for each grid cell to file"""
1110 if not filename:
1111 raise ValueError("Filename cannot be empty")
1112 lidar_wrapper.exportLiDARGtheta(self._cloud_ptr, filename)
1113
1114 def addTrianglesToContext(self, context: Context):
1115 """
1116 Add triangulated mesh to Context as triangle primitives.
1117
1118 Converts the triangulated point cloud mesh into Context triangle
1119 primitives that can be used for further analysis or visualization.
1120
1121 Args:
1122 context: Helios Context instance
1123
1124 Example:
1125 >>> with Context() as context:
1126 ... with LiDARCloud() as lidar:
1127 ... lidar.loadXML("scan.xml")
1128 ... lidar.triangulateHitPoints(Lmax=0.5, max_aspect_ratio=5)
1129 ... lidar.addTrianglesToContext(context)
1130 ... print(f"Added {context.getPrimitiveCount()} triangles to context")
1131 """
1132 if not isinstance(context, Context):
1133 raise TypeError("context must be a Context instance")
1134 lidar_wrapper.addLiDARTrianglesToContext(self._cloud_ptr, context.getNativePtr())
1135
1136 def initializeCollisionDetection(self, context: Context):
1137 """
1138 Initialize CollisionDetection plugin for ray tracing.
1139
1140 Required before performing synthetic scans.
1141
1142 Args:
1143 context: Helios Context instance containing geometry
1144 """
1145 if not isinstance(context, Context):
1146 raise TypeError("context must be a Context instance")
1147 lidar_wrapper.initializeLiDARCollisionDetection(self._cloud_ptr, context.getNativePtr())
1148
1149 def enableCDGPUAcceleration(self):
1150 """Enable GPU acceleration for collision detection ray tracing"""
1151 lidar_wrapper.enableLiDARCDGPUAcceleration(self._cloud_ptr)
1152
1153 def disableCDGPUAcceleration(self):
1154 """Disable GPU acceleration (use CPU ray tracing)"""
1155 lidar_wrapper.disableLiDARCDGPUAcceleration(self._cloud_ptr)
1156
1157 def is_available(self) -> bool:
1158 """
1159 Check if LiDAR is available in current build.
1160
1161 Returns:
1162 True if plugin is available, False otherwise
1163 """
1164 registry = get_plugin_registry()
1165 return registry.is_plugin_available('lidar')
1166
1167
1168# Convenience function
1169def create_lidar_cloud() -> LiDARCloud:
1170 """
1171 Create LiDARCloud instance.
1172
1173 Returns:
1174 LiDARCloud instance
1175 """
1176 return LiDARCloud()
High-level interface for LiDAR point cloud operations.
Definition LiDARCloud.py:63
__init__(self)
Initialize LiDARCloud.
Definition LiDARCloud.py:72
__enter__(self)
Context manager entry.
Definition LiDARCloud.py:90
__exit__(self, exc_type, exc_val, exc_tb)
Context manager exit - cleanup resources.
Definition LiDARCloud.py:94
__del__(self)
Fallback destructor for cleanup without context manage.
Exception raised for LiDAR-specific errors.
Definition LiDARCloud.py:23
Scan pattern returned by :meth:LiDARCloud.getScanPattern.
Definition LiDARCloud.py:33
Exception classes for PyHelios library.
Definition exceptions.py:10