pysdic.triangle_3_cast_rays#

triangle_3_cast_rays(vertices_coordinates, connectivity, ray_origins, ray_directions, nan_open3d_errors=True)[source]#

Cast rays into a triangular mesh and compute the intersection points.

This function uses Open3D to perform ray-mesh intersection tests.

Note

The inputs vertices_coordinates, ray_origins, and ray_directions will be converted to numpy.float64 for computation. The input connectivity will be converted to numpy.int64 for computation.

Warning

No tests are performed to check if the inputs array content are consistent (e.g., if the connectivity contains invalid vertex indices, if the rays are valid, …). Only the shapes and types of the input arrays are validated. The behavior of the function is undefined if the input arrays are not consistent.

Parameters:
  • vertices_coordinates (ArrayLike) – An array of shape \((N_v, 3)\) representing the coordinates of the vertices.

  • connectivity (ArrayLike) – An array of shape \((N_e, 3)\) representing the connectivity of the triangular elements.

  • ray_origins (ArrayLike) – An array of shape \((N_r, 3)\) representing the origins of the rays.

  • ray_directions (ArrayLike) – An array of shape \((N_r, 3)\) representing the directions of the rays.

  • nan_open3d_errors (bool, optional) – If True, handle NaN errors due to Open3D float32 precision issues by setting invalid intersections to NaN and -1, by default True.

Returns:

  • natural_coordinates (numpy.ndarray) – An array of shape \((N_r, 2)\) representing the natural coordinates (\(\xi, \eta\)) of the intersection points.

  • element_indices (numpy.ndarray) – An array of shape \((N_r,)\) representing the indices of the intersected elements. If a ray does not intersect any element, the index is -1.

Raises:
  • TypeError – If the inputs are not numpy arrays.

  • ValueError – If the input arrays do not have the correct shapes.

Return type:

Tuple[ndarray, ndarray]

See also

pysdic.triangle_3_mesh_to_open3d

For converting a Mesh instance to an Open3D TriangleMesh object.

Examples

Casting rays into a triangular mesh:

 1import numpy as np
 2from pysdic import triangle_3_cast_rays
 3
 4vertices_coordinates = np.array([[0, 0, 0],
 5                     [1, 0, 0],
 6                     [0, 1, 0],
 7                     [0, 0, 1]])  # shape (Nv, E) = (4, 3)
 8
 9connectivity = np.array([[0, 1, 2],
10                         [0, 1, 3],
11                         [0, 2, 3],
12                         [1, 2, 3]])  # shape (Ne, Nvpe) = (4, 3)
13
14ray_origins = np.array([[0.1, 0.1, -1],
15                        [0.5, 0.5, -1]])  # shape (Nr, E) = (2, 3)
16
17ray_directions = np.array([[0, 0, 1],
18                           [0, 0, 1]])  # shape (Nr, E) = (2, 3)
19
20natural_coordinates, element_indices = triangle_3_cast_rays(vertices_coordinates, connectivity, ray_origins, ray_directions)
21print(natural_coordinates) # shape (Nr, K)
22print(element_indices) # shape (Nr,)
[[0.1 0.1]
 [0.5 0.5]]

[0 0]

Deduce the coordinates of the intersection points from the natural coordinates and the element connectivity.

 1from pysdic import triangle_3_shape_functions, interpolate_property
 2
 3shape_functions = triangle_3_shape_functions(natural_coordinates)  # shape (Nr, Nvpe) = (2, 3)
 4
 5points_coordinates = interpolate_property(
 6    vertices_coordinates,
 7    shape_functions,
 8    connectivity,
 9    element_indices
10)  # shape (Nr, E) = (2, 3)