89     integer :: object_index = -1
 
  104     procedure, pass(this), 
public :: get_object_index => &
 
  106     procedure, pass(this), 
public :: get_parent_index => &
 
  108     procedure, pass(this), 
public :: get_left_index => &
 
  110     procedure, pass(this), 
public :: get_right_index => &
 
  121     generic :: 
operator(.lt.) => less
 
  122     generic :: 
operator(.gt.) => greater
 
 
  134     integer :: allocated_node_count = 0
 
  136     integer :: node_capacity = 0
 
  137     integer :: growth_size = 1
 
  145     procedure, pass(this), 
public :: insert_object => &
 
  147     generic :: build => build_generic
 
  152     procedure, pass(this), 
public :: get_root_index => &
 
  154     procedure, pass(this), 
public :: get_parent_index => &
 
  156     procedure, pass(this), 
public :: get_left_index => &
 
  158     procedure, pass(this), 
public :: get_right_index => &
 
  162     procedure, pass(this), 
public :: get_root_node => &
 
  164     procedure, pass(this), 
public :: get_parent_node => &
 
  166     procedure, pass(this), 
public :: get_left_node => &
 
  168     procedure, pass(this), 
public :: get_right_node => &
 
  173     procedure, pass(this), 
public :: query_overlaps => &
 
 
  202    this%object_index = -1
 
 
  223    integer :: object_index
 
  225    object_index = this%object_index
 
 
  231    integer :: parent_index
 
  233    parent_index = this%parent_node_index
 
 
  239    integer :: left_index
 
  241    left_index = this%left_node_index
 
 
  247    integer :: right_index
 
  249    right_index = this%right_node_index
 
 
  255    real(kind=
dp), 
dimension(3), 
intent(in) :: p
 
  256    real(kind=
dp) :: distance
 
  258    distance = 0.5_rp * this%aabb%get_diameter() &
 
  259         - norm2(this%aabb%get_center() - p)
 
 
  279    if (this%is_leaf()) 
then 
  283       & this%object_index .gt. 0
 
  288       & this%object_index .eq. -1
 
 
  302    res = this%aabb .lt. other%aabb
 
 
  312    res = this%aabb .gt. other%aabb
 
 
  323    integer, 
intent(in) :: initial_capacity
 
  324    integer :: i, nonzero_capacity
 
  326    if (initial_capacity < 1) 
then 
  329       nonzero_capacity = initial_capacity
 
  333    this%allocated_node_count = 0
 
  334    this%next_free_node_index = 1
 
  335    this%node_capacity = nonzero_capacity
 
  336    this%growth_size = nonzero_capacity
 
  338    if (
allocated(this%nodes)) 
deallocate(this%nodes)
 
  339    allocate(this%nodes(nonzero_capacity))
 
  341    do i = 1, nonzero_capacity
 
  342       this%nodes(i)%next_node_index = i + 1
 
 
  350    type(
aabb_t), 
intent(in) :: objects(:)
 
  351    real(kind=
dp), 
optional, 
intent(in) :: padding
 
  353    integer :: i_obj, i_node, i
 
  356    integer :: start_layer, end_layer
 
  358    type(
aabb_t), 
allocatable :: box_list(:)
 
  359    integer, 
dimension(:), 
allocatable :: sorted_indices
 
  361    real(kind=
dp) :: aabb_padding
 
  363    if (
allocated(box_list)) 
deallocate(box_list)
 
  364    allocate(box_list(
size(objects)))
 
  366    call this%init(
size(objects) * 2)
 
  367    if (
size(objects) .eq. 0) 
then 
  376    if (
present(padding)) 
then 
  377       aabb_padding = padding
 
  379       aabb_padding = 0.0_dp
 
  382    do i_obj = 1, 
size(objects)
 
  383       box_list(i_obj) = 
get_aabb(objects(i_obj), aabb_padding)
 
  385    call sort(box_list, sorted_indices)
 
  387    do i = 1, 
size(sorted_indices)
 
  388       i_obj = sorted_indices(i)
 
  389       i_node = this%allocate_node()
 
  390       this%nodes(i_node)%aabb = box_list(i_obj)
 
  391       this%nodes(i_node)%object_index = i_obj
 
  396    end_layer = 
size(objects)
 
  398    do while (.not. done)
 
  401       do i = start_layer, end_layer - 1, 2
 
  402          i_node = this%allocate_node()
 
  404          this%nodes(i_node)%aabb = 
merge(this%nodes(i)%aabb, &
 
  405               this%nodes(i + 1)%aabb)
 
  407          this%nodes(i_node)%left_node_index = i
 
  408          this%nodes(i_node)%right_node_index = i + 1
 
  410          this%nodes(i)%parent_node_index = i_node
 
  411          this%nodes(i + 1)%parent_node_index = i_node
 
  416       if (mod(end_layer - start_layer, 2) .eq. 0) 
then 
  417          i_node = this%allocate_node()
 
  418          this%nodes(i_node)%aabb = this%nodes(end_layer)%aabb
 
  419          this%nodes(i_node)%left_node_index = end_layer
 
  422          this%nodes(end_layer)%parent_node_index = i_node
 
  426       start_layer = end_layer + 1
 
  427       end_layer = this%allocated_node_count
 
  430       done = start_layer .eq. end_layer
 
  434    this%root_node_index = this%allocated_node_count
 
  436    if (this%get_size() .ne. 
size(objects)) 
then 
  437       print *, 
"this%get_size() = ", this%get_size()
 
  438       print *, 
"size(objects) = ", 
size(objects)
 
 
  449    class(*), 
target, 
intent(in) :: objects(:)
 
  450    real(kind=
dp), 
optional, 
intent(in) :: padding
 
  452    integer :: i_obj, i_node, i
 
  455    integer :: start_layer, end_layer
 
  457    type(
aabb_t), 
allocatable :: box_list(:)
 
  458    integer, 
dimension(:), 
allocatable :: sorted_indices
 
  460    real(kind=
dp) :: aabb_padding
 
  462    if (
allocated(box_list)) 
deallocate(box_list)
 
  463    allocate(box_list(
size(objects)))
 
  465    call this%init(
size(objects) * 2)
 
  466    if (
size(objects) .eq. 0) 
then 
  475    if (
present(padding)) 
then 
  476       aabb_padding = padding
 
  478       aabb_padding = 0.0_dp
 
  481    do i_obj = 1, 
size(objects)
 
  482       box_list(i_obj) = 
get_aabb(objects(i_obj), aabb_padding)
 
  484    call sort(box_list, sorted_indices)
 
  486    do i = 1, 
size(sorted_indices)
 
  487       i_obj = sorted_indices(i)
 
  488       i_node = this%allocate_node()
 
  489       this%nodes(i_node)%aabb = box_list(i_obj)
 
  490       this%nodes(i_node)%object_index = i_obj
 
  495    end_layer = 
size(objects)
 
  497    do while (.not. done)
 
  500       do i = start_layer, end_layer - 1, 2
 
  501          i_node = this%allocate_node()
 
  503          this%nodes(i_node)%aabb = 
merge(this%nodes(i)%aabb, &
 
  504               this%nodes(i + 1)%aabb)
 
  506          this%nodes(i_node)%left_node_index = i
 
  507          this%nodes(i_node)%right_node_index = i + 1
 
  509          this%nodes(i)%parent_node_index = i_node
 
  510          this%nodes(i + 1)%parent_node_index = i_node
 
  515       if (mod(end_layer - start_layer, 2) .eq. 0) 
then 
  516          i_node = this%allocate_node()
 
  517          this%nodes(i_node)%aabb = this%nodes(end_layer)%aabb
 
  518          this%nodes(i_node)%left_node_index = end_layer
 
  521          this%nodes(end_layer)%parent_node_index = i_node
 
  525       start_layer = end_layer + 1
 
  526       end_layer = this%allocated_node_count
 
  529       done = start_layer .eq. end_layer
 
  533    this%root_node_index = this%allocated_node_count
 
  535    if (this%get_size() .ne. 
size(objects)) 
then 
  536       print *, 
"this%get_size() = ", this%get_size()
 
  537       print *, 
"size(objects) = ", 
size(objects)
 
 
  544  subroutine sort(array, indices)
 
  545    type(
aabb_t), 
dimension(:), 
intent(in) :: array
 
  546    integer, 
intent(inout), 
dimension(:), 
allocatable :: indices
 
  547    logical, 
dimension(:), 
allocatable :: visited
 
  552    allocate(indices(
size(array)))
 
  553    allocate(visited(
size(array)))
 
  557    do i = 1, 
size(array)
 
  559       do imin = 1, 
size(array)
 
  560          if (.not. visited(imin) .and. minidx .eq. -1) minidx = imin
 
  561          if (minidx .gt. -1) 
then 
  562             if (visited(imin) .and. array(imin) .lt. array(minidx)) minidx = imin
 
  567       visited(minidx) = .true.
 
 
  583    call simple_stack%init(this%allocated_node_count)
 
  585    tmp = this%get_root_index()
 
  587       call simple_stack%push(tmp)
 
  590    do while (.not. simple_stack%is_empty())
 
  591       idx = simple_stack%pop()
 
  594       if (this%nodes(idx)%is_leaf()) 
then 
  597          tmp = this%get_left_index(idx)
 
  598          call simple_stack%push(tmp)
 
  599          tmp = this%get_right_index(idx)
 
  600          call simple_stack%push(tmp)
 
 
  612    integer :: root_index
 
  614    root_index = this%root_node_index
 
 
  621    integer, 
intent(in) :: node_index
 
  622    integer :: parent_index
 
  624    parent_index = this%nodes(node_index)%parent_node_index
 
 
  631    integer, 
intent(in) :: node_index
 
  632    integer :: left_index
 
  634    left_index = this%nodes(node_index)%left_node_index
 
 
  641    integer, 
intent(in) :: node_index
 
  642    integer :: right_index
 
  644    right_index = this%nodes(node_index)%right_node_index
 
 
  653    integer, 
intent(in) :: node_index
 
  656    node = this%nodes(node_index)
 
 
  664    root_node = this%nodes(this%root_node_index)
 
 
  671    integer, 
intent(in) :: node_index
 
  674    parent_node = this%nodes(this%nodes(node_index)%parent_node_index)
 
 
  680    integer, 
intent(in) :: node_index
 
  683    left_node = this%nodes(this%nodes(node_index)%left_node_index)
 
 
  690    integer, 
intent(in) :: node_index
 
  693    right_node = this%nodes(this%nodes(node_index)%right_node_index)
 
 
  698    integer, 
intent(in) :: node_index
 
  701    out_box = this%nodes(node_index)%aabb
 
 
  709    class(*), 
intent(in) :: object
 
  710    integer, 
intent(in) :: object_index
 
  712    integer :: node_index
 
  714    node_index = this%allocate_node()
 
  715    this%nodes(node_index)%aabb = 
get_aabb(object)
 
  716    this%nodes(node_index)%object_index = object_index
 
  718    call this%insert_leaf(node_index)
 
 
  724    class(*), 
intent(in) :: object
 
  725    integer, 
intent(in) :: object_index
 
  729    type(
aabb_t) :: object_box
 
  731    integer :: root_index, left_index, right_index
 
  732    integer :: node_index, tmp_index
 
  735    root_index = this%get_root_index()
 
  737    call simple_stack%init()
 
  738    call simple_stack%push(root_index)
 
  740    do while (.not. simple_stack%is_empty())
 
  741       node_index = simple_stack%pop()
 
  745       if (this%nodes(node_index)%aabb%overlaps(object_box)) 
then 
  746          if (this%nodes(node_index)%is_leaf()) 
then 
  747             if (this%nodes(node_index)%object_index .ne. object_index) 
then 
  748                tmp_index = this%nodes(node_index)%object_index
 
  749                call overlaps%push(tmp_index)
 
  752             left_index = this%get_left_index(node_index)
 
  754                call simple_stack%push(left_index)
 
  756             right_index = this%get_right_index(node_index)
 
  758                call simple_stack%push(right_index)
 
  763    call simple_stack%free()
 
 
  772    integer :: node_index
 
  775       call this%resize_node_pool(this%node_capacity + this%growth_size)
 
  778    node_index = this%next_free_node_index
 
  780    associate(new_node => this%nodes(node_index))
 
  781      this%next_free_node_index = new_node%next_node_index
 
  787      this%next_free_node_index = new_node%next_node_index
 
  788      this%allocated_node_count = this%allocated_node_count + 1
 
 
  796    integer, 
intent(in) :: node_index
 
  798    this%nodes(node_index)%next_node_index = this%next_free_node_index
 
  799    this%next_free_node_index = node_index
 
  800    this%allocated_node_count = this%allocated_node_count - 1
 
 
  806    integer, 
intent(in) :: leaf_node_index
 
  808    integer :: tree_node_index
 
  810    real(kind=rp) :: cost_left
 
  811    real(kind=rp) :: cost_right
 
  818    type(aabb_t) :: combined_aabb
 
  819    real(kind=rp) :: new_parent_node_cost
 
  820    real(kind=rp) :: minimum_push_down_cost
 
  821    type(aabb_t) :: new_left_aabb
 
  822    type(aabb_t) :: new_right_aabb
 
  823    integer :: leaf_sibling_index
 
  825    integer :: old_parent_index
 
  826    integer :: new_parent_index
 
  831    leaf_node = this%nodes(leaf_node_index)
 
  835       this%root_node_index = leaf_node_index
 
  845    tree_node_index = this%root_node_index
 
  846    tree_node = this%get_node(tree_node_index)
 
  847    do while (.not. tree_node%is_leaf())
 
  851       left_node = this%get_left_node(tree_node_index)
 
  852       right_node = this%get_right_node(tree_node_index)
 
  856       combined_aabb = merge(tree_node%aabb, leaf_node%get_aabb())
 
  858       new_parent_node_cost = 2.0_rp * combined_aabb%get_surface_area()
 
  859       minimum_push_down_cost = 2.0_rp * ( &
 
  860       & combined_aabb%get_surface_area() &
 
  861       & - tree_node%aabb%get_surface_area()&
 
  866       if (left_node%is_leaf()) 
then 
  867          new_left_aabb = merge(leaf_node%aabb, left_node%get_aabb())
 
  868          cost_left = new_left_aabb%get_surface_area() + minimum_push_down_cost
 
  870          new_left_aabb = merge(leaf_node%aabb, left_node%get_aabb())
 
  872          & new_left_aabb%get_surface_area() &
 
  873          & - left_node%aabb%get_surface_area()&
 
  874          & ) + minimum_push_down_cost
 
  877       if (right_node%is_leaf()) 
then 
  878          new_right_aabb = merge(leaf_node%aabb, right_node%aabb)
 
  879          cost_right = new_right_aabb%get_surface_area() + &
 
  880               minimum_push_down_cost
 
  882          new_right_aabb = merge(leaf_node%aabb, right_node%aabb)
 
  884          & new_right_aabb%get_surface_area() &
 
  885          & - right_node%aabb%get_surface_area() &
 
  886          & ) + minimum_push_down_cost
 
  892       if (new_parent_node_cost < cost_left .and. &
 
  893            new_parent_node_cost < cost_right) 
then 
  898       if (cost_left .lt. cost_right) 
then 
  899          tree_node_index = tree_node%get_left_index()
 
  901          tree_node_index = tree_node%get_right_index()
 
  906       tree_node = this%get_node(tree_node_index)
 
  911    leaf_sibling_index = tree_node_index
 
  912    leaf_sibling = this%nodes(leaf_sibling_index)
 
  913    old_parent_index = this%get_parent_index(leaf_sibling_index)
 
  914    new_parent_index = this%allocate_node()
 
  915    new_parent = this%nodes(new_parent_index)
 
  916    new_parent%parent_node_index = old_parent_index
 
  917    new_parent%aabb = merge(leaf_node%aabb, leaf_sibling%aabb)
 
  919    if (leaf_node .lt. leaf_sibling) 
then 
  920       new_parent%left_node_index = leaf_node_index
 
  921       new_parent%right_node_index = leaf_sibling_index
 
  923       new_parent%left_node_index = leaf_sibling_index
 
  924       new_parent%right_node_index = leaf_node_index
 
  927    leaf_node%parent_node_index = new_parent_index
 
  928    leaf_sibling%parent_node_index = new_parent_index
 
  932       this%root_node_index = new_parent_index
 
  936       old_parent = this%nodes(old_parent_index)
 
  937       if (old_parent%left_node_index .eq. leaf_sibling_index) 
then 
  938          old_parent%left_node_index = new_parent_index
 
  940          old_parent%right_node_index = new_parent_index
 
  942       this%nodes(old_parent_index) = old_parent
 
  945    this%nodes(leaf_node_index) = leaf_node
 
  946    this%nodes(leaf_sibling_index) = leaf_sibling
 
  947    this%nodes(new_parent_index) = new_parent
 
  950    tree_node_index = leaf_node%parent_node_index
 
  952    call this%fix_upwards_tree(tree_node_index)
 
 
  961    type(stack_i4_t) :: simple_stack
 
  962    integer :: current_index
 
  963    integer :: root_index, left_index, right_index
 
  970    root_index = this%get_root_index()
 
  972    call simple_stack%init(this%node_capacity)
 
  973    call simple_stack%push(root_index)
 
  975    do while (.not. simple_stack%is_empty())
 
  976       current_index = simple_stack%pop()
 
  979       valid = valid .and. this%nodes(current_index)%is_valid()
 
  981       if (.not. this%nodes(current_index)%is_leaf()) 
then 
  982          left_index = this%get_left_index(current_index)
 
  983          right_index = this%get_right_index(current_index)
 
  985          call simple_stack%push(left_index)
 
  986          call simple_stack%push(right_index)
 
 
  996    integer, 
intent(in) :: tree_start_index
 
 1000    integer :: tree_node_index
 
 1002    tree_node_index = tree_start_index
 
 1004       left_node = this%get_left_node(tree_node_index)
 
 1005       right_node = this%get_right_node(tree_node_index)
 
 1007       this%nodes(tree_node_index)%aabb = merge(left_node%aabb, right_node%aabb)
 
 1009       tree_node_index = this%get_parent_index(tree_node_index)
 
 
 1016    type(stack_i4_t) :: simple_stack
 
 1018    integer :: current_index
 
 1019    integer :: root_index, left_index, right_index
 
 1021    root_index = this%get_root_index()
 
 1022    call simple_stack%init(this%node_capacity)
 
 1023    call simple_stack%push(root_index)
 
 1025    do while (.not. simple_stack%is_empty())
 
 1026       current_index = simple_stack%pop()
 
 1029       left_index = this%get_left_index(current_index)
 
 1030       right_index = this%get_right_index(current_index)
 
 1032       call simple_stack%push(this%nodes(current_index)%left_node_index)
 
 1033       call simple_stack%push(this%nodes(current_index)%right_node_index)
 
 1035       write(*, *) 
"i = ", current_index
 
 1036       write(*, *) 
"  Parent  : ", this%get_parent_index(current_index)
 
 1037       write(*, *) 
"  Children: ", this%get_left_index(current_index), &
 
 1038            this%get_right_index(current_index)
 
 1040       write(*, *) 
"  object_index = ", this%nodes(current_index)%object_index
 
 
 1048    integer, 
intent(in) :: new_capacity
 
 1050    type(
aabb_node_t), 
dimension(:), 
allocatable :: temp
 
 1053    allocate(temp(new_capacity))
 
 1054    temp(:this%node_capacity) = this%nodes(:this%node_capacity)
 
 1056    do i = this%allocated_node_count, new_capacity
 
 1057       temp(i)%next_node_index = i + 1
 
 1061    call move_alloc(temp, this%nodes)
 
 1063    this%node_capacity = new_capacity
 
 1064    this%next_free_node_index = this%allocated_node_count + 1
 
 
Axis Aligned Bounding Box (aabb) Tree data structure.
 
pure type(aabb_t) function aabb_tree_get_aabb(this, node_index)
 
subroutine aabb_tree_fix_upwards_tree(this, tree_start_index)
Fixes the tree upwards.
 
pure logical function aabb_node_is_valid(this)
Returns true if the node is a valid node.
 
subroutine aabb_node_init(this)
Initializes the AABB node.
 
pure type(aabb_node_t) function aabb_tree_get_right_node(this, node_index)
Returns the right node of the node at the given index.
 
pure type(aabb_node_t) function aabb_tree_get_parent_node(this, node_index)
Returns the parent node of the node at the given index.
 
subroutine aabb_tree_build_tree_aabb(this, objects, padding)
Builds the tree.
 
pure type(aabb_node_t) function aabb_tree_get_root_node(this)
Returns the root node of the tree.
 
subroutine aabb_tree_insert_object(this, object, object_index)
Inserts an object into the tree.
 
pure type(aabb_node_t) function aabb_tree_get_node(this, node_index)
Returns the node at the given index.
 
subroutine aabb_tree_print(this)
Prints the tree.
 
pure logical function aabb_node_greater(this, other)
Returns true if the node is greater than the other node.
 
integer function aabb_tree_allocate_node(this)
Allocates a new node in the tree.
 
pure type(aabb_t) function aabb_node_get_aabb(this)
Returns the Axis Aligned Bounding Box (aabb) of the node.
 
pure integer function aabb_tree_get_root_index(this)
Returns the index of the root node.
 
pure integer function aabb_node_get_right_index(this)
Returns the right index of the node.
 
subroutine aabb_tree_deallocate_node(this, node_index)
Deallocates a node in the tree.
 
pure integer function aabb_node_get_left_index(this)
Returns the left index of the node.
 
pure logical function aabb_node_is_leaf(this)
Returns true if the node is a leaf node.
 
real(kind=dp) function aabb_node_min_distance(this, p)
Get the minimum possible distance from the aabb to a point.
 
subroutine aabb_tree_init(this, initial_capacity)
Initializes the AABB tree.
 
subroutine aabb_tree_insert_leaf(this, leaf_node_index)
Inserts a leaf into the tree.
 
subroutine aabb_tree_query_overlaps(this, object, object_index, overlaps)
Queries the tree for overlapping objects.
 
subroutine sort(array, indices)
Return a list of sorted indices of the aabb nodes.
 
pure logical function aabb_node_less(this, other)
Returns true if the node is less than the other node.
 
pure integer function aabb_tree_get_parent_index(this, node_index)
Returns the index of the parent node of the node at the given index.
 
pure integer function aabb_node_get_parent_index(this)
Returns the parent index of the node.
 
pure integer function aabb_tree_get_left_index(this, node_index)
Returns the index of the left node of the node at the given index.
 
subroutine aabb_tree_build_tree(this, objects, padding)
Builds the tree.
 
integer, parameter, public aabb_null_node
 
pure integer function aabb_node_get_object_index(this)
Returns the object index of the node.
 
logical function aabb_tree_valid_tree(this)
Validates the tree.
 
pure type(aabb_node_t) function aabb_tree_get_left_node(this, node_index)
Returns the left node of the node at the given index.
 
integer function aabb_tree_get_size(this)
Returns the size of the tree, in number of leaves.
 
subroutine aabb_tree_resize_node_pool(this, new_capacity)
Resizes the node pool.
 
pure integer function aabb_tree_get_right_index(this, node_index)
Returns the index of the right node of the node at the given index.
 
Axis Aligned Bounding Box (aabb) implementation in Fortran.
 
type(aabb_t) function, public get_aabb(object, padding)
Construct the aabb of a predefined object.
 
integer, parameter, public dp
 
integer, parameter, public rp
Global precision used in computations.
 
Implements a dynamic stack ADT.
 
Defines a triangular element.
 
Axis Aligned Bounding Box (aabb) data structure.
 
Node type for the Axis Aligned Bounding Box (aabb) Tree.
 
Axis Aligned Bounding Box (aabb) Tree.