Neko  0.9.0
A portable framework for high-order spectral element flow simulations
aabb Module Reference

Axis Aligned Bounding Box (aabb) implementation in Fortran. More...

Data Types

interface  merge
 Merge two aabbs. More...
 
interface  intersection
 Intersect two aabbs. More...
 
type  aabb_t
 Axis Aligned Bounding Box (aabb) data structure. More...
 

Functions/Subroutines

type(aabb_t) function, public get_aabb (object, padding)
 Construct the aabb of a predefined object. More...
 
subroutine add_padding (this, padding)
 Add padding to the aabb. More...
 
type(aabb_t) function get_aabb_element (object)
 Get the aabb of an arbitrary element. More...
 
type(aabb_t) function get_aabb_mesh (object)
 Get the aabb of a mesh. More...
 
type(aabb_t) function get_aabb_tri_mesh (object)
 Get the aabb of a triangular mesh. More...
 
type(aabb_t) function get_aabb_tet_mesh (object)
 Get the aabb of a tetrahedral mesh. More...
 
subroutine aabb_init (this, lower_left_front, upper_right_back)
 Initialize the aabb. More...
 
pure real(kind=dp) function, dimension(3) aabb_get_min (this)
 Get the minimum point of the aabb. More...
 
pure real(kind=dp) function, dimension(3) aabb_get_max (this)
 Get the maximum point of the aabb. More...
 
pure real(kind=dp) function aabb_get_width (this)
 Get the width of the aabb. Also known as the x-axis length. More...
 
pure real(kind=dp) function aabb_get_depth (this)
 Get the depth of the aabb. Also known as the y-axis length. More...
 
pure real(kind=dp) function aabb_get_height (this)
 Get the height of the aabb. Also known as the z-axis length. More...
 
pure real(kind=dp) function aabb_get_diameter (this)
 Get the diameter length of the aabb. More...
 
pure real(kind=dp) function aabb_get_surface_area (this)
 Get the surface area of the aabb. More...
 
pure real(kind=dp) function, dimension(3) aabb_get_center (this)
 Get the center of the aabb. More...
 
pure real(kind=dp) function, dimension(3) aabb_get_diagonal (this)
 Get the diagonal of the aabb. More...
 
logical function aabb_overlaps (this, other)
 Check if two aabbs are overlapping. More...
 
logical function aabb_contains_other (this, other)
 Check if this aabb contains another aabb. More...
 
logical function aabb_contains_point (this, p)
 Check if this aabb contains a point. More...
 
type(aabb_t) function merge_aabb (box1, box2)
 Merge two aabbs. More...
 
type(aabb_t) function intersection_aabb (box1, box2)
 Get the intersection of two aabbs. More...
 
pure real(kind=dp) function calculate_surface_area (this)
 Calculate the surface area of the aabb. More...
 
pure logical function aabb_less (this, other)
 Less than comparison operator. More...
 
pure logical function aabb_greater (this, other)
 Greater than comparison operator. More...
 

Detailed Description

This is a Fortran implementation of an Axis Aligned Bounding Box (aabb) data structure. The aabb is a box that is aligned to the x, y and z axes. It is defined by two points, the lower left front corner and the upper right back corner. This is the base data structure for the aabb_Tree, which is used to accelerate a Signed Distance Function.

Function/Subroutine Documentation

◆ aabb_contains_other()

logical function aabb::aabb_contains_other ( class(aabb_t), intent(in)  this,
class(aabb_t), intent(in)  other 
)
private

Definition at line 447 of file aabb.f90.

◆ aabb_contains_point()

logical function aabb::aabb_contains_point ( class(aabb_t), intent(in)  this,
real(kind=dp), dimension(3), intent(in)  p 
)
private

Definition at line 462 of file aabb.f90.

◆ aabb_get_center()

pure real(kind=dp) function, dimension(3) aabb::aabb_get_center ( class(aabb_t), intent(in)  this)
private

Definition at line 410 of file aabb.f90.

◆ aabb_get_depth()

pure real(kind=dp) function aabb::aabb_get_depth ( class(aabb_t), intent(in)  this)
private

Definition at line 378 of file aabb.f90.

◆ aabb_get_diagonal()

pure real(kind=dp) function, dimension(3) aabb::aabb_get_diagonal ( class(aabb_t), intent(in)  this)
private

Definition at line 418 of file aabb.f90.

◆ aabb_get_diameter()

pure real(kind=dp) function aabb::aabb_get_diameter ( class(aabb_t), intent(in)  this)
private

Definition at line 394 of file aabb.f90.

◆ aabb_get_height()

pure real(kind=dp) function aabb::aabb_get_height ( class(aabb_t), intent(in)  this)
private

Definition at line 386 of file aabb.f90.

◆ aabb_get_max()

pure real(kind=dp) function, dimension(3) aabb::aabb_get_max ( class(aabb_t), intent(in)  this)
private

Definition at line 362 of file aabb.f90.

◆ aabb_get_min()

pure real(kind=dp) function, dimension(3) aabb::aabb_get_min ( class(aabb_t), intent(in)  this)
private

Definition at line 354 of file aabb.f90.

◆ aabb_get_surface_area()

pure real(kind=dp) function aabb::aabb_get_surface_area ( class(aabb_t), intent(in)  this)
private

Definition at line 402 of file aabb.f90.

◆ aabb_get_width()

pure real(kind=dp) function aabb::aabb_get_width ( class(aabb_t), intent(in)  this)
private

Definition at line 370 of file aabb.f90.

◆ aabb_greater()

pure logical function aabb::aabb_greater ( class(aabb_t), intent(in)  this,
class(aabb_t), intent(in)  other 
)
private

Definition at line 564 of file aabb.f90.

◆ aabb_init()

subroutine aabb::aabb_init ( class(aabb_t), intent(inout)  this,
real(kind=dp), dimension(3), intent(in)  lower_left_front,
real(kind=dp), dimension(3), intent(in)  upper_right_back 
)
private
Parameters
lower_left_frontThe lower left front corner of the aabb.
upper_right_backThe upper right back corner of the aabb.

Definition at line 334 of file aabb.f90.

◆ aabb_less()

pure logical function aabb::aabb_less ( class(aabb_t), intent(in)  this,
class(aabb_t), intent(in)  other 
)
private

Definition at line 527 of file aabb.f90.

◆ aabb_overlaps()

logical function aabb::aabb_overlaps ( class(aabb_t), intent(in)  this,
class(aabb_t), intent(in)  other 
)
private

Definition at line 430 of file aabb.f90.

◆ add_padding()

subroutine aabb::add_padding ( class(aabb_t), intent(inout)  this,
real(kind=dp), intent(in)  padding 
)

This function adds padding to the aabb. The padding is a multiple of the diameter of the aabb. This is used to avoid numerical issues when the object itself it axis aligned.

Parameters
[in]paddingThe padding of the aabb.

Definition at line 221 of file aabb.f90.

◆ calculate_surface_area()

pure real(kind=dp) function aabb::calculate_surface_area ( class(aabb_t), intent(in)  this)
private

Definition at line 511 of file aabb.f90.

◆ get_aabb()

type(aabb_t) function, public aabb::get_aabb ( class(*), intent(in)  object,
real(kind=dp), intent(in), optional  padding 
)

This function is used to get the aabb of a predefined object. Optionally, the user can define the padding of the aabb, which is a multiple of the diameter of the aabb. This is used to avoid numerical issues when the object itself it axis aligned.

Current support:

  • Triangle (tri_t)
  • Quadrilateral (quad_t)
  • Tetrahedron (tet_t)
  • Hexahedron (hex_t)
  • Mesh (mesh_t)
  • Triangular mesh (tri_mesh_t)
  • Tetrahedral mesh (tet_mesh_t)
Parameters
[in]objectThe object to get the aabb of.
[in]paddingThe padding of the aabb.
Returns
The aabb of the object.

Definition at line 180 of file aabb.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_aabb_element()

type(aabb_t) function aabb::get_aabb_element ( class(element_t), intent(in)  object)
private

This function calculates the aabb of an element. The aabb is defined by the lower left front corner and the upper right back corner. The aabb is calculated by finding the minimum and maximum x, y and z coordinate for all points in the arbitrary element type.

Parameters
objectThe arbitrary element to get the aabb of.
Returns
The aabb of the element.

Definition at line 241 of file aabb.f90.

Here is the caller graph for this function:

◆ get_aabb_mesh()

type(aabb_t) function aabb::get_aabb_mesh ( type(mesh_t), intent(in)  object)
private

This function calculates the aabb of a mesh. The aabb is defined by the lower left front corner and the upper right back corner. The aabb is calculated by merging the aabb of all elements in the mesh.

Parameters
objectThe mesh to get the aabb of.
Returns
The aabb of the mesh.

Definition at line 268 of file aabb.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_aabb_tet_mesh()

type(aabb_t) function aabb::get_aabb_tet_mesh ( type(tet_mesh_t), intent(in)  object)
private

This function calculates the aabb of a mesh. The aabb is defined by the lower left front corner and the upper right back corner. The aabb is calculated by merging the aabb of all elements in the mesh.

Parameters
objectThe tetrahedral mesh to get the aabb of.
Returns
The aabb of the mesh.

Definition at line 312 of file aabb.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_aabb_tri_mesh()

type(aabb_t) function aabb::get_aabb_tri_mesh ( type(tri_mesh_t), intent(in)  object)
private

This function calculates the aabb of a mesh. The aabb is defined by the lower left front corner and the upper right back corner. The aabb is calculated by merging the aabb of all elements in the mesh.

Parameters
objectThe triangular mesh to get the aabb of.
Returns
The aabb of the mesh.

Definition at line 290 of file aabb.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ intersection_aabb()

type(aabb_t) function aabb::intersection_aabb ( class(aabb_t), intent(in)  box1,
class(aabb_t), intent(in)  box2 
)
private

Definition at line 493 of file aabb.f90.

◆ merge_aabb()

type(aabb_t) function aabb::merge_aabb ( class(aabb_t), intent(in)  box1,
class(aabb_t), intent(in)  box2 
)
private

Definition at line 479 of file aabb.f90.