CellModules
HertzianPhysics::PluginHertzianPhysics< cell_t > Class Template Reference

Class for managing Delaunay triangulation and Hertzian physics. More...

#include <PluginHertzianPhysics.hpp>

Public Member Functions

 PluginHertzianPhysics ()
 Constructor. More...
 
void setCoherenceCoeff (double c)
 Sets the coherence coefficient. More...
 
template<typename world_t >
void beginUpdate (world_t *w)
 Begin update hook for MecaCell. More...
 
template<typename world_t >
void onAddCell (world_t *w)
 On add cell hook for MecaCell. More...
 
template<typename world_t >
void preDeleteDeadCellsUpdate (world_t *w)
 Pre-delete dead cells update hook for MecaCell. More...
 
template<typename world_t >
void preBehaviorUpdate (world_t *w)
 Pre-behavior update hook for MecaCell. More...
 
void setMinConstraint (MecaCell::Vec min)
 Sets the minimum constraint. More...
 
void setMaxConstraint (MecaCell::Vec max)
 Sets the maximum constraint. More...
 
PhysicsConstraints getConstraints ()
 Gets the physics constraints. More...
 
void addPlane (MecaCell::Vec _n=MecaCell::Vec(0, 1, 0), MecaCell::Vec _p=MecaCell::Vec(0, 0, 0), bool _sticky=false)
 Adds a plane constraint. More...
 
void addSymAxis (std::function< double(double, MecaCell::Vec, MecaCell::Vec)> _radiusFunction, MecaCell::Vec _axis=MecaCell::Vec(0., 1., 0.), MecaCell::Vec _origin=MecaCell::Vec(0., 0., 0.), double _length=100, double _maxEffectRadius=100, bool _fitShape=true, double _velocity_attraction=0.001)
 Adds a symmetrical axis constraint. More...
 
void addVelocity (MecaCell::Vec _v=MecaCell::Vec(0, 1, 0))
 Adds a velocity constraint. More...
 
void setCurvedBox (double radiusBoxPlane, double curvedDistToMaxSpeed, double maxSpeed)
 Sets the curved box parameters. More...
 

Private Member Functions

bool isDelaunayCoherent ()
 Checks the coherence of the Delaunay triangulation. More...
 
void moveDelaunay ()
 Moves the points in the triangulation to match their real positions. More...
 
template<typename world_t >
void computeDelaunay (world_t *w)
 Computes the Delaunay triangulation from the list of cells. More...
 
template<typename world_t >
void addCells (world_t *w)
 Adds new cells to the triangulation. More...
 
void removeCells ()
 Removes dead cells from the triangulation. More...
 
template<typename world_t >
void updateConnections (world_t *w)
 Updates the connections between cells. More...
 
template<typename world_t >
void customConstraint (world_t *w)
 Applies custom spatial constraints to the cells. More...
 
template<typename world_t >
double updateForces (world_t *w)
 Updates the forces acting on the cells. More...
 
void updateInteractionForces (cell_t *c1)
 Updates the interaction forces between cells. More...
 
template<typename world_t >
void updatePositions (world_t *w, double dt)
 Updates the positions of the cells. More...
 
template<typename world_t >
void resetForces (world_t *w)
 Resets the forces acting on the cells. More...
 
template<typename world_t >
void applyPhysics (world_t *w)
 Applies the physics with a coherent time step. More...
 

Private Attributes

Delaunay triangulation
 
double coherenceCoeff = 0.2
 
bool areForcesNegligible = false
 
bool areMovementsNegligible = false
 
MecaCell::Vec minConstraint
 
MecaCell::Vec maxConstraint
 
PhysicsConstraints constraints
 

Detailed Description

template<typename cell_t>
class HertzianPhysics::PluginHertzianPhysics< cell_t >

Class for managing Delaunay triangulation and Hertzian physics.

Template Parameters
cell_tType of the cell.

Definition at line 41 of file PluginHertzianPhysics.hpp.

Constructor & Destructor Documentation

◆ PluginHertzianPhysics()

template<typename cell_t >
HertzianPhysics::PluginHertzianPhysics< cell_t >::PluginHertzianPhysics ( )
inline

Constructor.

Initializes the custom constraints.

Definition at line 366 of file PluginHertzianPhysics.hpp.

Member Function Documentation

◆ addCells()

template<typename cell_t >
template<typename world_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::addCells ( world_t *  w)
inlineprivate

Adds new cells to the triangulation.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.

Definition at line 133 of file PluginHertzianPhysics.hpp.

◆ addPlane()

template<typename cell_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::addPlane ( MecaCell::Vec  _n = MecaCell::Vec(0, 1, 0),
MecaCell::Vec  _p = MecaCell::Vec(0, 0, 0),
bool  _sticky = false 
)
inline

Adds a plane constraint.

Parameters
_nNormal vector of the plane.
_pPoint on the plane.
_stickyIndicates if the plane is sticky.

Definition at line 464 of file PluginHertzianPhysics.hpp.

◆ addSymAxis()

template<typename cell_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::addSymAxis ( std::function< double(double, MecaCell::Vec, MecaCell::Vec)>  _radiusFunction,
MecaCell::Vec  _axis = MecaCell::Vec(0., 1., 0.),
MecaCell::Vec  _origin = MecaCell::Vec(0., 0., 0.),
double  _length = 100,
double  _maxEffectRadius = 100,
bool  _fitShape = true,
double  _velocity_attraction = 0.001 
)
inline

Adds a symmetrical axis constraint.

Parameters
_radiusFunctionRadius function.
_axisAxis vector.
_originOrigin point.
_lengthLength of the axis.
_maxEffectRadiusMaximum effect radius.
_fitShapeIndicates if the shape should fit.
_velocity_attractionAttraction velocity.

Definition at line 477 of file PluginHertzianPhysics.hpp.

◆ addVelocity()

template<typename cell_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::addVelocity ( MecaCell::Vec  _v = MecaCell::Vec(0, 1, 0))
inline

Adds a velocity constraint.

Parameters
_vVelocity vector.

Definition at line 486 of file PluginHertzianPhysics.hpp.

◆ applyPhysics()

template<typename cell_t >
template<typename world_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::applyPhysics ( world_t *  w)
inlineprivate

Applies the physics with a coherent time step.

Applies the physics with a coherent physics dt. Stops when movements or forces are too weak.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.

Definition at line 343 of file PluginHertzianPhysics.hpp.

◆ beginUpdate()

template<typename cell_t >
template<typename world_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::beginUpdate ( world_t *  w)
inline

Begin update hook for MecaCell.

Recomputes or readjusts the triangulation and regenerates the connections.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.

Definition at line 386 of file PluginHertzianPhysics.hpp.

◆ computeDelaunay()

template<typename cell_t >
template<typename world_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::computeDelaunay ( world_t *  w)
inlineprivate

Computes the Delaunay triangulation from the list of cells.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.

Definition at line 113 of file PluginHertzianPhysics.hpp.

◆ customConstraint()

template<typename cell_t >
template<typename world_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::customConstraint ( world_t *  w)
inlineprivate

Applies custom spatial constraints to the cells.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.

Definition at line 197 of file PluginHertzianPhysics.hpp.

◆ getConstraints()

template<typename cell_t >
PhysicsConstraints HertzianPhysics::PluginHertzianPhysics< cell_t >::getConstraints ( )
inline

Gets the physics constraints.

Returns
PhysicsConstraints object.

Definition at line 455 of file PluginHertzianPhysics.hpp.

◆ isDelaunayCoherent()

template<typename cell_t >
bool HertzianPhysics::PluginHertzianPhysics< cell_t >::isDelaunayCoherent ( )
inlineprivate

Checks the coherence of the Delaunay triangulation.

Returns
False if the triangulation needs to be recomputed.

Definition at line 58 of file PluginHertzianPhysics.hpp.

◆ moveDelaunay()

template<typename cell_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::moveDelaunay ( )
inlineprivate

Moves the points in the triangulation to match their real positions.

Definition at line 86 of file PluginHertzianPhysics.hpp.

◆ onAddCell()

template<typename cell_t >
template<typename world_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::onAddCell ( world_t *  w)
inline

On add cell hook for MecaCell.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.

Definition at line 398 of file PluginHertzianPhysics.hpp.

◆ preBehaviorUpdate()

template<typename cell_t >
template<typename world_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::preBehaviorUpdate ( world_t *  w)
inline

Pre-behavior update hook for MecaCell.

Applies the physics to each cell.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.

Definition at line 429 of file PluginHertzianPhysics.hpp.

◆ preDeleteDeadCellsUpdate()

template<typename cell_t >
template<typename world_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::preDeleteDeadCellsUpdate ( world_t *  w)
inline

Pre-delete dead cells update hook for MecaCell.

Removes the dead cells from the triangulation.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.

Definition at line 411 of file PluginHertzianPhysics.hpp.

◆ removeCells()

template<typename cell_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::removeCells ( )
inlineprivate

Removes dead cells from the triangulation.

Definition at line 145 of file PluginHertzianPhysics.hpp.

◆ resetForces()

template<typename cell_t >
template<typename world_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::resetForces ( world_t *  w)
inlineprivate

Resets the forces acting on the cells.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.

Definition at line 327 of file PluginHertzianPhysics.hpp.

◆ setCoherenceCoeff()

template<typename cell_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::setCoherenceCoeff ( double  c)
inline

Sets the coherence coefficient.

Parameters
cCoherence coefficient.

Definition at line 375 of file PluginHertzianPhysics.hpp.

◆ setCurvedBox()

template<typename cell_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::setCurvedBox ( double  radiusBoxPlane,
double  curvedDistToMaxSpeed,
double  maxSpeed 
)
inline

Sets the curved box parameters.

Parameters
radiusBoxPlaneRadius of the box plane.
curvedDistToMaxSpeedCurved distance to maximum speed.
maxSpeedMaximum speed.

Definition at line 495 of file PluginHertzianPhysics.hpp.

◆ setMaxConstraint()

template<typename cell_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::setMaxConstraint ( MecaCell::Vec  max)
inline

Sets the maximum constraint.

Parameters
maxMaximum constraint vector.

Definition at line 447 of file PluginHertzianPhysics.hpp.

◆ setMinConstraint()

template<typename cell_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::setMinConstraint ( MecaCell::Vec  min)
inline

Sets the minimum constraint.

Parameters
minMinimum constraint vector.

Definition at line 438 of file PluginHertzianPhysics.hpp.

◆ updateConnections()

template<typename cell_t >
template<typename world_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::updateConnections ( world_t *  w)
inlineprivate

Updates the connections between cells.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.

Definition at line 165 of file PluginHertzianPhysics.hpp.

◆ updateForces()

template<typename cell_t >
template<typename world_t >
double HertzianPhysics::PluginHertzianPhysics< cell_t >::updateForces ( world_t *  w)
inlineprivate

Updates the forces acting on the cells.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.
Returns
Minimum time step for the simulation.

Definition at line 239 of file PluginHertzianPhysics.hpp.

◆ updateInteractionForces()

template<typename cell_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::updateInteractionForces ( cell_t *  c1)
inlineprivate

Updates the interaction forces between cells.

Parameters
c1Pointer to the first cell.

Definition at line 265 of file PluginHertzianPhysics.hpp.

◆ updatePositions()

template<typename cell_t >
template<typename world_t >
void HertzianPhysics::PluginHertzianPhysics< cell_t >::updatePositions ( world_t *  w,
double  dt 
)
inlineprivate

Updates the positions of the cells.

Updates positions and checks if the movements are too weak to be considered.

Template Parameters
world_tType of the world.
Parameters
wPointer to the world.
dtTime step.

Definition at line 305 of file PluginHertzianPhysics.hpp.

Member Data Documentation

◆ areForcesNegligible

template<typename cell_t >
bool HertzianPhysics::PluginHertzianPhysics< cell_t >::areForcesNegligible = false
private

Flag indicating if forces are negligible

Definition at line 46 of file PluginHertzianPhysics.hpp.

◆ areMovementsNegligible

template<typename cell_t >
bool HertzianPhysics::PluginHertzianPhysics< cell_t >::areMovementsNegligible = false
private

Flag indicating if movements are negligible

Definition at line 47 of file PluginHertzianPhysics.hpp.

◆ coherenceCoeff

template<typename cell_t >
double HertzianPhysics::PluginHertzianPhysics< cell_t >::coherenceCoeff = 0.2
private

Coefficient used to determine the spatial coherence of each cell

Definition at line 45 of file PluginHertzianPhysics.hpp.

◆ constraints

template<typename cell_t >
PhysicsConstraints HertzianPhysics::PluginHertzianPhysics< cell_t >::constraints
private

Object managing physical constraints

Definition at line 52 of file PluginHertzianPhysics.hpp.

◆ maxConstraint

template<typename cell_t >
MecaCell::Vec HertzianPhysics::PluginHertzianPhysics< cell_t >::maxConstraint
private

Maximum constraint vector

Definition at line 50 of file PluginHertzianPhysics.hpp.

◆ minConstraint

template<typename cell_t >
MecaCell::Vec HertzianPhysics::PluginHertzianPhysics< cell_t >::minConstraint
private

Minimum constraint vector

Definition at line 49 of file PluginHertzianPhysics.hpp.

◆ triangulation

template<typename cell_t >
Delaunay HertzianPhysics::PluginHertzianPhysics< cell_t >::triangulation
private

Contains and computes the Delaunay triangulation of the cells

Definition at line 44 of file PluginHertzianPhysics.hpp.


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