Skip to content

initialize_aabbtree

initialize_aabbtree(V, F=None, vmin=None, vmax=None)

Axis-Aligned Bounding Box hierarchy for efficient computation

Computes an AABB tree by recursively dividing cells along the biggest dimension. May be called for any dimension and simplex size

Parameters:

Name Type Description Default
V numpy double array

Matrix of mesh vertex coordinates

required
F numpy int array, optional (default None)

Matrix of element indices into V (if None, treated as point cloud)

None

Returns:

Name Type Description
C numpy double array

Matrix of box centers

W numpy double array

Matrix of box widths

CH numpy int array

Matrix of child indeces (-1 if leaf node)

PAR numpy int array

Vector of immediate parent indeces (to traverse upwards)

D numpy int array

Vector of tree depths

tri_indices numpy int array

Vector of element indices (-1 if not leaf node)

See Also

traverse_aabbtree, initialize_quadtree.

Notes

This code is purposefully not optimized beyond asymptotics for simplicity in understanding its functionality and translating it to other programming languages beyond prototyping.

Examples:

from gpytoolbox import initialize_aabbtree
P = np.array([[0,0,0],[0.1,0,0],[0,0,0],[0.02,0.1,0],[1,1,0.9],[1,1,1]])
F = np.array([[0,1],[2,3],[4,5]],dtype=int)
C,W,CH,PAR,D,tri_ind = gpytoolbox.initialize_aabbtree(P,F)
Source code in src/gpytoolbox/initialize_aabbtree.py
  3
  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
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
def initialize_aabbtree(V,F=None,vmin=None,vmax=None):
    """Axis-Aligned Bounding Box hierarchy for efficient computation

    Computes an AABB tree by recursively dividing cells along the biggest dimension. May be called for any dimension and simplex size

    Parameters
    ----------
    V : numpy double array
        Matrix of mesh vertex coordinates
    F : numpy int array, optional (default None)
        Matrix of element indices into V (if None, treated as point cloud)

    Returns
    -------
    C : numpy double array 
        Matrix of box centers
    W : numpy double array 
        Matrix of box widths
    CH : numpy int array
        Matrix of child indeces (-1 if leaf node)
    PAR : numpy int array 
        Vector of immediate parent indeces (to traverse upwards)
    D : numpy int array
        Vector of tree depths
    tri_indices : numpy int array
        Vector of element indices (-1 if *not* leaf node)

    See Also
    --------
    traverse_aabbtree, initialize_quadtree.

    Notes
    -----
    This code is *purposefully* not optimized beyond asymptotics for simplicity in understanding its functionality and translating it to other programming languages beyond prototyping.

    Examples
    --------
    ```python
    from gpytoolbox import initialize_aabbtree
    P = np.array([[0,0,0],[0.1,0,0],[0,0,0],[0.02,0.1,0],[1,1,0.9],[1,1,1]])
    F = np.array([[0,1],[2,3],[4,5]],dtype=int)
    C,W,CH,PAR,D,tri_ind = gpytoolbox.initialize_aabbtree(P,F)
    ```
    """

    if (F is None):
        F = np.linspace(0,V.shape[0]-1,V.shape[0],dtype=int)[:,None]

    # We start with a bounding box
    dim = V.shape[1]
    simplex_size = F.shape[1]
    num_s = F.shape[0]
    if (vmin is None):
        vmin = np.amin(V,axis=0)
    if (vmax is None):
        vmax = np.amax(V,axis=0)
    C = (vmin + vmax)/2.0
    C = C[None,:]
    #print(C)

    # We need to compute this once, we'll need it for the subdivision:
    vmin_big = 10000*np.ones((num_s,dim))
    vmax_big = -10000*np.ones((num_s,dim))
    for j in range(dim):
        for i in range(simplex_size):
            vmax_big[:,j] = np.amax(np.hstack((V[F[:,i],j][:,None],vmax_big[:,j][:,None])),axis=1)[:]
            vmin_big[:,j] = np.amin(np.hstack((V[F[:,i],j][:,None],vmin_big[:,j][:,None])),axis=1)[:]


    W = np.reshape(vmax-vmin,(1,dim))
    CH = np.tile(np.array([[-1]],dtype=int),(1,2)) # for now it's leaf node
    D = np.array([1],dtype=int)
    PAR = np.array([-1],dtype=int) # supreme Neanderthal ancestral node
    tri_indices = np.array([-1]) # this will have the index in leaf nodes (at least temporarily)
    tri_index_list = [ list(range(F.shape[0])) ]
    # Now, we will loop
    box_ind = -1
    while True:
        box_ind = box_ind + 1
        if box_ind>=C.shape[0]:
            break
        is_child = (CH[box_ind,1]==-1)
        assert(is_child) # This check should be superfluous I think
        #tris_in_box = is_in_box(V,F,C[box_ind,:],W[box_ind,:])
        tris_in_box = tri_index_list[box_ind]
        # print(box_ind)
        # print(tris_in_box)
        # print(tri_indices)
        # print(tris_in_box)
        num_tris_in_box = len(tris_in_box)
        # print(num_tris_in_box)
        assert(num_tris_in_box>0) # There can't be a node with zero triangles...
        if (is_child and num_tris_in_box>=2):
            # Does this quad contain more than one triangle? Then subdivide it
            C,W,CH,PAR,D,tri_indices,tri_index_list = subdivide_box(box_ind,V,F,tris_in_box,C,W,CH,PAR,D,tri_indices,tri_index_list,vmin_big,vmax_big)
        if (is_child and num_tris_in_box==1):
            tri_indices[box_ind] = tris_in_box[0] # Check this??

    return C,W,CH,PAR,D,tri_indices