pyRTX.core.utils_rt
Functions
|
Extract and process intersection results from an Embree RayHit structure. |
|
Initialize mesh geometry for the Embree 3 ray tracing kernel. |
|
Initialize Embree 3 RayHit data structure. |
|
Main ray tracing kernel wrapper. |
|
Perform an element-wise dot product between two arrays of vectors. |
Compute unit vectors for a block of vectors. |
|
|
Initialize mesh geometry for the CGAL ray tracing kernel. |
|
Divide an iterator or array into approximately equal-sized chunks for parallel processing. |
|
Prepare ray origins and directions for subsequent bounces. |
|
Compute multiple diffuse reflection directions for an array of surface normals. |
|
Efficiently build a pixel array coordinate grid using Numba. |
|
Return a cached EmbreeTrimeshShapeModel for this mesh, creating if needed. |
|
Compute the geometric centroids of all triangular faces in a mesh. |
|
Compute cross products of edge vectors for all triangular faces in a mesh. |
|
Compute the areas of all triangular faces in a mesh. |
Generate a unit vector orthogonal to the input vector. |
|
|
Compute outward-pointing unit normal vectors for all triangular faces. |
Efficiently compute both surface normals and face areas simultaneously. |
|
|
Generate a rectangular grid of rays for ray tracing. |
|
Generate a rectangular pixel array (optimized version). |
|
Convert a 3x3 SPICE rotation matrix to a 4x4 homogeneous transformation matrix. |
|
Compute reflected ray directions. |
|
Generate direction vectors following the Lambert cosine distribution. |
|
Save ray tracing results to a pickled dictionary for visualization. |
Classes
|
A triangle mesh shape model that uses the CGAL AABB tree for ray tracing. |
|
A triangle mesh shape model that uses the Embree ray tracing kernel. |
An abstract base class for shape models. |
|
|
A shape model consisting of a single triangle mesh. |
- pyRTX.core.utils_rt.get_cached_embree_scene(mesh_obj)[source]
Return a cached EmbreeTrimeshShapeModel for this mesh, creating if needed.
- Parameters:
mesh_obj (
trimesh.Trimesh) – The mesh object to get a scene for.- Returns:
The cached or newly created Embree shape model.
- Return type:
- pyRTX.core.utils_rt.chunker(iterator, chunks)[source]
Divide an iterator or array into approximately equal-sized chunks for parallel processing.
- pyRTX.core.utils_rt.pxform_convert(pxform)[source]
Convert a 3x3 SPICE rotation matrix to a 4x4 homogeneous transformation matrix.
- Parameters:
pxform (
array_like,shape (3,3)) – The 3x3 rotation matrix from SPICE.- Returns:
The 4x4 homogeneous transformation matrix.
- Return type:
ndarray,shape (4,4)
- pyRTX.core.utils_rt.block_normalize(V)[source]
Compute unit vectors for a block of vectors.
- Parameters:
V (
ndarray,shape (N,3)or(3,)) – Array of vectors to normalize.- Returns:
The normalized vectors.
- Return type:
ndarray
- pyRTX.core.utils_rt.block_dot(a, b)[source]
Perform an element-wise dot product between two arrays of vectors.
- Parameters:
a (
ndarray,shape (N,m)) – First array of vectors.b (
ndarray,shape (N,m)) – Second array of vectors.
- Returns:
Array of dot products.
- Return type:
ndarray,shape (N,)
- pyRTX.core.utils_rt.pixel_plane(d0, lon, lat, width=1, height=1, ray_spacing=0.1)[source]
Generate a rectangular grid of rays for ray tracing.
- Parameters:
d0 (
float) – Distance of the pixel plane from the origin.lon (
float) – Longitude of the pixel plane’s center direction in radians.lat (
float) – Latitude of the pixel plane’s center direction in radians.width (
float, default1) – Width of the plane.height (
float, default1) – Height of the plane.ray_spacing (
float, default0.1) – Spacing between adjacent rays.
- Returns:
A tuple containing: - locs (ndarray): Ray origin positions. - dirs (ndarray): Ray direction unit vectors.
- Return type:
- pyRTX.core.utils_rt.fast_vector_build(linsp1, linsp2, dim1, dim2)[source]
Efficiently build a pixel array coordinate grid using Numba.
- pyRTX.core.utils_rt.pixel_plane_opt(d0, lon, lat, width=1, height=1, ray_spacing=0.1, packets=1)[source]
Generate a rectangular pixel array (optimized version).
- Parameters:
d0 (
float) – Distance of the pixel plane from the origin.lon (
float) – Longitude of the pixel plane’s center direction in radians.lat (
float) – Latitude of the pixel plane’s center direction in radians.width (
float, default1) – Width of the plane.height (
float, default1) – Height of the plane.ray_spacing (
float, default0.1) – Spacing between adjacent rays.packets (
int, default1) – Number of ray packets to subdivide the rays into.
- Returns:
A tuple containing: - locs (ndarray): Ray origin positions. - dirs (ndarray): Ray direction unit vectors.
- Return type:
- pyRTX.core.utils_rt.reflected(incoming, normal)[source]
Compute reflected ray directions.
- Parameters:
incoming (
ndarray) – Incoming ray direction vectors.normal (
ndarray) – Surface normal vectors.
- Returns:
Reflected ray direction vectors.
- Return type:
ndarray
- pyRTX.core.utils_rt.get_orthogonal(v)[source]
Generate a unit vector orthogonal to the input vector.
- Parameters:
v (
ndarray) – Input vector.- Returns:
Orthogonal unit vector.
- Return type:
ndarray
- pyRTX.core.utils_rt.sample_lambert_dist(normal, num=100)[source]
Generate direction vectors following the Lambert cosine distribution.
- Parameters:
normal (
ndarray) – Surface normal vector.num (
int, default100) – Number of sample directions to generate.
- Returns:
Array of sampled direction vectors.
- Return type:
ndarray
- pyRTX.core.utils_rt.diffuse(normals, num=10)[source]
Compute multiple diffuse reflection directions for an array of surface normals.
- Parameters:
normals (
ndarray) – Array of surface normal unit vectors.num (
int, default10) – Number of diffuse samples to generate for each normal.
- Returns:
Array of sampled diffuse direction vectors.
- Return type:
ndarray
- pyRTX.core.utils_rt.compute_secondary_bounce(location, index_tri, mesh_obj, ray_directions, index_ray, diffusion=False, num_diffuse=None)[source]
Prepare ray origins and directions for subsequent bounces.
- Parameters:
location (
ndarray) – 3D coordinates of intersection points.index_tri (
ndarray) – Indices of the intersected mesh faces.mesh_obj (
trimesh.Trimesh) – The mesh object.ray_directions (
ndarray) – Direction vectors of the incident rays.index_ray (
ndarray) – Indices of the rays that intersected the mesh.diffusion (
bool, defaultFalse) – If True, compute diffuse reflection directions.num_diffuse (
int, optional) – Number of diffuse samples per intersection point.
- Returns:
A tuple containing: - location (ndarray): Intersection point coordinates. - reflect_dirs (ndarray): Specularly reflected ray directions. - diffuse_dirs (ndarray or int): Diffusely reflected directions or -1.
- Return type:
- pyRTX.core.utils_rt.save_for_visualization(outputFilePath, mesh, ray_origins, ray_directions, location, index_tri, diffusion_pack)[source]
Save ray tracing results to a pickled dictionary for visualization.
- Parameters:
outputFilePath (
str) – Path to the output file.mesh (
trimesh.Trimesh) – The mesh object.ray_origins (
listofndarrays) – List of ray origin arrays for each bounce.ray_directions (
listofndarrays) – List of ray direction arrays for each bounce.location (
listofndarrays) – List of intersection point coordinates for each bounce.index_tri (
listofndarrays) – List of indices of intersected triangles for each bounce.
- pyRTX.core.utils_rt.Embree3_init_geometry(mesh_obj)[source]
Initialize mesh geometry for the Embree 3 ray tracing kernel.
- Parameters:
mesh_obj (
trimesh.Trimesh) – Input mesh object.- Returns:
The Embree shape model.
- Return type:
- pyRTX.core.utils_rt.Embree3_init_rayhit(ray_origins, ray_directions)[source]
Initialize Embree 3 RayHit data structure.
- Parameters:
ray_origins (
ndarray) – Starting positions of rays.ray_directions (
ndarray) – Direction vectors of rays.
- Returns:
Initialized RayHit structure.
- Return type:
embree.RayHit1M
- pyRTX.core.utils_rt.Embree3_dump_solution(rayhit, V, F)[source]
Extract and process intersection results from an Embree RayHit structure.
- Parameters:
rayhit (
embree.RayHit1M) – The RayHit structure populated by Embree.V (
ndarray) – Vertex coordinates of the mesh.F (
ndarray) – Face indices of the mesh.
- Returns:
A tuple containing: - hits (ndarray or int): Indices of intersected triangles. - nhits (int): Number of intersected rays. - idh (ndarray or int): Indices of rays that hit the mesh. - Ph (ndarray or int): 3D coordinates of intersection points.
- Return type:
- pyRTX.core.utils_rt.cgal_init_geometry(mesh_obj)[source]
Initialize mesh geometry for the CGAL ray tracing kernel.
- Parameters:
mesh_obj (
trimesh.Trimesh) – Input mesh object.- Returns:
The CGAL shape model.
- Return type:
- pyRTX.core.utils_rt.get_centroids(V, F)[source]
Compute the geometric centroids of all triangular faces in a mesh.
- Parameters:
V (
ndarray) – Vertex coordinates of the mesh.F (
ndarray) – Face indices.
- Returns:
Centroid coordinates for each face.
- Return type:
ndarray
- pyRTX.core.utils_rt.get_cross_products(V, F)[source]
Compute cross products of edge vectors for all triangular faces in a mesh.
- Parameters:
V (
ndarray) – Vertex coordinates of the mesh.F (
ndarray) – Face indices.
- Returns:
Cross product vectors for each face.
- Return type:
ndarray
- pyRTX.core.utils_rt.get_face_areas(V, F)[source]
Compute the areas of all triangular faces in a mesh.
- Parameters:
V (
ndarray) – Vertex coordinates of the mesh.F (
ndarray) – Face indices.
- Returns:
Area of each face.
- Return type:
ndarray
- pyRTX.core.utils_rt.get_surface_normals(V, F)[source]
Compute outward-pointing unit normal vectors for all triangular faces.
- Parameters:
V (
ndarray) – Vertex coordinates of the mesh.F (
ndarray) – Face indices.
- Returns:
Unit normal vectors for each face.
- Return type:
ndarray
- pyRTX.core.utils_rt.get_surface_normals_and_face_areas(V, F)[source]
Efficiently compute both surface normals and face areas simultaneously.
- Parameters:
V (
ndarray) – Vertex coordinates of the mesh.F (
ndarray) – Face indices.
- Returns:
A tuple containing: - N (ndarray): Unit normal vectors for each face. - A (ndarray): Area of each face.
- Return type:
- class pyRTX.core.utils_rt.TrimeshShapeModel(V, F, N=None, P=None, A=None)[source]
A shape model consisting of a single triangle mesh.
- __init__(V, F, N=None, P=None, A=None)[source]
Initialize a triangle mesh shape model.
- Parameters:
V (
array_like) – An array with shape (num_verts, 3) of vertex coordinates.F (
array_like) – An array with shape (num_faces, 3) of face indices.N (
array_like, optional) – An array with shape (num_faces, 3) of face normals.P (
array_like, optional) – An array with shape (num_faces, 3) of triangle centroids.A (
array_like, optional) – An array of shape (num_faces,) of triangle areas.
- property num_faces
The number of faces in the mesh.
- property num_verts
The number of vertices in the mesh.
- intersect1(x, d)[source]
Trace a single ray.
- Parameters:
x (
array_like) – The origin of the ray.d (
array_like) – The direction of the ray.
- Returns:
A tuple containing the index of the hit and the distance to the hit.
- Return type:
- class pyRTX.core.utils_rt.CgalTrimeshShapeModel(V, F, N=None, P=None, A=None)[source]
A triangle mesh shape model that uses the CGAL AABB tree for ray tracing.
- class pyRTX.core.utils_rt.EmbreeTrimeshShapeModel(V, F, N=None, P=None, A=None)[source]
A triangle mesh shape model that uses the Embree ray tracing kernel.
- pyRTX.core.utils_rt.RTXkernel(mesh_obj, ray_origins, ray_directions, bounces=1, kernel='Embree3', diffusion=False, num_diffuse=None, errorMsg=True)[source]
Main ray tracing kernel wrapper.
- Parameters:
mesh_obj (
trimesh.Trimesh) – The mesh geometry to ray trace.ray_origins (
ndarray) – Starting positions of rays.ray_directions (
ndarray) – Direction vectors of rays.bounces (
int, default1) – Number of reflection bounces to simulate.kernel (
str, default'Embree3') – Ray tracing backend to use (‘Embree3’, ‘CGAL’, or ‘Native’).diffusion (
bool, defaultFalse) – If True, compute diffuse reflections.num_diffuse (
int, optional) – Number of diffuse samples per intersection.errorMsg (
bool, defaultTrue) – If True, print warnings when no intersections are found.
- Returns:
A tuple containing the results of the ray tracing.
- Return type: