BaseInteractable Class Referenceabstract

Defines the basic properties that a interactable object can have. More...

#include <BaseInteractable.h>

+ Inheritance diagram for BaseInteractable:

Public Member Functions

 BaseInteractable ()
 Default BaseInteractable constructor. More...
 
 BaseInteractable (const BaseInteractable &p)
 Copy constructor. More...
 
 ~BaseInteractable () override
 Destructor, it simply destructs the BaseInteractable and all the objects it contains. More...
 
void read (std::istream &is) override
 Reads a BaseInteractable from an input stream. More...
 
void write (std::ostream &os) const override
 Write a BaseInteractable to an output stream. More...
 
unsigned int getIndSpecies () const
 Returns the index of the species associated with the interactable object. More...
 
virtual void setIndSpecies (unsigned int indSpecies)
 Sets the index of the Species of this BaseInteractable. More...
 
const ParticleSpeciesgetSpecies () const
 Returns a pointer to the species of this BaseInteractable. More...
 
void setSpecies (const ParticleSpecies *species)
 Sets the species of this BaseInteractable. More...
 
const Vec3DgetForce () const
 Returns the force on this BaseInteractable. More...
 
const Vec3DgetTorque () const
 Returns the torque on this BaseInteractable. More...
 
void setForce (const Vec3D &force)
 Sets the force on this BaseInteractable. More...
 
void setTorque (const Vec3D &torque)
 Sets the torque on this BaseInteractable. More...
 
void addForce (const Vec3D &addForce)
 Adds an amount to the force on this BaseInteractable. More...
 
void addTorque (const Vec3D &addTorque)
 Adds an amount to the torque on this BaseInteractable. More...
 
virtual void resetForceTorque (int numberOfOMPthreads)
 
void sumForceTorqueOMP ()
 
const Vec3DgetPosition () const
 Returns the position of this BaseInteractable. More...
 
const QuaterniongetOrientation () const
 Returns the orientation of this BaseInteractable. More...
 
virtual void setPosition (const Vec3D &position)
 Sets the position of this BaseInteractable. More...
 
void setOrientationViaNormal (Vec3D normal)
 Sets the orientation of this BaseInteractable by defining the vector that results from the rotation of the (1,0,0) vector. More...
 
void setOrientationViaEuler (Vec3D eulerAngle)
 Sets the orientation of this BaseInteractable by defining the euler angles. More...
 
virtual void setOrientation (const Quaternion &orientation)
 Sets the orientation of this BaseInteractable. More...
 
virtual void move (const Vec3D &move)
 Moves this BaseInteractable by adding an amount to the position. More...
 
virtual void rotate (const Vec3D &angularVelocityDt)
 Rotates this BaseInteractable. More...
 
const std::vector< BaseInteraction * > & getInteractions () const
 Returns a list of interactions which belong to this interactable. More...
 
void addInteraction (BaseInteraction *I)
 Adds an interaction to this BaseInteractable. More...
 
bool removeInteraction (BaseInteraction *I)
 Removes an interaction from this BaseInteractable. More...
 
void copyInteractionsForPeriodicParticles (const BaseInteractable &p)
 Copies interactions to this BaseInteractable whenever a periodic copy made. More...
 
void setVelocity (const Vec3D &velocity)
 set the velocity of the BaseInteractable. More...
 
void setAngularVelocity (const Vec3D &angularVelocity)
 set the angular velocity of the BaseInteractble. More...
 
void addVelocity (const Vec3D &velocity)
 adds an increment to the velocity. More...
 
void addAngularVelocity (const Vec3D &angularVelocity)
 add an increment to the angular velocity. More...
 
virtual const Vec3DgetVelocity () const
 Returns the velocity of this interactable. More...
 
virtual const Vec3DgetAngularVelocity () const
 Returns the angular velocity of this interactable. More...
 
void setPrescribedPosition (const std::function< Vec3D(double)> &prescribedPosition)
 Allows the position of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedPosition (double time)
 Computes the position from the user defined prescribed position function. More...
 
void setPrescribedVelocity (const std::function< Vec3D(double)> &prescribedVelocity)
 Allows the velocity of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedVelocity (double time)
 Computes the velocity from the user defined prescribed velocity function. More...
 
void setPrescribedOrientation (const std::function< Quaternion(double)> &prescribedOrientation)
 Allows the orientation of the infinite mass interactbale to be prescribed. More...
 
void applyPrescribedOrientation (double time)
 Computes the orientation from the user defined prescribed orientation function. More...
 
void setPrescribedAngularVelocity (const std::function< Vec3D(double)> &prescribedAngularVelocity)
 Allows the angular velocity of the infinite mass interactable to be prescribed. More...
 
void applyPrescribedAngularVelocity (double time)
 Computes the angular velocity from the user defined prescribed angular velocity. More...
 
virtual BaseInteractiongetInteractionWith (BaseParticle *P, unsigned timeStamp, InteractionHandler *interactionHandler)=0
 Returns the interaction between this object and a given BaseParticle. More...
 
virtual const Vec3D getVelocityAtContact (const Vec3D &contact) const
 Returns the velocity at the contact point, use by many force laws. More...
 
void integrateBeforeForceComputation (double time, double timeStep)
 This is part of integrate routine for objects with infinite mass. More...
 
void integrateAfterForceComputation (double time, double timeStep)
 This is part of the integration routine for objects with infinite mass. More...
 
virtual bool isFixed () const =0
 
virtual Mdouble getInvMass () const
 
virtual Mdouble getCurvature (const Vec3D &labFixedCoordinates) const
 
virtual bool isFaceContact (const Vec3D &normal) const
 
- Public Member Functions inherited from BaseObject
 BaseObject ()=default
 Default constructor. More...
 
 BaseObject (const BaseObject &p)=default
 Copy constructor, copies all the objects BaseObject contains. More...
 
virtual ~BaseObject ()=default
 virtual destructor More...
 
virtual std::string getName () const =0
 A purely virtual function. More...
 
virtual void moveInHandler (unsigned int index)
 Except that it is virtual, it does the same thing as setIndex() does. More...
 
void setIndex (unsigned int index)
 Allows one to assign an index to an object in the handler/container. More...
 
void setId (unsigned long id)
 Assigns a unique identifier to each object in the handler (container) which remains constant even after the object is deleted from the container/handler. More...
 
unsigned int getIndex () const
 Returns the index of the object in the handler. More...
 
unsigned int getId () const
 Returns the unique identifier of any particular object. More...
 
void setGroupId (unsigned groupId)
 
unsigned getGroupId () const
 

Private Attributes

std::function< Vec3D(double)> prescribedPosition_
 
std::function< Vec3D(double)> prescribedVelocity_
 
std::function< Quaternion(double)> prescribedOrientation_
 
std::function< Vec3D(double)> prescribedAngularVelocity_
 
Vec3D position_
 
Quaternion orientation_
 
Vec3D angularVelocity_
 
Vec3D force_
 
Vec3D torque_
 
std::vector< Vec3DforceOMP_
 
std::vector< Vec3DtorqueOMP_
 
const ParticleSpeciesspecies_
 
unsigned int indSpecies_
 
Vec3D velocity_
 
std::vector< BaseInteraction * > interactions_
 

Detailed Description

Defines the basic properties that a interactable object can have.

Inherits from class BaseObject (public) Also it includes a lot of code to deal with interactable objects that have a prescibed motion. Most of the code in here is MercuryDPM internal. The only place an user will interface with this code is for setting the lambda functions that prescribe the motion of infinite mass particles.

Todo:
Check prescribed objects have infinite mass.

Constructor & Destructor Documentation

◆ BaseInteractable() [1/2]

BaseInteractable::BaseInteractable ( )

Default BaseInteractable constructor.

Simply creates an empty BaseInteractable, with all vectors relating to the positions and motions of the current object initialised to zero and all pointers to null.

Note that the function also sets the species index (indSpecies_) to zero by default, so any objects created will, by default, possess the properties associated with species 0.

42  :
43  BaseObject()
44 {
45  //setting all vectors to zero
50  force_.setZero();
51  torque_.setZero();
52  //setting the default species to species 0
53  indSpecies_ = 0;
54  //setting all relevant pointers to nullptr
55  species_ = nullptr;
56  prescribedPosition_ = nullptr;
57  prescribedVelocity_ = nullptr;
58  prescribedOrientation_ = nullptr;
60  logger(DEBUG, "BaseInteractable::BaseInteractable() finished");
61 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ DEBUG
unsigned int indSpecies_
Definition: BaseInteractable.h:497
std::function< Vec3D(double)> prescribedAngularVelocity_
Definition: BaseInteractable.h:442
std::function< Vec3D(double)> prescribedVelocity_
Definition: BaseInteractable.h:430
const ParticleSpecies * species_
Definition: BaseInteractable.h:492
std::function< Vec3D(double)> prescribedPosition_
Definition: BaseInteractable.h:424
Vec3D position_
Definition: BaseInteractable.h:448
Vec3D velocity_
Definition: BaseInteractable.h:502
Quaternion orientation_
Definition: BaseInteractable.h:454
Vec3D torque_
Definition: BaseInteractable.h:469
Vec3D angularVelocity_
Definition: BaseInteractable.h:459
std::function< Quaternion(double)> prescribedOrientation_
Definition: BaseInteractable.h:436
Vec3D force_
Definition: BaseInteractable.h:464
BaseObject()=default
Default constructor.
void setUnity()
Sets quaternion value to (1,0,0,0)
Definition: Quaternion.cc:53
void setZero()
Sets all elements to zero.
Definition: Vector.cc:43

References angularVelocity_, DEBUG, force_, indSpecies_, logger, orientation_, position_, prescribedAngularVelocity_, prescribedOrientation_, prescribedPosition_, prescribedVelocity_, Quaternion::setUnity(), Vec3D::setZero(), species_, torque_, and velocity_.

◆ BaseInteractable() [2/2]

BaseInteractable::BaseInteractable ( const BaseInteractable p)

Copy constructor.

Copies an existing BaseInteractable, p, and (almost) all objects it contains.

Note, that the interactions are not copied as these often require extra work. Rather, the interactions list is simply emptied using the standard C++ clear function, destroying all contents and leaving the interactions_ vector with a size of zero.

All the other properties are copied normally.

Please use this copy with care.

76  : BaseObject(p)
77 {
78  interactions_.clear();
79  position_ = p.position_;
81  velocity_ = p.velocity_;
83  force_ = p.force_;
84  torque_ = p.torque_;
85  species_ = p.species_;
91  logger(DEBUG, "BaseInteractable::BaseInteractable(const BaseInteractable &p finished");
92 }
std::vector< BaseInteraction * > interactions_
Definition: BaseInteractable.h:507

References angularVelocity_, DEBUG, force_, indSpecies_, interactions_, logger, orientation_, position_, prescribedAngularVelocity_, prescribedOrientation_, prescribedPosition_, prescribedVelocity_, species_, torque_, and velocity_.

◆ ~BaseInteractable()

BaseInteractable::~BaseInteractable ( )
override

Destructor, it simply destructs the BaseInteractable and all the objects it contains.

Removes all the interactions from the interactable.

98 {
99  logger(VERBOSE, "Deleting BaseInteractable with index= % and id = %, size = %", getIndex(), getId(),
100  interactions_.size());
101  while (!interactions_.empty())
102  {
103  interactions_.front()->removeFromHandler();
104  }
105  logger(DEBUG, "BaseInteractable::~BaseInteractable() finished");
106 }
@ VERBOSE
unsigned int getId() const
Returns the unique identifier of any particular object.
Definition: BaseObject.h:125
unsigned int getIndex() const
Returns the index of the object in the handler.
Definition: BaseObject.h:118

References DEBUG, BaseObject::getId(), BaseObject::getIndex(), interactions_, logger, and VERBOSE.

Member Function Documentation

◆ addAngularVelocity()

void BaseInteractable::addAngularVelocity ( const Vec3D angularVelocity)

add an increment to the angular velocity.

See also BaseInteractable::setAngularVelocity

Parameters
[in]angularVelocityincrement which to increase the angularVelocity by.
371 {
372  angularVelocity_ += angularVelocity;
373 }

References angularVelocity_.

Referenced by BaseParticle::angularAccelerate(), and TriangleMeshWall::getInteractionWith().

◆ addForce()

void BaseInteractable::addForce ( const Vec3D addForce)

Adds an amount to the force on this BaseInteractable.

Incremental version of BaseInteractable::setForce. Also see BaseInteraction::setForce for were this is used.

Parameters
[in]addForceVec3D incremental force which is added to the total force of the interactable.
117 {
118  if (OMP_THREAD_NUM==0) {
119  force_ += addForce;
120  } else {
122  }
123 }
#define OMP_THREAD_NUM
Definition: GeneralDefine.h:69
std::vector< Vec3D > forceOMP_
Definition: BaseInteractable.h:481
void addForce(const Vec3D &addForce)
Adds an amount to the force on this BaseInteractable.
Definition: BaseInteractable.cc:116

References force_, forceOMP_, and OMP_THREAD_NUM.

Referenced by MeshTriangle::checkInteractions(), DPMBase::computeExternalForces(), Contact::computeExternalForces(), VolumeCoupling::computeExternalForces(), AngledPeriodicBoundarySecondUnitTest::computeExternalForces(), ScaleCoupling< M, O >::computeExternalForces(), DPMBase::computeForcesDueToWalls(), Mercury3Dclump::computeForcesDueToWalls(), DPMBase::computeInternalForce(), Mercury3Dclump::computeInternalForce(), ChuteWithPeriodicInflow::computeInternalForces(), BaseWall::getInteractionWith(), and TriangleMeshWall::getInteractionWith().

◆ addInteraction()

void BaseInteractable::addInteraction ( BaseInteraction I)

Adds an interaction to this BaseInteractable.

Added a new interactions to the current interactable.

Parameters
[in]IPointer to the new interaction which is to be added to the list of interactions of this interactable.
293 {
294  interactions_.push_back(I);
295 }

References interactions_.

Referenced by InteractionHandler::addObject(), BaseInteraction::importI(), BaseInteraction::importP(), BaseInteraction::setI(), and BaseInteraction::setP().

◆ addTorque()

void BaseInteractable::addTorque ( const Vec3D addTorque)

Adds an amount to the torque on this BaseInteractable.

Incremental version of BaseInteractable::setTorque. Also see BaseInteraction::setTorque for were this is used.

Parameters
[in]addTorqueVec3D incremental force which is added to the total torque of the interactable.
133 {
134  if (OMP_THREAD_NUM==0) {
135  torque_ += addTorque;
136  } else {
138  }
139 }
std::vector< Vec3D > torqueOMP_
Definition: BaseInteractable.h:486
void addTorque(const Vec3D &addTorque)
Adds an amount to the torque on this BaseInteractable.
Definition: BaseInteractable.cc:132

References OMP_THREAD_NUM, torque_, and torqueOMP_.

Referenced by MeshTriangle::checkInteractions(), DPMBase::computeForcesDueToWalls(), Mercury3Dclump::computeForcesDueToWalls(), DPMBase::computeInternalForce(), Mercury3Dclump::computeInternalForce(), ChuteWithPeriodicInflow::computeInternalForces(), BaseWall::getInteractionWith(), and TriangleMeshWall::getInteractionWith().

◆ addVelocity()

void BaseInteractable::addVelocity ( const Vec3D velocity)
inline

adds an increment to the velocity.

See also BaseInteractable::setVelocity

Parameters
[in]velocityVec3D containing the velocity increment which to increase the velocity by.
313  { velocity_ += velocity; }

References velocity_.

Referenced by BaseParticle::accelerate(), TriangleMeshWall::getInteractionWith(), TimeDependentPeriodicBoundary::shiftAndBoostParticle(), ShearBoxBoundary::shiftHorizontalPosition(), LeesEdwardsBoundary::shiftVerticalPosition(), and ClumpParticle::updatePebblesVelPos().

◆ applyPrescribedAngularVelocity()

void BaseInteractable::applyPrescribedAngularVelocity ( double  time)

Computes the angular velocity from the user defined prescribed angular velocity.

This calls the prescribedAngularVelocity function if one has been defined. See also BaseInteractable::setPrescribedAngularVelocity

Parameters
[in]timedouble which is the current time of the simulation.
513 {
515  {
517  }
518 }
void setAngularVelocity(const Vec3D &angularVelocity)
set the angular velocity of the BaseInteractble.
Definition: BaseInteractable.cc:360

References prescribedAngularVelocity_, and setAngularVelocity().

Referenced by integrateAfterForceComputation(), and integrateBeforeForceComputation().

◆ applyPrescribedOrientation()

void BaseInteractable::applyPrescribedOrientation ( double  time)

Computes the orientation from the user defined prescribed orientation function.

This calls the prescribedOrientation function if one has been defined. See also BaseInteractable::setPrescribedOrientation

Parameters
[in]timedouble which is the current time of the simulation.
486 {
488  {
490  }
491 }
virtual void setOrientation(const Quaternion &orientation)
Sets the orientation of this BaseInteractable.
Definition: BaseInteractable.h:260

References prescribedOrientation_, and setOrientation().

Referenced by integrateBeforeForceComputation().

◆ applyPrescribedPosition()

void BaseInteractable::applyPrescribedPosition ( double  time)

Computes the position from the user defined prescribed position function.

This calls the prescribedPosition function if one has been defined. See also BaseInteractable::setPrescribedPosition

Parameters
[in]timedouble which is the current time of the simulation.
424 {
426  {
428  }
429 }
virtual void setPosition(const Vec3D &position)
Sets the position of this BaseInteractable.
Definition: BaseInteractable.h:239

References prescribedPosition_, and setPosition().

Referenced by integrateBeforeForceComputation().

◆ applyPrescribedVelocity()

void BaseInteractable::applyPrescribedVelocity ( double  time)

Computes the velocity from the user defined prescribed velocity function.

This calls the prescribedVeclocity function if one has been defined. See also BaseInteractable::setPrescribedVelocity

Parameters
[in]timedouble which is the current time of the simulation.
455 {
457  {
459  }
460 }
void setVelocity(const Vec3D &velocity)
set the velocity of the BaseInteractable.
Definition: BaseInteractable.cc:350

References prescribedVelocity_, and setVelocity().

Referenced by integrateAfterForceComputation(), and integrateBeforeForceComputation().

◆ copyInteractionsForPeriodicParticles()

void BaseInteractable::copyInteractionsForPeriodicParticles ( const BaseInteractable pOriginal)

Copies interactions to this BaseInteractable whenever a periodic copy made.

This loops over all interactions of periodic (particle) and calls copySwitchPointer, which copies the interactions.

Parameters
[in]pOriginalReference to the BaseInteractable which is to be copied to create the ghost particles.
387 {
388  for (BaseInteraction* interaction : pOriginal.interactions_)
389  {
390  //So here this is the ghost and it is the interaction of the ghost/
391  interaction->copySwitchPointer(&pOriginal, this);
392  }
393 }
Stores information about interactions between two interactable objects; often particles but could be ...
Definition: BaseInteraction.h:60
void copySwitchPointer(const BaseInteractable *original, BaseInteractable *ghost) const
This copies the interactions of the original particle and replaces the original with the ghost copy.
Definition: BaseInteraction.cc:341

References BaseInteraction::copySwitchPointer(), and interactions_.

Referenced by ConstantMassFlowMaserBoundary::createGhostCopy(), SubcriticalMaserBoundary::createGhostCopy(), PeriodicBoundary::createGhostParticle(), TimeDependentPeriodicBoundary::createGhostParticle(), LeesEdwardsBoundary::createHorizontalPeriodicParticle(), ShearBoxBoundary::createHorizontalPeriodicParticle(), AngledPeriodicBoundary::createPeriodicParticle(), LeesEdwardsBoundary::createVerticalPeriodicParticle(), and ShearBoxBoundary::createVerticalPeriodicParticle().

◆ getAngularVelocity()

const Vec3D & BaseInteractable::getAngularVelocity ( ) const
virtual

◆ getCurvature()

virtual Mdouble BaseInteractable::getCurvature ( const Vec3D labFixedCoordinates) const
inlinevirtual

returns the inverse radius, or curvature, of the surface. This value is zero for walls and gets overridden for particles that have finite radius

Todo:
should be wall-type dependent

Reimplemented in BaseParticle, and SuperQuadricParticle.

415  { return 0.0; }

Referenced by BaseInteraction::getEffectiveRadius(), and BaseInteraction::getOverlapVolume().

◆ getForce()

const Vec3D& BaseInteractable::getForce ( ) const
inline

Returns the force on this BaseInteractable.

Return the current force being to the BaseInteractable. Note, the code works by first computing the forces of each interaction and then it loops over all BaseInteracables applying forces to them from the interactions they are involved in.

Returns
const Vec3D reference that is the total force applied to this interactable.
127  { return force_; }

References force_.

Referenced by T_protectiveWall::actionsAfterTimeStep(), ClosedCSCWalls::actionsAfterTimeStep(), MercuryProblem::actionsAfterTimeStep(), protectiveWall::actionsAfterTimeStep(), Slide::actionsBeforeTimeStep(), SphericalIndenter::getForceOnIndenter(), BaseParticle::integrateAfterForceComputation(), ClumpParticle::integrateAfterForceComputation(), BaseParticle::integrateBeforeForceComputation(), ClumpParticle::integrateBeforeForceComputation(), main(), BaseWall::setForceControl(), BaseWall::setVelocityControl(), ClumpParticle::updateExtraQuantities(), and Slide::writeEneTimeStep().

◆ getIndSpecies()

◆ getInteractions()

◆ getInteractionWith()

virtual BaseInteraction* BaseInteractable::getInteractionWith ( BaseParticle P,
unsigned  timeStamp,
InteractionHandler interactionHandler 
)
pure virtual

Returns the interaction between this object and a given BaseParticle.

Todo:

TW make sure this function sets normal, distance, overlap, contact point

AT why is this a BaseParticle and not a BaseInteratable.

Implemented in VChute, TriangulatedWall, TriangleMeshWall, SineWall, RestrictedWall, ParabolaChute, MeshTriangle, Combtooth, BaseWall, ArcWall, SuperQuadricParticle, and BaseParticle.

Referenced by PeriodicBoundaryHandler::processReceivedInteractionData(), and Domain::processReceivedInteractionData().

◆ getInvMass()

virtual Mdouble BaseInteractable::getInvMass ( ) const
inlinevirtual

returns the inverse mass. This value is zero for walls and gets overridden for particles that have finite mass

Reimplemented in MeshTriangle, and BaseParticle.

409  { return 0.0; }

Referenced by BaseInteraction::getEffectiveMass().

◆ getOrientation()

const Quaternion& BaseInteractable::getOrientation ( ) const
inline

Returns the orientation of this BaseInteractable.

Returns the reference to a Vec3D which contains the orientation of the interactionable. Please note the interpretation of this depends on which interactable. Please see derived objects for details.

Returns
Returns a reference to a Vec3D returns the position of the interactable.
231  { return orientation_; }

References orientation_.

Referenced by MarbleRun::actionsAfterTimeStep(), NautaMixer::addScrew(), compareParticles(), ScrewsymmetricIntersectionOfWalls::computeDeltaZ(), SuperQuadricParticle::computeHessianLabFixed(), SuperQuadricParticle::computeShape(), SuperQuadricParticle::computeShapeGradientLabFixed(), WearableTriangleMeshWall::computeWear(), AxisymmetricIntersectionOfWalls::convertLimits(), HorizontalBaseScrew::convertLimits(), ScrewsymmetricIntersectionOfWalls::convertLimits(), MPISphericalParticle::copyDataFromParticleToMPIParticle(), copyPositionFrom(), InfiniteWall::createVTK(), BaseWall::getAxis(), SuperQuadricParticle::getContactPointPlanB(), InfiniteWall::getDistance(), AxisymmetricIntersectionOfWalls::getDistanceAndNormal(), HorizontalBaseScrew::getDistanceAndNormal(), NurbsWall::getDistanceAndNormal(), Screw::getDistanceAndNormal(), BasicIntersectionOfWalls::getDistanceAndNormal(), BasicUnionOfWalls::getDistanceAndNormal(), InfiniteWall::getDistanceAndNormal(), IntersectionOfWalls::getDistanceAndNormal(), LevelSetWall::getDistanceAndNormal(), TriangleMeshWall::getDistanceAndNormal(), InfiniteWall::getDistanceNormalOverlapSuperquadric(), BaseWall::getInteractionWith(), TriangleMeshWall::getInteractionWith(), InfiniteWall::getNormal(), BasicIntersectionOfWalls::getVTK(), BasicUnionOfWalls::getVTK(), integrateAfterForceComputation(), BaseParticle::integrateAfterForceComputation(), integrateBeforeForceComputation(), BaseParticle::integrateBeforeForceComputation(), objectivenessTest(), ChuteWithPeriodicInflow::outputXBallsDataParticlee(), CGFields::OrientationField::setFields(), EllipticalSuperQuadricCollision::setupInitialConditions(), WearableNurbsWall::storeDebris(), ShapeGradientHessianTester::testCushion(), ShapeGradientHessianTester::testEllipsoid(), ShapeGradientHessianTester::testRoundedBeam(), TriangleMeshWall::updateBoundingBoxGlobal(), PeriodicBoundaryHandler::updateParticles(), ClumpParticle::updatePebblesVelPos(), TriangleWall::updateVertexAndNormal(), SuperQuadricParticle::writeDebugMessageMiddleOfLoop(), SuperQuadricParticle::writeDebugMessageStep1(), AxisymmetricIntersectionOfWalls::writeVTK(), BasicUnionOfWalls::writeVTK(), IntersectionOfWalls::writeVTK(), LevelSetWall::writeVTK(), NurbsWall::writeVTK(), Screw::writeVTK(), ScrewsymmetricIntersectionOfWalls::writeVTK(), TriangleMeshWall::writeVTK(), NurbsWall::writeWallDetailsVTK(), and WearableNurbsWall::writeWallDetailsVTK().

◆ getPosition()

const Vec3D& BaseInteractable::getPosition ( ) const
inline

Returns the position of this BaseInteractable.

Returns the reference to a Vec3D which contains the position of the interactionable. Please note the interpretation of this depends on which interactable. For particles this is the centre of the particle; where for walls it is one point of the wall given \(r.n=p\)

Returns
Returns a reference to a Vec3D returns the position of the interactable.
219  { return position_; }

References position_.

Referenced by ClosedCSCRestart::actionsAfterTimeStep(), ClosedCSCRun::actionsAfterTimeStep(), multiParticleT1::actionsAfterTimeStep(), BouncingSuperQuadric::actionsAfterTimeStep(), DrivenParticleClass::actionsAfterTimeStep(), AngleOfRepose::actionsBeforeTimeStep(), Slide::actionsBeforeTimeStep(), Chutebelt::actionsBeforeTimeStep(), NautaMixer::addParticlesAtWall(), NautaMixer::addScrew(), statistics_while_running< T >::auto_set_domain(), statistics_while_running< T >::auto_set_z(), ConstantMassFlowMaserBoundary::checkBoundaryAfterParticleMoved(), DeletionBoundary::checkBoundaryAfterParticleMoved(), FluxBoundary::checkBoundaryAfterParticleMoved(), HeaterBoundary::checkBoundaryAfterParticleMoved(), SubcriticalMaserBoundaryTEST::checkBoundaryAfterParticleMoved(), Mercury3Dclump::checkClumpForInteractionPeriodic(), DPMBase::checkParticleForInteractionLocal(), DPMBase::checkParticleForInteractionLocalPeriodic(), ChuteWithPeriodicInflow::cleanChute(), ChuteWithContraction::cleanChute(), Funnel::cleanChute(), Chute::cleanChute(), compareParticles(), Contact::computeExternalForces(), AngledPeriodicBoundarySecondUnitTest::computeExternalForces(), ScaleCoupling< M, O >::computeExternalForces(), DPMBase::computeForcesDueToWalls(), Mercury3Dclump::computeForcesDueToWalls(), SuperQuadricParticle::computeHessianLabFixed(), DPMBase::computeInternalForce(), Mercury3Dclump::computeInternalForce(), Mercury2D::computeInternalForces(), Mercury3D::computeInternalForces(), ChuteWithPeriodicInflow::computeInternalForces(), DPM::computeLocalCGHGrid(), DPM::computeLocalVolumeFractionHGrid(), SuperQuadricParticle::computeShape(), SuperQuadricParticle::computeShapeGradientLabFixed(), WearableTriangleMeshWall::computeWear(), Domain::containsParticle(), AxisymmetricIntersectionOfWalls::convertLimits(), HorizontalBaseScrew::convertLimits(), ScrewsymmetricIntersectionOfWalls::convertLimits(), MPISphericalParticle::copyDataFromParticleToMPIParticle(), copyPositionFrom(), Funnel::create_funnel(), Funnel::create_inflow_particle(), CircularPeriodicBoundary::createPeriodicParticle(), SubcriticalMaserBoundaryTEST::createPeriodicParticle(), InfiniteWall::createVTK(), ChuteWithWedge::extendBottom(), Domain::findNearbyBoundaries(), Domain::findNewMPIInteractions(), PeriodicBoundaryHandler::findNewParticle(), Domain::flushParticlesFromList(), BaseInteraction::gatherContactStatistics(), BaseParticle::getAngularMomentum(), SuperQuadricParticle::getContactPointPlanB(), BaseInteraction::getCP(), BaseParticle::getDisplacement2(), ConstantMassFlowMaserBoundary::getDistance(), SubcriticalMaserBoundary::getDistance(), PeriodicBoundary::getDistance(), TimeDependentPeriodicBoundary::getDistance(), SphericalWall::getDistance(), InfiniteWall::getDistance(), ArcWall::getDistanceAndNormal(), AxisymmetricIntersectionOfWalls::getDistanceAndNormal(), HorizontalBaseScrew::getDistanceAndNormal(), HorizontalScrew::getDistanceAndNormal(), NurbsWall::getDistanceAndNormal(), Screw::getDistanceAndNormal(), ScrewsymmetricIntersectionOfWalls::getDistanceAndNormal(), BasicIntersectionOfWalls::getDistanceAndNormal(), BasicUnionOfWalls::getDistanceAndNormal(), Coil::getDistanceAndNormal(), Combtooth::getDistanceAndNormal(), InfiniteWall::getDistanceAndNormal(), IntersectionOfWalls::getDistanceAndNormal(), LevelSetWall::getDistanceAndNormal(), ParabolaChute::getDistanceAndNormal(), RestrictedWall::getDistanceAndNormal(), SimpleDrumSuperquadrics::getDistanceAndNormal(), SineWall::getDistanceAndNormal(), SphericalWall::getDistanceAndNormal(), TriangleMeshWall::getDistanceAndNormal(), TriangleWall::getDistanceAndNormal(), VChute::getDistanceAndNormal(), TriangulatedWall::Face::getDistanceAndNormal(), InfiniteWall::getDistanceNormalOverlapSuperquadric(), SimpleDrumSuperquadrics::getDistanceNormalOverlapSuperquadric(), MeshTriangle::getDistanceNormalOverlapType(), LeesEdwardsBoundary::getHorizontalDistance(), ShearBoxBoundary::getHorizontalDistance(), BaseInteraction::getIC(), SphericalIndenter::getIndenterHeight(), SuperQuadricParticle::getInitialGuessForContact(), CGCoordinates::R::getINormal(), CGCoordinates::RZ::getINormal(), CGCoordinates::X::getINormal(), CGCoordinates::XY::getINormal(), CGCoordinates::XYZ::getINormal(), CGCoordinates::XZ::getINormal(), CGCoordinates::Y::getINormal(), CGCoordinates::YZ::getINormal(), CGCoordinates::Z::getINormal(), BaseParticle::getInteractionWith(), SuperQuadricParticle::getInteractionWith(), ArcWall::getInteractionWith(), BaseWall::getInteractionWith(), Combtooth::getInteractionWith(), MeshTriangle::getInteractionWith(), RestrictedWall::getInteractionWith(), SineWall::getInteractionWith(), TriangleMeshWall::getInteractionWith(), TriangulatedWall::getInteractionWith(), VChute::getInteractionWith(), SuperQuadricParticle::getInteractionWithSuperQuad(), DomainHandler::getParticleDomainGlobalIndex(), ClumpParticle::getPebblePositions(), CGCoordinates::R::getPNormal(), CGCoordinates::RZ::getPNormal(), CGCoordinates::X::getPNormal(), CGCoordinates::XY::getPNormal(), CGCoordinates::XZ::getPNormal(), CGCoordinates::Y::getPNormal(), CGCoordinates::YZ::getPNormal(), CGCoordinates::Z::getPNormal(), CGCoordinates::RZ::getTangentialSquared(), CGCoordinates::XY::getTangentialSquared(), CGCoordinates::XYZ::getTangentialSquared(), CGCoordinates::XZ::getTangentialSquared(), CGCoordinates::YZ::getTangentialSquared(), getVelocityAtContact(), LeesEdwardsBoundary::getVerticalDistance(), ShearBoxBoundary::getVerticalDistance(), Mercury3D::hGridFindContactsWithTargetCell(), Mercury2D::hGridFindParticleContacts(), Mercury3D::hGridFindParticleContacts(), Mercury2D::hGridGetInteractingParticleList(), Mercury3D::hGridGetInteractingParticleList(), Mercury2D::hGridHasParticleContacts(), Mercury3D::hGridHasParticleContacts(), Mercury2D::hGridUpdateParticle(), Mercury3D::hGridUpdateParticle(), BaseParticle::integrateBeforeForceComputation(), ClumpParticle::integrateBeforeForceComputation(), HorizontalMixer::introduceParticlesAtWall(), ContactDetectionIntersectionOfWallsTest::introduceParticlesAtWall(), PeriodicBoundary::isClosestToLeftBoundary(), TimeDependentPeriodicBoundary::isClosestToLeftBoundary(), ConstantMassFlowMaserBoundary::isClosestToRightBoundary(), SubcriticalMaserBoundary::isClosestToRightBoundary(), BaseParticle::isInContactWith(), LawinenBox::LawinenBox(), load(), main(), SubcriticalMaserBoundaryTEST::modifyPeriodicComplexity(), objectivenessTest(), SphericalIndenter::outputXBallsData(), ChuteWithPeriodicInflow::outputXBallsDataParticlee(), CubeDeletionBoundarySelfTest::printTime(), DeletionBoundarySelfTest::printTime(), PeriodicBoundaryHandler::processLocalGhostParticles(), Domain::processReceivedBoundaryParticleData(), PeriodicBoundaryHandler::processReceivedGhostParticleData(), FileReader::read(), regimeForceUnitTest::regimeForceUnitTest(), save(), Membrane::saveVertexPositions(), ClosedCSCWalls::saveWalls(), CGFields::GradVelocityField::setCylindricalFields(), CGFields::StandardFields::setCylindricalFields(), IntersectionOfWalls::setPointsAndLines(), MarbleRun::setupInitialConditions(), AngleOfRepose::setupInitialConditions(), NewtonsCradleSelfTest::setupInitialConditions(), CoilSelfTest::setupInitialConditions(), Packing::setupInitialConditions(), TriangleWall::setVertices(), ShearBoxBoundary::shiftHorizontalPosition(), SinterPair::SinterPair(), WearableNurbsWall::storeDebris(), FlowFrontChute::stretch(), TriangleMeshWall::updateBoundingBoxGlobal(), PeriodicBoundaryHandler::updateMaserParticle(), PeriodicBoundaryHandler::updateParticles(), PeriodicBoundaryHandler::updateParticleStatus(), ClumpParticle::updatePebblesVelPos(), SCoupling< M, O >::updateTriangleWall(), TriangleWall::updateVertexAndNormal(), MeshTriangle::updateVerticesFromParticles(), viscoElasticUnitTest::viscoElasticUnitTest(), SuperQuadricParticle::writeDebugMessageMiddleOfLoop(), SuperQuadricParticle::writeDebugMessageStep1(), Penetration::writeEneTimeStep(), Slide::writeEneTimeStep(), SingleParticle< SpeciesType >::writeEneTimeStep(), AxisymmetricIntersectionOfWalls::writeVTK(), BasicUnionOfWalls::writeVTK(), HorizontalBaseScrew::writeVTK(), IntersectionOfWalls::writeVTK(), LevelSetWall::writeVTK(), NurbsWall::writeVTK(), Screw::writeVTK(), ScrewsymmetricIntersectionOfWalls::writeVTK(), TriangleMeshWall::writeVTK(), NurbsWall::writeWallDetailsVTK(), and WearableNurbsWall::writeWallDetailsVTK().

◆ getSpecies()

const ParticleSpecies* BaseInteractable::getSpecies ( ) const
inline

Returns a pointer to the species of this BaseInteractable.

This function return a ParticleSpecies* for the current interacable. Please note, this is a ParticleSpecies; not, a BaseSpecies as interactables must have physically properties as well.

Returns
constant ParticleSpecies* pointer to the species storing the physical properties of this interactable.
109  {
110  return species_;
111  }

References species_.

Referenced by ClumpParticle::actionsAfterAddObject(), WallHandler::addObject(), IntersectionOfWalls::addObject(), ConstantMassFlowMaserBoundary::addParticleToMaser(), SubcriticalMaserBoundary::addParticleToMaser(), NautaMixer::addScrew(), ConstantMassFlowMaserBoundary::checkBoundaryAfterParticleMoved(), SubcriticalMaserBoundary::checkBoundaryAfterParticleMoved(), ChargedBondedInteraction::computeAdhesionForce(), CurvyChute::createBottom(), NurbsWall::getDistanceAndNormal(), Screw::getDistanceAndNormal(), BasicIntersectionOfWalls::getDistanceAndNormal(), BasicUnionOfWalls::getDistanceAndNormal(), IntersectionOfWalls::getDistanceAndNormal(), TriangleMeshWall::getDistanceAndNormal(), ChargedBondedInteraction::getElasticEnergy(), BaseParticle::getInfo(), BaseParticle::getInteractionDistance(), SuperQuadricParticle::getInteractionRadius(), BaseParticle::getMaxInteractionRadius(), BasicIntersectionOfWalls::getVTK(), BasicUnionOfWalls::getVTK(), ConstantMassFlowMaserBoundary::isMaserParticle(), SubcriticalMaserBoundary::isMaserParticle(), ConstantMassFlowMaserBoundary::isNormalParticle(), SubcriticalMaserBoundary::isNormalParticle(), DPMBase::readNextDataFile(), TriangleMeshWall::refineTriangle(), ConstantMassFlowMaserBoundary::removeParticleFromMaser(), SubcriticalMaserBoundary::removeParticleFromMaser(), RestrictedWall::set(), TriangleMeshWall::set(), SuperQuadricParticle::setAxes(), LinearViscoelasticNormalSpecies::setCollisionTimeAndRestitutionCoefficient(), SPHNormalSpecies::setCollisionTimeAndRestitutionCoefficient(), SuperQuadricParticle::setExponents(), SuperQuadricParticle::setInertia(), BaseParticle::setRadius(), WallParticleCollision::setupInitialConditions(), BaseParticle::unfix(), and BaseInteraction::writeInteraction().

◆ getTorque()

const Vec3D& BaseInteractable::getTorque ( ) const
inline

Returns the torque on this BaseInteractable.

Return the current torque being to the BaseInteractable. Note, the code works by first computing the forces of each interaction and then it loops over all BaseInteracables applying forces to them from the interactions they are involved in.

Returns
const Vec3D reference that is the total force applied to this interactable.
139  { return torque_; }

References torque_.

Referenced by ClumpParticle::angularAccelerateClumpIterative(), BaseParticle::integrateAfterForceComputation(), BaseParticle::integrateBeforeForceComputation(), and ClumpParticle::updateExtraQuantities().

◆ getVelocity()

const Vec3D & BaseInteractable::getVelocity ( ) const
virtual

Returns the velocity of this interactable.

Returns the velocity of the BaseInteractbale. Note, this is the same for all BaseInteractables; it is the direction the object is moving in.

Returns
Vec3D reference which contains the current velocity of the current BaseInteractable.
330 {
331  return velocity_;
332 }

References velocity_.

Referenced by HertzianBSHPInteractionTwoParticleElasticCollision::actionsAfterSolve(), SphericalSuperQuadricCollision::actionsAfterSolve(), SpeciesTest::actionsAfterSolve(), ClosedCSCRestart::actionsAfterTimeStep(), ClosedCSCRun::actionsAfterTimeStep(), ClosedCSCWalls::actionsAfterTimeStep(), SphericalIndenter::actionsAfterTimeStep(), BouncingSuperQuadric::actionsAfterTimeStep(), LawinenBox::actionsBeforeTimeStep(), Slide::actionsBeforeTimeStep(), DeletionBoundary::checkBoundaryAfterParticleMoved(), HeaterBoundary::checkBoundaryAfterParticleMoved(), ChuteWithContraction::ChuteWithContraction(), compareParticles(), DPMBase::computeExternalForces(), FrictionInteraction::computeFrictionForce(), ChuteWithPeriodicInflow::computeInternalForces(), ClosedCSCWalls::continueSolve(), MPISphericalParticle::copyDataFromParticleToMPIParticle(), copyVelocityFrom(), BaseParticle::getAngularMomentum(), SphericalIndenter::getIndenterVelocity(), WallSpecies::getInfo(), TriangleMeshWall::getInteractionWith(), BaseParticle::getKineticEnergy(), ClumpParticle::getKineticEnergy(), BaseParticle::getMomentum(), ParticleParticleCollision::getRelativeVelocity(), WallParticleCollision::getRelativeVelocity(), getVelocityAtContact(), ClumpParticle::integrateAfterForceComputation(), integrateBeforeForceComputation(), BaseParticle::integrateBeforeForceComputation(), ClumpParticle::integrateBeforeForceComputation(), load(), main(), InfiniteWallWithHole::move_time(), objectivenessTest(), SphericalIndenter::outputXBallsData(), ChuteWithPeriodicInflow::outputXBallsDataParticlee(), LawinenBox::printTime(), ClosedCSCRestart::printTime(), ClosedCSCRun::printTime(), ClosedCSCWalls::printTime(), ChuteWithPeriodicInflow::printTime(), save(), CGFields::GradVelocityField::setFields(), CGFields::StandardFields::setFields(), FreeCooling2DinWalls::setupInitialConditions(), FreeCooling2DinWallsDemo::setupInitialConditions(), FreeCooling3DDemoProblem::setupInitialConditions(), FreeCooling3DinWallsDemo::setupInitialConditions(), FreeCoolingDemoProblem::setupInitialConditions(), ContactDetectionNormalSpheresTest::test(), HertzContactRestitutionUnitTest::test(), PeriodicBoundaryHandler::updateParticles(), ClumpParticle::updatePebblesVelPos(), MeshTriangle::updateVerticesFromParticles(), and viscoElasticUnitTest::viscoElasticUnitTest().

◆ getVelocityAtContact()

const Vec3D BaseInteractable::getVelocityAtContact ( const Vec3D contact) const
virtual

Returns the velocity at the contact point, use by many force laws.

Reimplemented in MeshTriangle.

376 {
377  return getVelocity() - Vec3D::cross(contact - getPosition(), getAngularVelocity());
378 }
virtual const Vec3D & getAngularVelocity() const
Returns the angular velocity of this interactable.
Definition: BaseInteractable.cc:341
virtual const Vec3D & getVelocity() const
Returns the velocity of this interactable.
Definition: BaseInteractable.cc:329
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:218
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:163

References Vec3D::cross(), getAngularVelocity(), getPosition(), and getVelocity().

◆ integrateAfterForceComputation()

void BaseInteractable::integrateAfterForceComputation ( double  time,
double  timeStep 
)

This is part of the integration routine for objects with infinite mass.

This is the last part of time integration for interactable objects which have an infinite mass.

Parameters
[in]timedouble which is the current simulation time
[in]timeStepdouble which is the current delta time of the simulation i.e. the size of each time step.
612 {
614  {
616  {
617  applyPrescribedVelocity(time + timeStep);
618  }
619  else
620  {
621  setVelocity((prescribedPosition_(time + 1.1 * timeStep) - prescribedPosition_(time + 0.9 * timeStep)) /
622  (0.2 * timeStep));
623  }
624  }
625  else
626  {
628  {
629  applyPrescribedVelocity(time + 0.5 * timeStep);
630  }
631  }
633  {
635  {
636  applyPrescribedAngularVelocity(time + timeStep);
637  }
638  else
639  {
640  setAngularVelocity(getOrientation().applyCInverse(
641  (prescribedOrientation_(time + 1.1 * timeStep) - prescribedOrientation_(time + 0.9 * timeStep)) /
642  (0.2 * timeStep)
643  ));
644  }
645  }
646  else
647  {
649  {
650  applyPrescribedAngularVelocity(time + 0.5 * timeStep);
651  }
652  }
653 }
const Quaternion & getOrientation() const
Returns the orientation of this BaseInteractable.
Definition: BaseInteractable.h:230
void applyPrescribedAngularVelocity(double time)
Computes the angular velocity from the user defined prescribed angular velocity.
Definition: BaseInteractable.cc:512
void applyPrescribedVelocity(double time)
Computes the velocity from the user defined prescribed velocity function.
Definition: BaseInteractable.cc:454

References applyPrescribedAngularVelocity(), applyPrescribedVelocity(), getOrientation(), prescribedAngularVelocity_, prescribedOrientation_, prescribedPosition_, prescribedVelocity_, setAngularVelocity(), and setVelocity().

Referenced by DPMBase::integrateAfterForceComputation(), BaseParticle::integrateAfterForceComputation(), and ClumpParticle::integrateAfterForceComputation().

◆ integrateBeforeForceComputation()

void BaseInteractable::integrateBeforeForceComputation ( double  time,
double  timeStep 
)

This is part of integrate routine for objects with infinite mass.

This does not first part of verlet integration but for objects with an infinite mass i.e. there motion is prescribed and not calculated from the applied forces. First it deals with the translation degrees of freedom and then secondly if deals with the angular degrees of freedom. Note, in both cases if the user has prescribed both positions and velocity these are used. If only only position is prescribed the velocity is computed from a finite difference. If only the velocity is prescribed the position is computed from integrating the velocity.

In the weird case they neither is set. The objects computed velocity is used to update its position.

Parameters
[in]timedouble which is the current simulation time
[in]timeStepdouble which is the current delta time of the simulation i.e. the size of each time step.
539 {
541  {
543  {
544  //Both the velocity and position are defined; as we are using leap
545  //frog method so the velocity is evaluated half a time later.
546  applyPrescribedPosition(time + timeStep);
547  applyPrescribedVelocity(time + 0.5 * timeStep);
548  }
549  else
550  {
551  //Only the position is defined.
552  //Velocity is evaluated from a finite different of the Position
553  //Note, we use 0.5 +- 0.1 timeStep for the velocity eval.
554  applyPrescribedPosition(time + timeStep);
555  setVelocity((prescribedPosition_(time + 0.6 * timeStep) - prescribedPosition_(time + 0.4 * timeStep)) /
556  (0.2 * timeStep));
557  }
558  }
559  else
560  {
562  {
563  //Only the velocity is set. The position is calculated from the
564  //the integral of velocity.
565  applyPrescribedVelocity(time + 0.5 * timeStep);
566  move(getVelocity() * timeStep);
567  }
568  else
569  {
570  //Neither is set move based on the computed velocity of the object.
571  move(getVelocity() * timeStep);
572  }
573  }
575  {
577  {
578  applyPrescribedOrientation(time + timeStep);
579  applyPrescribedAngularVelocity(time + 0.5 * timeStep);
580  }
581  else
582  {
583  applyPrescribedOrientation(time + timeStep);
584  setAngularVelocity(getOrientation().applyCInverse(
585  (prescribedOrientation_(time + 0.6 * timeStep) - prescribedOrientation_(time + 0.4 * timeStep)) /
586  (0.2 * timeStep)
587  ));
588  }
589  }
590  else
591  {
593  {
594  applyPrescribedAngularVelocity(time + 0.5 * timeStep);
595  rotate(getAngularVelocity() * timeStep);
596  }
597  else
598  {
599  rotate(getAngularVelocity() * timeStep);
600  }
601  }
602 }
void applyPrescribedPosition(double time)
Computes the position from the user defined prescribed position function.
Definition: BaseInteractable.cc:423
void applyPrescribedOrientation(double time)
Computes the orientation from the user defined prescribed orientation function.
Definition: BaseInteractable.cc:485
virtual void move(const Vec3D &move)
Moves this BaseInteractable by adding an amount to the position.
Definition: BaseInteractable.cc:215
virtual void rotate(const Vec3D &angularVelocityDt)
Rotates this BaseInteractable.
Definition: BaseInteractable.cc:230

References applyPrescribedAngularVelocity(), applyPrescribedOrientation(), applyPrescribedPosition(), applyPrescribedVelocity(), getAngularVelocity(), getOrientation(), getVelocity(), move(), prescribedAngularVelocity_, prescribedOrientation_, prescribedPosition_, prescribedVelocity_, rotate(), setAngularVelocity(), and setVelocity().

Referenced by DPMBase::integrateBeforeForceComputation(), BaseParticle::integrateBeforeForceComputation(), and ClumpParticle::integrateBeforeForceComputation().

◆ isFaceContact()

virtual bool BaseInteractable::isFaceContact ( const Vec3D normal) const
inlinevirtual

Reimplemented in TriangleWall.

417 {return true;}

Referenced by BaseWall::getInteractionWith().

◆ isFixed()

virtual bool BaseInteractable::isFixed ( ) const
pure virtual

used to distinguish particles which belong to the flow and fixed particles/walls

Implemented in BaseWall, and BaseParticle.

◆ move()

void BaseInteractable::move ( const Vec3D move)
virtual

Moves this BaseInteractable by adding an amount to the position.

Moves (displaces) the interacable a given distance. Note, this just updates the position by the move.

Parameters
[in]moveReference to Vec3D which is the distance to move the interactable.

Reimplemented in TriangulatedWall, TriangleWall, TriangleMeshWall, and MeshTriangle.

216 {
217  position_ += move;
218 }

References position_.

Referenced by SphericalIndenter::actionsAfterTimeStep(), ChuteWithPeriodicInflow::AddContinuingBottom(), ConstantMassFlowMaserBoundary::addParticleToMaser(), ConstantMassFlowMaserBoundary::checkBoundaryAfterParticleMoved(), ChuteWithContraction::ChuteWithContraction(), ChuteWithPeriodicInflowAndContinuingBottom::ChuteWithPeriodicInflowAndContinuingBottom(), ChuteWithPeriodicInflowAndContraction::ChuteWithPeriodicInflowAndContraction(), ChuteWithPeriodicInflowAndVariableBottom::ChuteWithPeriodicInflowAndVariableBottom(), ContractionWithPeriodicInflow::ContractionWithPeriodicInflow(), ConstantMassFlowMaserBoundary::createPeriodicParticle(), ChuteWithPeriodicInflow::ExtendInWidth(), integrateBeforeForceComputation(), BaseParticle::integrateBeforeForceComputation(), ClumpParticle::integrateBeforeForceComputation(), MeshTriangle::move(), TriangleMeshWall::move(), TriangleWall::move(), TriangulatedWall::move(), HopperInsertionBoundary::placeParticle(), ConstantMassFlowMaserBoundary::removeParticleFromMaser(), ScalingTestRun::setupInitialConditions(), TimeDependentPeriodicBoundary::shiftAndBoostParticle(), LeesEdwardsBoundary::shiftHorizontalPosition(), ShearBoxBoundary::shiftHorizontalPosition(), ConstantMassFlowMaserBoundary::shiftPosition(), SubcriticalMaserBoundary::shiftPosition(), PeriodicBoundary::shiftPosition(), LeesEdwardsBoundary::shiftVerticalPosition(), and ShearBoxBoundary::shiftVerticalPosition().

◆ read()

void BaseInteractable::read ( std::istream &  is)
overridevirtual

Reads a BaseInteractable from an input stream.

BaseInteacatable read functions. Reads in all the information about an interacatable. Note, this can be from any istream but would normally be a file See also BaseInteractable::write

Parameters
[in]isstd::istream to which the information is read from.

Implements BaseObject.

Reimplemented in WearableTriangulatedWall, WearableTriangleMeshWall, WearableNurbsWall, VChute, TriangulatedWall, TriangleWall, TriangleMeshWall, SphericalWall, SineWall, Screw, RestrictedWall, ParabolaChute, NurbsWall, MeshTriangle, LevelSetWall, IntersectionOfWalls, InfiniteWallWithHole, InfiniteWall, HorizontalScrew, CylindricalWall, Combtooth, Coil, BasicUnionOfWalls, BasicIntersectionOfWalls, BaseWall, ArcWall, SuperQuadricParticle, ClumpParticle, BaseParticle, SimpleDrumSuperquadrics, ScrewsymmetricIntersectionOfWalls, HorizontalBaseScrew, and AxisymmetricIntersectionOfWalls.

245 {
246  BaseObject::read(is);
247  std::string dummy;
248  is >> dummy >> indSpecies_;
249  is >> dummy >> position_;
250  is >> dummy;
251  //this if-statement is added to read Kasper van der Vaart's data files, which contain an additional variable named positionAbsolute
252  if (dummy == "positionAbsolute")
253  {
254  is >> dummy >> dummy >> dummy >> dummy;
255  }
256  is >> orientation_;
257  is >> dummy >> velocity_;
258  is >> dummy >> angularVelocity_;
259  is >> dummy;
260  if (dummy == "0")
261  is >> dummy;
262  is >> force_;
263  is >> dummy >> torque_;
264 }
virtual void read(std::istream &is)=0
Definition: BaseObject.cc:81

References angularVelocity_, force_, indSpecies_, orientation_, position_, BaseObject::read(), torque_, and velocity_.

Referenced by BaseParticle::read(), and BaseWall::read().

◆ removeInteraction()

bool BaseInteractable::removeInteraction ( BaseInteraction I)

Removes an interaction from this BaseInteractable.

Removes a given interaction form the list of interactions belonging to the current interacable. This functions returns true to the interaction was found and returns false if the given interaction did not exist and the interaction was not removed. Note that this cannot be done with a range-based for, since we need the iterator to erase the interaction.

Parameters
[in]IBaseInteraction pointer which is the interaction to be removed.
Returns
bool True if the interaction was found and removed; false if the interaction did not exist for that interactable.
309 {
310  for (std::vector<BaseInteraction*>::iterator it = interactions_.begin(); it != interactions_.end(); ++it)
311  {
312  if (I == (*it))
313  {
314  interactions_.erase(it);
315  return true;
316  }
317  }
318  logger(WARN, "Error in BaseInteractable::removeInteraction: Interaction could not be removed");
319  return false;
320 }
@ WARN

References interactions_, logger, and WARN.

Referenced by BaseInteraction::importI(), BaseInteraction::importP(), BaseInteraction::setI(), BaseInteraction::setP(), and BaseInteraction::~BaseInteraction().

◆ resetForceTorque()

void BaseInteractable::resetForceTorque ( int  numberOfOMPthreads)
virtual

Reimplemented in TriangleMeshWall.

142 {
143  #ifdef MERCURYDPM_USE_OMP
144  forceOMP_.resize(numberOfOMPthreads-1);
145  torqueOMP_.resize(numberOfOMPthreads-1);
146  for (auto& force : forceOMP_)
147  {
148  force.setZero();
149  }
150 
151  for (auto& torque : torqueOMP_)
152  {
153  torque.setZero();
154  }
155  #endif
156 
157  force_.setZero();
158  torque_.setZero();
159 }

References force_, forceOMP_, Vec3D::setZero(), torque_, and torqueOMP_.

Referenced by SphericalIndenter::actionsBeforeTimeStep(), DPMBase::computeAllForces(), Mercury3Dclump::computeAllForces(), and TriangleMeshWall::resetForceTorque().

◆ rotate()

void BaseInteractable::rotate ( const Vec3D angularVelocityDt)
virtual

Rotates this BaseInteractable.

Rotates the interacable a given solid angle. Note, this just updates the orientation by the angle.

This function has been declared virtual, so it can be overridden for IntersectionOfWalls.

Parameters
[in]angularVelocityDtReference to Vec3D which is the solid angle through which the interactable is rotated.
Todo:
TW the move and rotate functions should only pass the time step, as teh velocity can be accessed directly by the object; this would simplify functions like Screw::rotate

Reimplemented in Screw, TriangleWall, TriangleMeshWall, and MeshTriangle.

231 {
232  //if (!angularVelocityDt.isZero()) {
233  orientation_.updateAngularDisplacement(angularVelocityDt);
234  //}
235 }
void updateAngularDisplacement(Vec3D angularVelocityDt)
Definition: Quaternion.cc:431

References orientation_, and Quaternion::updateAngularDisplacement().

Referenced by VerticalMixerStraightBlades::addBlades(), VerticalMixerAngledBlades::addBlades(), VerticalMixerAngledBlades::addPrettyBlades(), integrateBeforeForceComputation(), BaseParticle::integrateBeforeForceComputation(), ClumpParticle::integrateBeforeForceComputation(), MeshTriangle::rotate(), TriangleMeshWall::rotate(), TriangleWall::rotate(), and Screw::rotate().

◆ setAngularVelocity()

void BaseInteractable::setAngularVelocity ( const Vec3D angularVelocity)

set the angular velocity of the BaseInteractble.

See also BaseInteractable::getAngularVelocity

Parameters
[in]angularVelocityVec3D which is the angularVelocity of the interactable.
361 {
362  angularVelocity_ = angularVelocity;
363 }

References angularVelocity_.

Referenced by DrumRot::actionsBeforeTimeStep(), RotatingDrum::actionsBeforeTimeStep(), DrumRot::actionsOnRestart(), NautaMixer::addScrew(), ClumpParticle::angularAccelerateClumpIterative(), applyPrescribedAngularVelocity(), MPISphericalParticle::copyDataFromMPIParticleToParticle(), DPM::DPM(), BaseParticle::fixParticle(), GranuDrum::GranuDrum(), integrateAfterForceComputation(), integrateBeforeForceComputation(), loadingTest(), main(), normalAndTangentialLoadingTest(), objectivenessTest(), BaseParticle::oldRead(), DPMBase::readNextDataFile(), WallHandler::readTriangleWall(), HorizontalMixer::setScrewWalls(), multiParticleT1::setupInitialConditions(), RotatingDrumWet::setupInitialConditions(), MinimalExampleDrum::setupInitialConditions(), Drum::setupInitialConditions(), clumpTest::setupInitialConditions(), ParticleParticleCollision::setupInitialConditions(), WallParticleCollision::setupInitialConditions(), TriangulatedScrewSelfTest::setupInitialConditions(), RotatingDrum::setupInitialConditions(), Tutorial12::setupInitialConditions(), AngledPeriodicBoundarySecondUnitTest::setupInitialConditions(), AngledPeriodicBoundaryUnitTest::setupInitialConditions(), MovingWalls::setupInitialConditions(), PeriodicWallsWithSlidingFrictionUnitTest::setupInitialConditions(), TangentialSpringEnergyConservationUnitTest::setupInitialConditions(), ChuteBottom::setupInitialConditions(), PeriodicBoundaryHandler::updateParticles(), and ClumpParticle::updatePebblesVelPos().

◆ setForce()

void BaseInteractable::setForce ( const Vec3D force)
inline

Sets the force on this BaseInteractable.

This sets the force being applied to this interactable. Note, first the code computes all forces in the interactions and then loops over all interactable objects applying the forces from the interactions to the interactables involved in the interaction.

Parameters
[in]forceVec3D which is the force to be applied.
150  { force_ = force; }

References force_.

Referenced by SmoothChute::actionsBeforeTimeStep().

◆ setIndSpecies()

virtual void BaseInteractable::setIndSpecies ( unsigned int  indSpecies)
inlinevirtual

Sets the index of the Species of this BaseInteractable.

This set the species associated with this interactable. This function should not be used and BaseInteractable::setSpecies should be used instead. See also BaseInteractable::setSpecies

Reimplemented in BaseWall, and BaseParticle.

99  { indSpecies_ = indSpecies; }

References indSpecies_.

Referenced by BaseParticle::oldRead(), BaseParticle::setIndSpecies(), and BaseWall::setIndSpecies().

◆ setOrientation()

◆ setOrientationViaEuler()

void BaseInteractable::setOrientationViaEuler ( Vec3D  eulerAngle)

Sets the orientation of this BaseInteractable by defining the euler angles.

205 {
206  orientation_.setEuler(eulerAngle);
207 }
void setEuler(const Vec3D &e)
Convert Euler angles to a quaternion. See Wikipedia for details.
Definition: Quaternion.cc:473

References orientation_, and Quaternion::setEuler().

Referenced by DPMBase::readNextDataFile().

◆ setOrientationViaNormal()

void BaseInteractable::setOrientationViaNormal ( Vec3D  normal)

Sets the orientation of this BaseInteractable by defining the vector that results from the rotation of the (1,0,0) vector.

Sets the orientation of the interactable. interpretation depends on which interactable is being considered See also BaseInteractable::getOrientation.

Parameters
[in]orientationReference to Vec3D storing the orientation of the particle.
200 {
202 }
void setOrientationViaNormal(Vec3D normal)
Definition: Quaternion.cc:536

References orientation_, and Quaternion::setOrientationViaNormal().

Referenced by AxisymmetricIntersectionOfWalls::AxisymmetricIntersectionOfWalls(), HorizontalBaseScrew::HorizontalBaseScrew(), main(), Screw::Screw(), ScrewsymmetricIntersectionOfWalls::ScrewsymmetricIntersectionOfWalls(), AxisymmetricIntersectionOfWalls::setAxis(), HorizontalBaseScrew::setAxis(), ScrewsymmetricIntersectionOfWalls::setAxis(), SimpleDrumSuperquadrics::setAxis(), InfiniteWall::setNormal(), MinimalExampleDrum::setupInitialConditions(), Silo::setupInitialConditions(), BouncingSuperQuadric::setupInitialConditions(), EllipsoidsBouncingOnWallDemo::setupInitialConditions(), EllipticalSuperQuadricCollision::setupInitialConditions(), GranularCollapse::setupInitialConditions(), ContactDetectionNormalSpheresTest::setupInitialConditions(), ContactDetectionRotatedSpheresTest::setupInitialConditions(), MovingWalls::setupInitialConditions(), ContactDetectionWithWallTester::setupParticleAndWall(), ContactDetectionTester::setupParticles(), ShapeGradientHessianTester::testEllipsoid(), ContactDetectionTester::testEllipsoidsContact(), ContactDetectionWithWallTester::testEllipsoidsContact(), ContactDetectionTester::testSpheresContact(), and ContactDetectionWithWallTester::testSpheresContact().

◆ setPosition()

virtual void BaseInteractable::setPosition ( const Vec3D position)
inlinevirtual

Sets the position of this BaseInteractable.

Interpretation depends on which interactable is being considered See also BaseInteractable::getPosistion.

Parameters
[in]positionReference to Vec3D storing the position of the particle.

Reimplemented in TriangleWall, and TriangleMeshWall.

240  { position_ = position; }

References position_.

Referenced by DPM::actionsAfterSolve(), SmoothChute::actionsBeforeTimeStep(), Slide::actionsBeforeTimeStep(), Chutebelt::actionsBeforeTimeStep(), Chutebelt::actionsOnRestart(), VerticalMixerAngledBlades::addBlades(), IntersectionOfWalls::addObject(), NautaMixer::addParticles(), NautaMixer::addParticlesAtWall(), BaseWall::addParticlesAtWall(), VerticalMixerAngledBlades::addPrettyBlades(), NautaMixer::addScrew(), applyPrescribedPosition(), AxisymmetricIntersectionOfWalls::AxisymmetricIntersectionOfWalls(), DPMBase::checkParticleForInteractionLocalPeriodic(), BaseCluster::computeInternalStructure(), MercuryLogo::constructTextAsParticles(), MPISphericalParticle::copyDataFromMPIParticleToParticle(), Funnel::create_funnel(), LawinenBox::create_inflow_particle(), ChutePeriodic::create_inflow_particle(), ChuteWithContraction::create_inflow_particle(), Funnel::create_inflow_particle(), AngleOfRepose::create_inflow_particle(), FlowRule::create_inflow_particle(), SilbertPeriodic::create_inflow_particle(), SegregationWithHopper::create_inflow_particle(), Slide::create_rough_wall(), Chute::createBottom(), Chute::createFlowParticle(), ChuteWithWedge::createFlowParticle(), Membrane::createVertexParticles(), DPM::DPM(), ChuteWithWedge::extendBottom(), PeriodicBoundaryHandler::findTargetProcessor(), FluxAndPeriodicBoundarySelfTest::FluxAndPeriodicBoundarySelfTest(), SuperQuadricParticle::getContactPointPlanB(), BasicIntersectionOfWalls::getDistanceAndNormal(), BasicUnionOfWalls::getDistanceAndNormal(), TriangleMeshWall::getDistanceAndNormal(), TriangleMeshWall::getInteractionWith(), BasicIntersectionOfWalls::getVTK(), BasicUnionOfWalls::getVTK(), HorizontalBaseScrew::HorizontalBaseScrew(), InfiniteWall::InfiniteWall(), InitialConditions< SpeciesType >::InitialConditions(), HorizontalMixer::introduceParticlesAtWall(), ContactDetectionIntersectionOfWallsTest::introduceParticlesAtWall(), HorizontalMixer::introduceParticlesInDomain(), HorizontalMixer::introduceSingleParticle(), LawinenBox::LawinenBox(), load(), loadingTest(), main(), MercuryCGSelfTest::MercuryCGSelfTest(), normalAndTangentialLoadingTest(), objectivenessTest(), BaseParticle::oldRead(), InfiniteWall::oldRead(), ParticleBeam::ParticleBeam(), ParticleInclusion::ParticleInclusion(), BaseCluster::particleInsertionSuccessful(), FixedClusterInsertionBoundary::placeParticle(), ChuteInsertionBoundary::placeParticle(), CubeInsertionBoundary::placeParticle(), HopperInsertionBoundary::placeParticle(), RandomClusterInsertionBoundary::placeParticle(), FileReader::read(), DPMBase::readNextDataFile(), regimeForceUnitTest::regimeForceUnitTest(), ScrewsymmetricIntersectionOfWalls::ScrewsymmetricIntersectionOfWalls(), InfiniteWall::set(), Slide::set_Walls(), AxisymmetricWallSelfTest::setGeometry(), SphericalIndenter::setIndenterHeight(), MarbleRun::setParticlePosition(), TriangleMeshWall::setPosition(), TriangleWall::setPosition(), ExtremeOverlapUnitTest::setupInitialConditions(), ExtremeOverlapVolumeUnitTest::setupInitialConditions(), T_protectiveWall::setupInitialConditions(), ClosedCSCWalls::setupInitialConditions(), CSCInit::setupInitialConditions(), CSCWalls::setupInitialConditions(), AxisymmetricHopper::setupInitialConditions(), MercuryLogo::setupInitialConditions(), SmoothChute::setupInitialConditions(), multiParticleT1::setupInitialConditions(), MembraneDemo::setupInitialConditions(), Binary::setupInitialConditions(), my_problem::setupInitialConditions(), Chain::setupInitialConditions(), ForceLawsMPI2Test::setupInitialConditions(), MaserRepeatedOutInMPI2Test::setupInitialConditions(), PeriodicBounaryEnteringMPIDomainTest::setupInitialConditions(), SubcriticalMaserBoundaryTESTMPI2Test::setupInitialConditions(), TwoByTwoMPIDomainMPI4Test::setupInitialConditions(), FreeCooling2DinWallsDemo::setupInitialConditions(), FreeCooling3DDemoProblem::setupInitialConditions(), FreeCooling3DinWallsDemo::setupInitialConditions(), FreeCoolingDemoProblem::setupInitialConditions(), HeaterBoundaryTest::setupInitialConditions(), HourGlass2D::setupInitialConditions(), HourGlass::setupInitialConditions(), MinimalExampleDrum::setupInitialConditions(), ParameterStudy1DDemo::setupInitialConditions(), ParameterStudy2DDemo::setupInitialConditions(), ParameterStudy3DDemo::setupInitialConditions(), ShiftingConstantMassFlowMaserBoundarySelfTest::setupInitialConditions(), ShiftingMaserBoundarySelfTest::setupInitialConditions(), FiveParticles::setupInitialConditions(), Cstatic2d::setupInitialConditions(), Cstatic3D::setupInitialConditions(), statistics_while_running< T >::setupInitialConditions(), Drum::setupInitialConditions(), HertzSelfTest::setupInitialConditions(), MindlinSelfTest::setupInitialConditions(), Penetration::setupInitialConditions(), Silo::setupInitialConditions(), ConstantMassFlowMaserBoundaryMixedSpeciesSelfTest::setupInitialConditions(), CubeDeletionBoundarySelfTest::setupInitialConditions(), DeletionBoundarySelfTest::setupInitialConditions(), LeesEdwardsSelfTest::setupInitialConditions(), clumpTest::setupInitialConditions(), CGBasicSelfTest::setupInitialConditions(), CGHandlerSelfTest::setupInitialConditions(), CGStaticBalanceSelfTest::setupInitialConditions(), NewtonsCradleSelftest::setupInitialConditions(), NewtonsCradleSelfTest::setupInitialConditions(), DPM::setupInitialConditions(), ParticleCreation::setupInitialConditions(), ChargedBondedInteractionSelfTest::setupInitialConditions(), ParticleParticleCollision::setupInitialConditions(), WallParticleCollision::setupInitialConditions(), my_problem_HGRID::setupInitialConditions(), FreeFallInteractionSelfTest::setupInitialConditions(), FreeFallSelfTest::setupInitialConditions(), HertzianBSHPInteractionTwoParticleElasticCollision::setupInitialConditions(), ObliqueImpactSelfTest::setupInitialConditions(), TwoBondedParticleElasticCollision::setupInitialConditions(), CoilSelfTest::setupInitialConditions(), ContactDetectionIntersectionOfWallsTest::setupInitialConditions(), MembraneSelfTest::setupInitialConditions(), RollingOverTriangleWalls::setupInitialConditions(), TriangulatedScrewSelfTest::setupInitialConditions(), TriangulatedStepSelfTest::setupInitialConditions(), TriangulatedStepWallSelfTest::setupInitialConditions(), TriangulatedWallSelfTest::setupInitialConditions(), UnionOfWalls::setupInitialConditions(), DrumRot::setupInitialConditions(), RotatingDrum::setupInitialConditions(), ScalingTestInitialConditionsRelax::setupInitialConditions(), Contact::setupInitialConditions(), BouncingSuperQuadric::setupInitialConditions(), EllipsoidsBouncingOnWallDemo::setupInitialConditions(), EllipticalSuperQuadricCollision::setupInitialConditions(), GranularCollapse::setupInitialConditions(), SlidingSpheresUnitTest::setupInitialConditions(), SphericalSuperQuadricCollision::setupInitialConditions(), ContactDetectionNormalSpheresTest::setupInitialConditions(), ContactDetectionRotatedSpheresTest::setupInitialConditions(), ShapesDemo::setupInitialConditions(), VisualisationTest::setupInitialConditions(), MercuryProblem::setupInitialConditions(), protectiveWall::setupInitialConditions(), Tutorial11::setupInitialConditions(), Tutorial12::setupInitialConditions(), Tutorial1::setupInitialConditions(), Tutorial2::setupInitialConditions(), Tutorial3::setupInitialConditions(), Tutorial4::setupInitialConditions(), Tutorial5::setupInitialConditions(), Tutorial6::setupInitialConditions(), Tutorial7::setupInitialConditions(), Tutorial8::setupInitialConditions(), Tutorial9::setupInitialConditions(), AngledPeriodicBoundarySecondUnitTest::setupInitialConditions(), AngledPeriodicBoundaryUnitTest::setupInitialConditions(), Packing::setupInitialConditions(), CreateDataAndFStatFiles::setupInitialConditions(), ChargedBondedParticleUnitTest::setupInitialConditions(), ExtremeOverlapWithWallsUnitTest::setupInitialConditions(), FreeFallHertzMindlinUnitTest::setupInitialConditions(), FreeFall::setupInitialConditions(), FullRestartTest::setupInitialConditions(), HertzContactRestitutionUnitTest::setupInitialConditions(), HertzianSinterForceUnitTest::setupInitialConditions(), MD_demo::setupInitialConditions(), InclinedPlane::setupInitialConditions(), MpiMaserChuteTest::setupInitialConditions(), MovingIntersectionOfWallsUnitTest_Basic::setupInitialConditions(), MovingWalls::setupInitialConditions(), MovingWall::setupInitialConditions(), MovingWallReference::setupInitialConditions(), MovingWallSimpleIntegration::setupInitialConditions(), MovingWallPrescribedVelocity::setupInitialConditions(), MovingWallTangential::setupInitialConditions(), MovingWallTangentialReference::setupInitialConditions(), MovingWallTangentialSimpleIntegration::setupInitialConditions(), MovingWallTangentialPrescribedVelocity::setupInitialConditions(), MultiParticlesInsertion::setupInitialConditions(), MpiPeriodicBoundaryUnitTest::setupInitialConditions(), PeriodicWallsWithSlidingFrictionUnitTest::setupInitialConditions(), PlasticForceUnitTest::setupInitialConditions(), SeparateFilesSelfTest::setupInitialConditions(), SinterForceUnitTest::setupInitialConditions(), SpeciesTest::setupInitialConditions(), TangentialSpringEnergyConservationUnitTest::setupInitialConditions(), TangentialSpringUnitTest::setupInitialConditions(), WallSpecies::setupInitialConditions(), ChuteBottom::setupInitialConditions(), ScaleCoupledBeam::setupMercury(), CoupledBeam::setupMercury(), CoupledProblem::setupMercury(), ContactDetectionWithWallTester::setupParticleAndWall(), ContactDetectionTester::setupParticles(), TriangleWall::setVertices(), MeshTriangle::setVertices(), SinterPair::SinterPair(), Slide::Slide(), FlowFrontChute::stretch(), ContactDetectionTester::testEllipsoidsContact(), ContactDetectionWithWallTester::testEllipsoidsContact(), ContactDetectionTester::testSpheresContact(), ContactDetectionWithWallTester::testSpheresContact(), PeriodicBoundaryHandler::updateParticles(), ClumpParticle::updatePebblesVelPos(), and viscoElasticUnitTest::viscoElasticUnitTest().

◆ setPrescribedAngularVelocity()

void BaseInteractable::setPrescribedAngularVelocity ( const std::function< Vec3D(double)> &  prescribedAngularVelocity)

Allows the angular velocity of the infinite mass interactable to be prescribed.

502 {
503  prescribedAngularVelocity_ = prescribedAngularVelocity;
504 }

References prescribedAngularVelocity_.

◆ setPrescribedOrientation()

void BaseInteractable::setPrescribedOrientation ( const std::function< Quaternion(double)> &  prescribedOrientation)

Allows the orientation of the infinite mass interactbale to be prescribed.

This is similar to the BaseInteractable::setPrescribedPosition and
works the same way. See BaseInteractable::setPrescibedPosition and BaseInteractable::setPrescribedVelocity for more details how it works. Note, the rate of change of the orientation can also be set using the function BaseInteractable::setPrescribedAngularVelocity.

Parameters
[in]prescribedOrientationstd::function which is the lambda function and takes a double the time and returns a Vec 3D which is the orientation of the interactable for that time.
475 {
476  prescribedOrientation_ = prescribedOrientation;
477 }

References prescribedOrientation_.

◆ setPrescribedPosition()

void BaseInteractable::setPrescribedPosition ( const std::function< Vec3D(double)> &  prescribedPosition)

Allows the position of an infinite mass interactable to be prescribed.

This functions is used to give an interactable a prescribed motion, which is defined by a std::function. This is the new moving walls interface. A demo of the use would be: setPrescribedPosition([this] (double time) { return Vec3D(getXMin(),0.0,shaker_amp * std::sin(t * 2.0 * shaker_freq * constants::pi)); } ); This example moves the wall sinusoidally with time.

Parameters
[in]prescribedPositionstd::function which is the lambda function and takes a double the time and returns a Vec 3D which is the position of the interactable for that time. See also BaseInteractable::setPrescribedVelocity for more information.
414 {
415  prescribedPosition_ = prescribedPosition;
416 }

References prescribedPosition_.

Referenced by main(), Binary::setupInitialConditions(), my_problem::setupInitialConditions(), DrivenParticleClass::setupInitialConditions(), ExtremeOverlapWithWallsUnitTest::setupInitialConditions(), MovingWallPrescribedPosition::setupInitialConditions(), MovingWallPrescribedPositionPrescribedVelocity::setupInitialConditions(), MovingWallTangentialPrescribedPosition::setupInitialConditions(), MovingWallTangentialPrescribedPositionPrescribedVelocity::setupInitialConditions(), and SCoupling< M, O >::updateTriangleWall().

◆ setPrescribedVelocity()

void BaseInteractable::setPrescribedVelocity ( const std::function< Vec3D(double)> &  prescribedVelocity)

Allows the velocity of an infinite mass interactable to be prescribed.

In a similar manner to BaseInteractable::setPrescribedPosition this sets the velocity of the BaseInteactable. See also BaseInteractable::setPrescribedPosition Note, it is valid to set both the velocity and the position. No checking if these are consist is done at all. If you only set one of these function the other is automatically calculated using a by numerically differentiating or integrating the functions will is prescribed.

Parameters
[in]prescribedVelocitystd::function which is the lambda function and takes a double the time and returns a Vec 3D which is the velocity of the intertable for that time.
445 {
446  prescribedVelocity_ = prescribedVelocity;
447 }

References prescribedVelocity_.

Referenced by ParticleBeam::ParticleBeam(), BaseWall::setForceControl(), MovingWallPrescribedVelocity::setupInitialConditions(), MovingWallPrescribedPositionPrescribedVelocity::setupInitialConditions(), MovingWallTangentialPrescribedVelocity::setupInitialConditions(), MovingWallTangentialPrescribedPositionPrescribedVelocity::setupInitialConditions(), BaseWall::setVelocityControl(), and SCoupling< M, O >::updateTriangleWall().

◆ setSpecies()

void BaseInteractable::setSpecies ( const ParticleSpecies species)

Sets the species of this BaseInteractable.

This function sets the species associated with this interactable object. Again this must be a ParticleSpecies. Note, it also automatically sets the index on the indSpecies_ by working up the correct index. However, index should be carefully used

Parameters
[in]speciesParticleSpcies pointer which is species holding the physical properties.
186 {
187  logger.assert_debug(species->getHandler() != nullptr && species->getHandler()->getObject(species->getIndex())==species, "Error: Species is not part of any handler yet");
188  species_ = species;
189  indSpecies_ = species->getIndex();
190 }
T * getObject(const unsigned int id)
Gets a pointer to the Object at the specified index in the BaseHandler.
Definition: BaseHandler.h:613
SpeciesHandler * getHandler() const
Returns the pointer to the handler to which this species belongs.
Definition: BaseSpecies.cc:99

References BaseSpecies::getHandler(), BaseObject::getIndex(), BaseHandler< T >::getObject(), indSpecies_, logger, and species_.

Referenced by BaseParticle::setSpecies(), and BaseWall::setSpecies().

◆ setTorque()

void BaseInteractable::setTorque ( const Vec3D torque)
inline

Sets the torque on this BaseInteractable.

This sets the torque being applied to this interactable. Note, first the code computes all force/torques in the interactions and then loops over all interactable objects applying the torques from the interactions to the interactables involved in the interaction.

Parameters
[in]torqueVec3D which is the force to be applied.
162  { torque_ = torque; }

References torque_.

Referenced by SmoothChute::actionsBeforeTimeStep().

◆ setVelocity()

void BaseInteractable::setVelocity ( const Vec3D velocity)

set the velocity of the BaseInteractable.

See also BaseInteractable::getVelocity

Parameters
[in]velocityVec3D which is the velocity of the interactable.
351 {
352  velocity_ = velocity;
353 }

References velocity_.

Referenced by ClosedCSCRestart::actionsAfterTimeStep(), ClosedCSCRun::actionsAfterTimeStep(), ClosedCSCWalls::actionsAfterTimeStep(), SmoothChute::actionsBeforeTimeStep(), Slide::actionsBeforeTimeStep(), Chutebelt::actionsBeforeTimeStep(), Chutebelt::actionsOnRestart(), NautaMixer::addParticles(), NautaMixer::addParticlesAtWall(), applyPrescribedVelocity(), HeaterBoundary::checkBoundaryAfterParticleMoved(), ClosedCSCRun::ClosedCSCRun(), BaseCluster::computeInternalStructure(), MPISphericalParticle::copyDataFromMPIParticleToParticle(), LawinenBox::create_inflow_particle(), ChutePeriodic::create_inflow_particle(), ChuteWithContraction::create_inflow_particle(), Funnel::create_inflow_particle(), FlowRule::create_inflow_particle(), SilbertPeriodic::create_inflow_particle(), Chute::createBottom(), Chute::createFlowParticle(), ChuteWithWedge::createFlowParticle(), BaseParticle::fixParticle(), FluxAndPeriodicBoundarySelfTest::FluxAndPeriodicBoundarySelfTest(), integrateAfterForceComputation(), integrateBeforeForceComputation(), HorizontalMixer::introduceParticlesAtWall(), HorizontalMixer::introduceParticlesInDomain(), HorizontalMixer::introduceSingleParticle(), load(), loadingTest(), main(), normalAndTangentialLoadingTest(), objectivenessTest(), BaseParticle::oldRead(), InfiniteWall::oldRead(), InfiniteWallWithHole::oldRead(), BaseCluster::particleInsertionSuccessful(), FixedClusterInsertionBoundary::placeParticle(), ChuteInsertionBoundary::placeParticle(), CubeInsertionBoundary::placeParticle(), HopperInsertionBoundary::placeParticle(), RandomClusterInsertionBoundary::placeParticle(), FileReader::read(), DPMBase::readNextDataFile(), WallHandler::readTriangleWall(), ClosedCSCWalls::saveWalls(), Slide::set_Walls(), SphericalIndenter::setIndenterVelocity(), ExtremeOverlapUnitTest::setupInitialConditions(), MercuryLogo::setupInitialConditions(), SmoothChute::setupInitialConditions(), multiParticleT1::setupInitialConditions(), Binary::setupInitialConditions(), my_problem::setupInitialConditions(), ForceLawsMPI2Test::setupInitialConditions(), MaserRepeatedOutInMPI2Test::setupInitialConditions(), PeriodicBounaryEnteringMPIDomainTest::setupInitialConditions(), SubcriticalMaserBoundaryTESTMPI2Test::setupInitialConditions(), TwoByTwoMPIDomainMPI4Test::setupInitialConditions(), FreeCooling2DinWalls::setupInitialConditions(), FreeCooling2DinWallsDemo::setupInitialConditions(), FreeCooling3DDemoProblem::setupInitialConditions(), FreeCooling3DinWallsDemo::setupInitialConditions(), FreeCoolingDemoProblem::setupInitialConditions(), HeaterBoundaryTest::setupInitialConditions(), HourGlass2D::setupInitialConditions(), HourGlass::setupInitialConditions(), MinimalExampleDrum::setupInitialConditions(), ParameterStudy1DDemo::setupInitialConditions(), ParameterStudy2DDemo::setupInitialConditions(), ParameterStudy3DDemo::setupInitialConditions(), ShiftingConstantMassFlowMaserBoundarySelfTest::setupInitialConditions(), ShiftingMaserBoundarySelfTest::setupInitialConditions(), FiveParticles::setupInitialConditions(), Cstatic2d::setupInitialConditions(), Cstatic3D::setupInitialConditions(), statistics_while_running< T >::setupInitialConditions(), HertzSelfTest::setupInitialConditions(), Penetration::setupInitialConditions(), ConstantMassFlowMaserBoundaryMixedSpeciesSelfTest::setupInitialConditions(), clumpTest::setupInitialConditions(), CGBasicSelfTest::setupInitialConditions(), CGHandlerSelfTest::setupInitialConditions(), CGStaticBalanceSelfTest::setupInitialConditions(), NewtonsCradleSelftest::setupInitialConditions(), DPM::setupInitialConditions(), ChargedBondedInteractionSelfTest::setupInitialConditions(), ParticleParticleCollision::setupInitialConditions(), WallParticleCollision::setupInitialConditions(), my_problem_HGRID::setupInitialConditions(), FreeFallInteractionSelfTest::setupInitialConditions(), FreeFallSelfTest::setupInitialConditions(), HertzianBSHPInteractionTwoParticleElasticCollision::setupInitialConditions(), ObliqueImpactSelfTest::setupInitialConditions(), TwoBondedParticleElasticCollision::setupInitialConditions(), CoilSelfTest::setupInitialConditions(), RollingOverTriangleWalls::setupInitialConditions(), TriangulatedScrewSelfTest::setupInitialConditions(), TriangulatedStepSelfTest::setupInitialConditions(), UnionOfWalls::setupInitialConditions(), DrumRot::setupInitialConditions(), RotatingDrum::setupInitialConditions(), ScalingTestInitialConditionsRelax::setupInitialConditions(), Contact::setupInitialConditions(), Wall::setupInitialConditions(), BouncingSuperQuadric::setupInitialConditions(), EllipsoidsBouncingOnWallDemo::setupInitialConditions(), EllipticalSuperQuadricCollision::setupInitialConditions(), GranularCollapse::setupInitialConditions(), SlidingSpheresUnitTest::setupInitialConditions(), SphericalSuperQuadricCollision::setupInitialConditions(), ContactDetectionNormalSpheresTest::setupInitialConditions(), ContactDetectionRotatedSpheresTest::setupInitialConditions(), VisualisationTest::setupInitialConditions(), MercuryProblem::setupInitialConditions(), Tutorial11::setupInitialConditions(), Tutorial12::setupInitialConditions(), Tutorial1::setupInitialConditions(), Tutorial2::setupInitialConditions(), Tutorial3::setupInitialConditions(), Tutorial4::setupInitialConditions(), Tutorial5::setupInitialConditions(), Tutorial6::setupInitialConditions(), Tutorial7::setupInitialConditions(), Tutorial8::setupInitialConditions(), Tutorial9::setupInitialConditions(), AngledPeriodicBoundarySecondUnitTest::setupInitialConditions(), CreateDataAndFStatFiles::setupInitialConditions(), ChargedBondedParticleUnitTest::setupInitialConditions(), ExtremeOverlapWithWallsUnitTest::setupInitialConditions(), FreeFallHertzMindlinUnitTest::setupInitialConditions(), HertzContactRestitutionUnitTest::setupInitialConditions(), HertzianSinterForceUnitTest::setupInitialConditions(), MD_demo::setupInitialConditions(), InclinedPlane::setupInitialConditions(), MpiMaserChuteTest::setupInitialConditions(), MovingIntersectionOfWallsUnitTest_Basic::setupInitialConditions(), MovingIntersectionOfWallsUnitTest_MovingReferenceFrame::setupInitialConditions(), MovingWall::setupInitialConditions(), MovingWallReference::setupInitialConditions(), MovingWallSimpleIntegration::setupInitialConditions(), MovingWallTangential::setupInitialConditions(), MovingWallTangentialReference::setupInitialConditions(), MovingWallTangentialSimpleIntegration::setupInitialConditions(), MultiParticlesInsertion::setupInitialConditions(), MpiPeriodicBoundaryUnitTest::setupInitialConditions(), PeriodicWallsWithSlidingFrictionUnitTest::setupInitialConditions(), PlasticForceUnitTest::setupInitialConditions(), SeparateFilesSelfTest::setupInitialConditions(), SinterForceUnitTest::setupInitialConditions(), SpeciesTest::setupInitialConditions(), TangentialSpringEnergyConservationUnitTest::setupInitialConditions(), TangentialSpringUnitTest::setupInitialConditions(), WallSpecies::setupInitialConditions(), ChuteBottom::setupInitialConditions(), MembraneDemo::setUpMembrane(), MembraneSelfTest::setUpMembrane(), ScaleCoupledBeam::setupMercury(), CoupledBeam::setupMercury(), ContactDetectionWithWallTester::setupParticleAndWall(), ContactDetectionTester::setupParticles(), Slide::Slide(), ContactDetectionTester::testSpheresContact(), PeriodicBoundaryHandler::updateParticles(), ClumpParticle::updatePebblesVelPos(), and viscoElasticUnitTest::viscoElasticUnitTest().

◆ sumForceTorqueOMP()

void BaseInteractable::sumForceTorqueOMP ( )
163 {
164  #ifdef MERCURYDPM_USE_OMP
165  for (const auto force : forceOMP_)
166  {
167  force_ += force;
168  }
169 
170  for (const auto torque : torqueOMP_)
171  {
172  torque_ += torque;
173  }
174  #endif
175 }

References force_, forceOMP_, torque_, and torqueOMP_.

Referenced by SphericalIndenter::actionsAfterTimeStep(), and DPMBase::computeAllForces().

◆ write()

void BaseInteractable::write ( std::ostream &  os) const
overridevirtual

Write a BaseInteractable to an output stream.

Parameters
[in]osThe output stream to which the BaseInteractable is written.

BaseInteractable write function. Writes out all the information required to recreate this interactable. To write this interactable to the screen call write(std::cout). See also BaseInteractable::read

Parameters
[in]osstd::ostream to which the information is written. Note, is any ostream is can be file or screen.
Todo:
take the zero out

Implements BaseObject.

Reimplemented in WearableTriangulatedWall, WearableTriangleMeshWall, WearableNurbsWall, VChute, TriangulatedWall, TriangleWall, TriangleMeshWall, SphericalWall, SineWall, Screw, RestrictedWall, ParabolaChute, NurbsWall, MeshTriangle, IntersectionOfWalls, InfiniteWallWithHole, HorizontalScrew, CylindricalWall, Combtooth, Coil, BasicUnionOfWalls, BasicIntersectionOfWalls, BaseWall, ArcWall, SuperQuadricParticle, ClumpParticle, BaseParticle, SimpleDrumSuperquadrics, ScrewsymmetricIntersectionOfWalls, HorizontalBaseScrew, and AxisymmetricIntersectionOfWalls.

275 {
276  BaseObject::write(os);
277  os << " indSpecies " << indSpecies_
278  << " position " << position_
279  << " orientation " << orientation_
280  //<< " axis " << orientation_.getAxis()
281  << " velocity " << velocity_
282  << " angularVelocity " << angularVelocity_
283  << " force " << force_
284  << " torque " << torque_;
285 }
virtual void write(std::ostream &os) const =0
A purely virtual function which has an implementation which writes the name and the object id_ to the...
Definition: BaseObject.cc:91

References angularVelocity_, force_, indSpecies_, orientation_, position_, torque_, velocity_, and BaseObject::write().

Referenced by BaseParticle::write(), and BaseWall::write().

Member Data Documentation

◆ angularVelocity_

Vec3D BaseInteractable::angularVelocity_
private

Store the angular velocity of the interactable.

Referenced by addAngularVelocity(), BaseInteractable(), getAngularVelocity(), read(), setAngularVelocity(), and write().

◆ force_

Vec3D BaseInteractable::force_
private

Stores the force applied to the interactable.

Referenced by addForce(), BaseInteractable(), getForce(), read(), resetForceTorque(), setForce(), sumForceTorqueOMP(), and write().

◆ forceOMP_

std::vector<Vec3D> BaseInteractable::forceOMP_
private

◆ indSpecies_

unsigned int BaseInteractable::indSpecies_
private

Stores the index on the species associated with this interactable.

Referenced by BaseInteractable(), getIndSpecies(), read(), setIndSpecies(), setSpecies(), and write().

◆ interactions_

std::vector<BaseInteraction*> BaseInteractable::interactions_
private

List of interactions this interactable is involved with.

Referenced by addInteraction(), BaseInteractable(), copyInteractionsForPeriodicParticles(), getInteractions(), removeInteraction(), and ~BaseInteractable().

◆ orientation_

Quaternion BaseInteractable::orientation_
private

Stores the orientation of the interactable. Exactly what is stored depends on the type of interatable

Referenced by BaseInteractable(), getOrientation(), read(), rotate(), setOrientation(), setOrientationViaEuler(), setOrientationViaNormal(), and write().

◆ position_

Vec3D BaseInteractable::position_
private

Stores the position of the interactable. Exactly what is stored depends on the type of interactable.

Referenced by BaseInteractable(), getPosition(), move(), read(), setPosition(), and write().

◆ prescribedAngularVelocity_

std::function<Vec3D(double)> BaseInteractable::prescribedAngularVelocity_
private

User defined functions which if set describes the angular velocity of the interactable.

Referenced by applyPrescribedAngularVelocity(), BaseInteractable(), integrateAfterForceComputation(), integrateBeforeForceComputation(), and setPrescribedAngularVelocity().

◆ prescribedOrientation_

std::function<Quaternion(double)> BaseInteractable::prescribedOrientation_
private

User defined function which if set describes the orientation of the interactable

Referenced by applyPrescribedOrientation(), BaseInteractable(), integrateAfterForceComputation(), integrateBeforeForceComputation(), and setPrescribedOrientation().

◆ prescribedPosition_

std::function<Vec3D(double)> BaseInteractable::prescribedPosition_
private

User defined function which if set describes the position of the interactable

Referenced by applyPrescribedPosition(), BaseInteractable(), integrateAfterForceComputation(), integrateBeforeForceComputation(), and setPrescribedPosition().

◆ prescribedVelocity_

std::function<Vec3D(double)> BaseInteractable::prescribedVelocity_
private

User defined function which if set describes the velocity of the interactable.

Referenced by applyPrescribedVelocity(), BaseInteractable(), integrateAfterForceComputation(), integrateBeforeForceComputation(), and setPrescribedVelocity().

◆ species_

const ParticleSpecies* BaseInteractable::species_
private

Point to the ParticlesSpecies which stores density and other material properties of the interactable.

Referenced by BaseInteractable(), getSpecies(), and setSpecies().

◆ torque_

Vec3D BaseInteractable::torque_
private

Stores the torque applied to the interactable.

Referenced by addTorque(), BaseInteractable(), getTorque(), read(), resetForceTorque(), setTorque(), sumForceTorqueOMP(), and write().

◆ torqueOMP_

std::vector<Vec3D> BaseInteractable::torqueOMP_
private

See BaseInteractable::forceOMP_.

Referenced by addTorque(), resetForceTorque(), and sumForceTorqueOMP().

◆ velocity_

Vec3D BaseInteractable::velocity_
private

Stores the velocity of this interactable.

Referenced by addVelocity(), BaseInteractable(), getVelocity(), read(), setVelocity(), and write().


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