Skip to content

squared_distance_to_element

squared_distance_to_element(point, V, element)

Squared distance from a point to a mesh element (point, edge, triangle)

Parameters:

Name Type Description Default
point (dim,) numpy double array

Query point coordinates

required
V (v,dim) numpy double array

Matrix of mesh/polyline/pointcloud vertex coordinates

required
element (s,) numpy int array

Vector of element indices into V

required

Returns:

Name Type Description
sqrD double

Squared minimum distance from point to mesh element

lmb (s,) numpy double array

Barycentric coordinates into the element of the closest point to the query point

See Also

squared_distance.

Examples:

# Generate random mesh
V = np.random.rand(3,3)
F = np.array([0,1,2],dtype=int)
# Generate random query point
P = np.random.rand(3)
# Calculate distance from point to triangle
sqrD = gpytoolbox.squared_distance_to_element(P,V,F)
Source code in src/gpytoolbox/squared_distance_to_element.py
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
def squared_distance_to_element(point,V,element):
    """Squared distance from a point to a mesh element (point, edge, triangle)

    Parameters
    ----------
    point : (dim,) numpy double array
        Query point coordinates
    V : (v,dim) numpy double array
        Matrix of mesh/polyline/pointcloud vertex coordinates
    element : (s,) numpy int array
        Vector of element indices into V

    Returns
    -------
    sqrD : double
        Squared minimum distance from point to mesh element
    lmb : (s,) numpy double array
        Barycentric coordinates into the element of the closest point to the query point

    See Also
    --------
    squared_distance.

    Examples
    --------
    ```python
    # Generate random mesh
    V = np.random.rand(3,3)
    F = np.array([0,1,2],dtype=int)
    # Generate random query point
    P = np.random.rand(3)
    # Calculate distance from point to triangle
    sqrD = gpytoolbox.squared_distance_to_element(P,V,F)
    ```
    """
    dim = V.shape[1]
    if element.ndim>1:
        element = np.ravel(element)
    simplex_size = element.shape[0]
    if simplex_size==1:
        # Then this is just distance between two points
        sqrD = np.sum((V[element] - point)**2.0)
        lmb = 1
    elif simplex_size==2:
        if dim==2:
            V = np.hstack(( V,np.zeros((V.shape[0],1)) ))
            point = np.concatenate((point,np.array([0])))
        # Distance from point to segment
        start = V[element[0],:]
        end = V[element[1],:]
        line_vec = end - start
        pnt_vec = point - start
        line_len = np.linalg.norm(line_vec)
        line_unitvec = line_vec/line_len
        pnt_vec_scaled = pnt_vec/line_len
        t = np.dot(line_unitvec, pnt_vec_scaled)    
        if t < 0.0:
            t = 0.0
        elif t > 1.0:
            t = 1.0
        nearest = t*line_vec
        sqrD = np.sum((nearest - pnt_vec)**2.0)
        nearest = start + nearest
        lmb = [1-t,t]
    elif simplex_size==3:
        assert(dim==3)
        sqrD = 0
        tri = np.vstack((
            V[element[0],:],
            V[element[1],:],
            V[element[2],:]
        ))
        # print(tri)
        sqrD,nearest_point = pointTriangleDistance(tri,point)
        lmb = barycentric_coordinates(nearest_point,V[element[0],:],V[element[1],:],V[element[2],:])
    return sqrD,lmb