KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
List of all members
Kratos::GeometricalObjectsBins Class Reference

A bins container for 3 dimensional GeometricalObject entities. More...

#include <geometrical_objects_bins.h>

Collaboration diagram for Kratos::GeometricalObjectsBins:

Public Member Functions

Life Cycle
template<typename TIteratorType >
 GeometricalObjectsBins (TIteratorType GeometricalObjectsBegin, TIteratorType GeometricalObjectsEnd, const double Tolerance=1e-12)
 The constructor with all geometries to be stored. Please note that all of them should be available at construction time and cannot be modified after. More...
 
template<typename TContainer >
 GeometricalObjectsBins (TContainer &rGeometricalObjectsVector, const double Tolerance=1e-12)
 The constructor with all geometries to be stored. Please note that all of them should be available at construction time and cannot be modified after. More...
 
virtual ~GeometricalObjectsBins ()
 Destructor. More...
 
Operations
CellTypeGetCell (const std::size_t I, const std::size_t J, const std::size_t K)
 Accessing cell_ijk giving the 3 indices. More...
 
BoundingBox< PointTypeGetCellBoundingBox (const std::size_t I, const std::size_t J, const std::size_t K)
 Calculating the cell_ijk bounding box. More...
 
void SearchInRadius (const PointType &rPoint, const double Radius, std::vector< ResultType > &rResults)
 This method takes a point and finds all of the objects in the given radius to it. More...
 
template<typename TPointIteratorType >
void SearchInRadius (TPointIteratorType itPointBegin, TPointIteratorType itPointEnd, const double Radius, std::vector< std::vector< ResultType >> &rResults)
 This method takes a point and finds all of the objects in the given radius to it (iterative version). More...
 
ResultType SearchNearestInRadius (const PointType &rPoint, const double Radius)
 This method takes a point and finds the nearest object to it in a given radius. More...
 
template<typename TPointIteratorType >
std::vector< ResultTypeSearchNearestInRadius (TPointIteratorType itPointBegin, TPointIteratorType itPointEnd, const double Radius)
 This method takes a point and finds the nearest object to it in a given radius (iterative version). More...
 
ResultType SearchNearest (const PointType &rPoint)
 This method takes a point and finds the nearest object to it. More...
 
template<typename TPointIteratorType >
std::vector< ResultTypeSearchNearest (TPointIteratorType itPointBegin, TPointIteratorType itPointEnd)
 This method takes a point and finds the nearest object to it (iterative version). More...
 
ResultType SearchIsInside (const PointType &rPoint)
 This method takes a point and search if it's inside an geometrical object of the domain. More...
 
template<typename TPointIteratorType >
std::vector< ResultTypeSearchIsInside (TPointIteratorType itPointBegin, TPointIteratorType itPointEnd)
 This method takes a point and search if it's inside an geometrical object of the domain (iterative version). More...
 
Access
const BoundingBox< PointType > & GetBoundingBox () const
 Getting the bins bounding box. More...
 
const array_1d< double, 3 > & GetCellSizes ()
 return an array with the x,y and z size of the cube More...
 
const array_1d< std::size_t, 3 > & GetNumberOfCells ()
 returns a 3D array having the number of cells in direction x, y and z More...
 
std::size_t GetTotalNumberOfCells ()
 The total number of cells in the container. More...
 
Input and output
virtual std::string Info () const
 Turn back information as a string. More...
 
virtual void PrintInfo (std::ostream &rOStream) const
 Print information about this object. More...
 
virtual void PrintData (std::ostream &rOStream) const
 Print object's data. More...
 

Protected Member Functions

Protected Life Cycle
 GeometricalObjectsBins ()=default
 Default constructor protected. More...
 
Protected Operations
bool PointIsInsideBoundingBox (const array_1d< double, 3 > &rCoords)
 This method checks if a point is inside any bounding box of the global bounding boxes. More...
 
bool PointIsInsideBoundingBoxWithTolerance (const array_1d< double, 3 > &rCoords, const double Tolerance)
 This method checks if a point is inside any bounding box of the global bounding boxes considering a certain tolerance. More...
 
void CalculateCellSize (const std::size_t NumberOfCells)
 Calculate the cell sizes to be as equilateral as possible and tries to approximate (roughly) the given number of cells. More...
 
template<typename TIteratorType >
void AddObjectsToCells (TIteratorType GeometricalObjectsBegin, TIteratorType GeometricalObjectsEnd)
 Adding objects to the cells that intersecting with it. More...
 

Protected Attributes

Protected Member Variables
BoundingBox< PointTypemBoundingBox
 
array_1d< std::size_t, DimensionmNumberOfCells
 The bounding box of the domain. More...
 
array_1d< double, 3 > mCellSizes
 The number of cells in each direction. More...
 
array_1d< double, 3 > mInverseOfCellSize
 The size of each cell in each direction. More...
 
std::vector< CellTypemCells
 The inverse of the size of each cell in each direction. More...
 
double mTolerance
 The cells of the domain. More...
 

Static Protected Attributes

Protected Static Member Variables
static constexpr unsigned int Dimension = 3
 

Type Definitions

using PointType = Point
 The point type definition. More...
 
using ObjectType = GeometricalObject
 The type of geometrical object to be stored in the bins. More...
 
using CellType = std::vector< GeometricalObject * >
 The type of geometrical object to be stored in the bins. More...
 
using ResultType = SpatialSearchResult< GeometricalObject >
 
 KRATOS_CLASS_POINTER_DEFINITION (GeometricalObjectsBins)
 Pointer definition of GeometricalObjectsBins. More...
 

Detailed Description

A bins container for 3 dimensional GeometricalObject entities.

It provides efficient search in radius and search nearest methods. All of the geometries should be given at construction time. After constructing the bins the geometries cannot be modified. In case of any modification, the bins should be reconstructed.

Author
Pooyan Dadvand

Member Typedef Documentation

◆ CellType

The type of geometrical object to be stored in the bins.

◆ ObjectType

The type of geometrical object to be stored in the bins.

◆ PointType

The point type definition.

◆ ResultType

Constructor & Destructor Documentation

◆ GeometricalObjectsBins() [1/3]

template<typename TIteratorType >
Kratos::GeometricalObjectsBins::GeometricalObjectsBins ( TIteratorType  GeometricalObjectsBegin,
TIteratorType  GeometricalObjectsEnd,
const double  Tolerance = 1e-12 
)
inline

The constructor with all geometries to be stored. Please note that all of them should be available at construction time and cannot be modified after.

Parameters
GeometricalObjectsBeginThe begin iterator of the geometries to be stored
GeometricalObjectsEndThe end iterator of the geometries to be stored
Template Parameters
TIteratorTypeThe type of the iterator

◆ GeometricalObjectsBins() [2/3]

template<typename TContainer >
Kratos::GeometricalObjectsBins::GeometricalObjectsBins ( TContainer &  rGeometricalObjectsVector,
const double  Tolerance = 1e-12 
)
inline

The constructor with all geometries to be stored. Please note that all of them should be available at construction time and cannot be modified after.

Parameters
rGeometricalObjectsVectorThe geometries to be stored
Template Parameters
TContainerThe container type

◆ ~GeometricalObjectsBins()

virtual Kratos::GeometricalObjectsBins::~GeometricalObjectsBins ( )
inlinevirtual

Destructor.

◆ GeometricalObjectsBins() [3/3]

Kratos::GeometricalObjectsBins::GeometricalObjectsBins ( )
protecteddefault

Default constructor protected.

Member Function Documentation

◆ AddObjectsToCells()

template<typename TIteratorType >
void Kratos::GeometricalObjectsBins::AddObjectsToCells ( TIteratorType  GeometricalObjectsBegin,
TIteratorType  GeometricalObjectsEnd 
)
inlineprotected

Adding objects to the cells that intersecting with it.

This method takes a geometrical object and adds it to the cells that intersecting with it.

Template Parameters
TIteratorTypeThe type of the iterator of the geometrical objects
Parameters
GeometricalObjectsBeginThe begining of the geometrical objects
GeometricalObjectsEndThe end of the geometrical objects

◆ CalculateCellSize()

void Kratos::GeometricalObjectsBins::CalculateCellSize ( const std::size_t  NumberOfCells)
protected

Calculate the cell sizes to be as equilateral as possible and tries to approximate (roughly) the given number of cells.

This method calculates the cell sizes to be as equilateral as possible and tries to approximate (roughly) the given number of cells

Parameters
NumberOfCellsThe number of cells to be calculated

◆ GetBoundingBox()

const BoundingBox<PointType>& Kratos::GeometricalObjectsBins::GetBoundingBox ( ) const
inline

Getting the bins bounding box.

Returns
The bounding box of the bins

◆ GetCell()

GeometricalObjectsBins::CellType & Kratos::GeometricalObjectsBins::GetCell ( const std::size_t  I,
const std::size_t  J,
const std::size_t  K 
)

Accessing cell_ijk giving the 3 indices.

Parameters
IThe index in x direction
JThe index in y direction
KThe index in z direction
Returns
The cell_ijk

◆ GetCellBoundingBox()

BoundingBox< GeometricalObjectsBins::PointType > Kratos::GeometricalObjectsBins::GetCellBoundingBox ( const std::size_t  I,
const std::size_t  J,
const std::size_t  K 
)

Calculating the cell_ijk bounding box.

Parameters
IThe index in x direction
JThe index in y direction
KThe index in z direction
Returns
The bounding box of the cell

◆ GetCellSizes()

const array_1d<double, 3>& Kratos::GeometricalObjectsBins::GetCellSizes ( )
inline

return an array with the x,y and z size of the cube

Returns
The size of the cube

◆ GetNumberOfCells()

const array_1d<std::size_t, 3>& Kratos::GeometricalObjectsBins::GetNumberOfCells ( )
inline

returns a 3D array having the number of cells in direction x, y and z

Returns
The number of cells in each direction

◆ GetTotalNumberOfCells()

std::size_t Kratos::GeometricalObjectsBins::GetTotalNumberOfCells ( )
inline

The total number of cells in the container.

Returns
The total number of cells

◆ Info()

virtual std::string Kratos::GeometricalObjectsBins::Info ( ) const
inlinevirtual

Turn back information as a string.

◆ KRATOS_CLASS_POINTER_DEFINITION()

Kratos::GeometricalObjectsBins::KRATOS_CLASS_POINTER_DEFINITION ( GeometricalObjectsBins  )

Pointer definition of GeometricalObjectsBins.

◆ PointIsInsideBoundingBox()

bool Kratos::GeometricalObjectsBins::PointIsInsideBoundingBox ( const array_1d< double, 3 > &  rCoords)
protected

This method checks if a point is inside any bounding box of the global bounding boxes.

Parameters
rCoordsThe coordinates of the point
Returns
True if the point is inside the bounding box

◆ PointIsInsideBoundingBoxWithTolerance()

bool Kratos::GeometricalObjectsBins::PointIsInsideBoundingBoxWithTolerance ( const array_1d< double, 3 > &  rCoords,
const double  Tolerance 
)
protected

This method checks if a point is inside any bounding box of the global bounding boxes considering a certain tolerance.

Parameters
rCoordsThe coordinates of the point
ToleranceThe tolerance
Returns
True if the point is inside the bounding box

◆ PrintData()

virtual void Kratos::GeometricalObjectsBins::PrintData ( std::ostream &  rOStream) const
inlinevirtual

Print object's data.

◆ PrintInfo()

virtual void Kratos::GeometricalObjectsBins::PrintInfo ( std::ostream &  rOStream) const
inlinevirtual

Print information about this object.

◆ SearchInRadius() [1/2]

void Kratos::GeometricalObjectsBins::SearchInRadius ( const PointType rPoint,
const double  Radius,
std::vector< ResultType > &  rResults 
)

This method takes a point and finds all of the objects in the given radius to it.

The result contains the object and also its distance to the point.

Parameters
rPointThe point to be checked
RadiusThe radius to be checked
rResultsThe results of the search

◆ SearchInRadius() [2/2]

template<typename TPointIteratorType >
void Kratos::GeometricalObjectsBins::SearchInRadius ( TPointIteratorType  itPointBegin,
TPointIteratorType  itPointEnd,
const double  Radius,
std::vector< std::vector< ResultType >> &  rResults 
)
inline

This method takes a point and finds all of the objects in the given radius to it (iterative version).

The result contains the object and also its distance to the point.

Parameters
itPointBeginThe first point iterator
itPointEndThe last point iterator
RadiusThe radius to be checked
rResultsThe results of the search
Template Parameters
TPointIteratorTypeThe type of the point iterator

◆ SearchIsInside() [1/2]

GeometricalObjectsBins::ResultType Kratos::GeometricalObjectsBins::SearchIsInside ( const PointType rPoint)

This method takes a point and search if it's inside an geometrical object of the domain.

If it is inside an object, it returns it, and search distance is set to zero. If there is no object, the result will be set to not found. Result contains a flag is the object has been found or not. This method is a simplified and faster method of SearchNearest.

Parameters
rPointThe point to be checked
Returns
ResultType The result of the search

◆ SearchIsInside() [2/2]

template<typename TPointIteratorType >
std::vector<ResultType> Kratos::GeometricalObjectsBins::SearchIsInside ( TPointIteratorType  itPointBegin,
TPointIteratorType  itPointEnd 
)
inline

This method takes a point and search if it's inside an geometrical object of the domain (iterative version).

If it is inside an object, it returns it, and search distance is set to zero. If there is no object, the result will be set to not found. Result contains a flag is the object has been found or not. This method is a simplified and faster method of SearchNearest.

Parameters
itPointBeginThe first point iterator
itPointEndThe last point iterator
Returns
std::vector<ResultType> The result of the search
Template Parameters
TPointIteratorTypeThe type of the point iterator

◆ SearchNearest() [1/2]

GeometricalObjectsBins::ResultType Kratos::GeometricalObjectsBins::SearchNearest ( const PointType rPoint)

This method takes a point and finds the nearest object to it.

If there are more than one object in the same minimum distance only one is returned Result contains a flag is the object has been found or not.

Parameters
rPointThe point to be checked
Returns
ResultType The result of the search

◆ SearchNearest() [2/2]

template<typename TPointIteratorType >
std::vector<ResultType> Kratos::GeometricalObjectsBins::SearchNearest ( TPointIteratorType  itPointBegin,
TPointIteratorType  itPointEnd 
)
inline

This method takes a point and finds the nearest object to it (iterative version).

If there are more than one object in the same minimum distance only one is returned Result contains a flag is the object has been found or not.

Parameters
itPointBeginThe first point iterator
itPointEndThe last point iterator
Returns
std::vector<ResultType> The result of the search
Template Parameters
TPointIteratorTypeThe type of the point iterator

◆ SearchNearestInRadius() [1/2]

GeometricalObjectsBins::ResultType Kratos::GeometricalObjectsBins::SearchNearestInRadius ( const PointType rPoint,
const double  Radius 
)

This method takes a point and finds the nearest object to it in a given radius.

If there are more than one object in the same minimum distance only one is returned If there are no objects in that radius the result will be set to not found. Result contains a flag is the object has been found or not.

Parameters
rPointThe point to be checked
RadiusThe radius to be checked
Returns
ResultType The result of the search

◆ SearchNearestInRadius() [2/2]

template<typename TPointIteratorType >
std::vector<ResultType> Kratos::GeometricalObjectsBins::SearchNearestInRadius ( TPointIteratorType  itPointBegin,
TPointIteratorType  itPointEnd,
const double  Radius 
)
inline

This method takes a point and finds the nearest object to it in a given radius (iterative version).

If there are more than one object in the same minimum distance only one is returned If there are no objects in that radius the result will be set to not found. Result contains a flag is the object has been found or not.

Parameters
itPointBeginThe first point iterator
itPointEndThe last point iterator
RadiusThe radius to be checked
Returns
std::vector<ResultType> The result of the search
Template Parameters
TPointIteratorTypeThe type of the point iterator

Member Data Documentation

◆ Dimension

constexpr unsigned int Kratos::GeometricalObjectsBins::Dimension = 3
staticconstexprprotected

◆ mBoundingBox

BoundingBox<PointType> Kratos::GeometricalObjectsBins::mBoundingBox
protected

◆ mCells

std::vector<CellType> Kratos::GeometricalObjectsBins::mCells
protected

The inverse of the size of each cell in each direction.

◆ mCellSizes

array_1d<double, 3> Kratos::GeometricalObjectsBins::mCellSizes
protected

The number of cells in each direction.

◆ mInverseOfCellSize

array_1d<double, 3> Kratos::GeometricalObjectsBins::mInverseOfCellSize
protected

The size of each cell in each direction.

◆ mNumberOfCells

array_1d<std::size_t, Dimension> Kratos::GeometricalObjectsBins::mNumberOfCells
protected

The bounding box of the domain.

◆ mTolerance

double Kratos::GeometricalObjectsBins::mTolerance
protected

The cells of the domain.


The documentation for this class was generated from the following files: