InfiniteWall Class Referencefinal

A infinite wall fills the half-space {point: (position_-point)*normal_<=0}. More...

#include <InfiniteWall.h>

+ Inheritance diagram for InfiniteWall:

Public Member Functions

 InfiniteWall ()
 Default constructor, the normal is infinitely long. More...
 
 InfiniteWall (const InfiniteWall &w)
 Copy constructor, copy the given wall. More...
 
 InfiniteWall (const ParticleSpecies *species)
 Constructor setting species. More...
 
 InfiniteWall (Vec3D normal, Vec3D point, const ParticleSpecies *species)
 Constructor setting values. More...
 
 InfiniteWall (Vec3D PointA, Vec3D PointB, Vec3D PointC, const ParticleSpecies *species)
 Constructor setting values if 3 coordinates are given. More...
 
 ~InfiniteWall () override
 Default destructor. More...
 
InfiniteWallcopy () const override
 Wall copy method. It calls the copy constructor of this Wall, useful for polymorphism. More...
 
void set (Vec3D normal, Vec3D point)
 Defines a standard wall, given an outward normal vector s.t. normal*x=normal*point for all x of the wall. More...
 
void setNormal (Vec3D normal)
 Changes the normal of the InfiniteWall. More...
 
MERCURYDPM_DEPRECATED void set (Vec3D normal, Mdouble position)
 Defines a standard wall by computing normal*position = point and using the overloaded function set(Vec3D, vec3D). More...
 
Mdouble getDistance (Vec3D position) const
 Returns the distance of the wall to the particle. More...
 
bool getDistanceAndNormal (const BaseParticle &p, Mdouble &distance, Vec3D &normal_return) const override
 Compute the distance from the wall for a given BaseParticle and return if there is a collision. If there is a collision, also return the normal vector. More...
 
bool getDistanceNormalOverlapSuperquadric (const SuperQuadricParticle &p, Mdouble &distance, Vec3D &normal_return, Mdouble &overlap) const override
 Compute the distance from the wall for a given BaseParticle and return if there is a collision. If there is a collision, also return the normal vector. More...
 
void read (std::istream &is) override
 Reads InfiniteWall from a restart file. More...
 
void oldRead (std::istream &is)
 Reads InfiniteWall from an old-style restart file. More...
 
std::string getName () const override
 Writes the InfiniteWall to an output stream, usually a restart file. More...
 
Vec3D getNormal () const
 Access function for normal. More...
 
void createVTK (std::vector< Vec3D > &myPoints) const
 
void createVTK (std::vector< Vec3D > &myPoints, Vec3D max, Vec3D min) const
 
void writeVTK (VTKContainer &vtk) const override
 
Vec3D getFurthestPointSuperQuadric (const Vec3D &normalBodyFixed, const Vec3D &axes, Mdouble eps1, Mdouble eps2) const override
 Largely untested, use at your own risk for anything other than ellipsoids. More...
 
- Public Member Functions inherited from BaseWall
 BaseWall ()
 Default constructor. More...
 
 BaseWall (const BaseWall &w)
 Copy constructor. More...
 
 ~BaseWall () override
 Default destructor. More...
 
void write (std::ostream &os) const override
 Function that writes a BaseWall to an output stream, usually a restart file. More...
 
virtual bool getDistanceNormalOverlap (const BaseParticle &P, Mdouble &distance, Vec3D &normal_return, Mdouble &overlap) const
 
virtual void setHandler (WallHandler *handler)
 A function which sets the WallHandler for this BaseWall. More...
 
WallHandlergetHandler () const
 A function which returns the WallHandler that handles this BaseWall. More...
 
void setIndSpecies (unsigned int indSpecies) override
 Define the species of this wall using the index of the species in the SpeciesHandler in this DPMBase. More...
 
void setSpecies (const ParticleSpecies *species)
 Defines the species of the current wall. More...
 
bool isFixed () const override
 
void setForceControl (Vec3D forceGoal, Vec3D gainFactor, Vec3D baseVelocity={0, 0, 0})
 Slowly adjusts the force on a wall towards a specified goal, by adjusting (prescribing) the velocity of the wall. More...
 
virtual bool isLocal (Vec3D &min, Vec3D &max) const
 
bool getLinePlaneIntersect (Vec3D &intersect, const Vec3D &p0, const Vec3D &p1, const Vec3D &n, const Vec3D &p)
 
bool isInsideWallVTK (const Vec3D &point, const Vec3D &normal, const Vec3D &position) const
 
void projectOntoWallVTK (Vec3D &point0, const Vec3D &point1, const Vec3D &normal, const Vec3D &position) const
 
void intersectVTK (std::vector< Vec3D > &points, Vec3D normal, Vec3D position) const
 
virtual BaseInteractiongetInteractionWithSuperQuad (SuperQuadricParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler)
 
void getVTK (std::vector< Vec3D > &points, std::vector< std::vector< double >> &triangleStrips)
 
const Vec3D getAxis () const
 
BaseInteractiongetInteractionWith (BaseParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler) override
 Returns the interaction between this wall and a given particle, nullptr if there is no interaction. More...
 
virtual void actionsOnRestart ()
 No implementation but can be overidden in its derived classes. More...
 
virtual void actionsAfterParticleGhostUpdate ()
 No implementation but can be overidden in its derived classes. More...
 
virtual void handleParticleAddition (unsigned int id, BaseParticle *p)
 Handles the addition of particles to the particleHandler. More...
 
virtual void handleParticleRemoval (unsigned int id)
 Handles the addition of particles to the particleHandler. More...
 
virtual void checkInteractions (InteractionHandler *interactionHandler, unsigned int timeStamp)
 Check if all interactions are valid. More...
 
bool getVTKVisibility () const
 
void setVTKVisibility (bool vtkVisibility)
 
void addRenderedWall (BaseWall *w)
 
BaseWallgetRenderedWall (size_t i) const
 
std::vector< BaseWall * > getRenderedWalls () const
 
void removeRenderedWalls ()
 
void renderWall (VTKContainer &vtk)
 
void addParticlesAtWall (unsigned numElements=50)
 
void setVelocityControl (Vec3D forceGoal, Vec3D gainFactor, Vec3D baseVelocity)
 
virtual void writeWallDetailsVTK (VTKData &data) const
 
virtual void computeWear ()
 
- Public Member Functions inherited from BaseInteractable
 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...
 
unsigned int getIndSpecies () const
 Returns the index of the species associated with the interactable object. 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 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 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 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
 

Additional Inherited Members

- Static Public Member Functions inherited from BaseWall
static void addToVTK (const std::vector< Vec3D > &points, VTKContainer &vtk)
 Takes the points provided and adds a triangle strip connecting these points to the vtk container. More...
 

Detailed Description

A infinite wall fills the half-space {point: (position_-point)*normal_<=0}.

This is a class defining walls. It defines the interaction of regular walls and periodic walls with particles as defined in Particle Modifications:

Thus, the surface of the wall is a plane through position position_ with normal_ the outward unit normal vector of the wall (pointing away from the particles, into the wall). Please note that this wall is infinite and straight.

A particle touches an infinite wall if (position_-point)*normal_<=radius.

Constructor & Destructor Documentation

◆ InfiniteWall() [1/5]

InfiniteWall::InfiniteWall ( )

Default constructor, the normal is infinitely long.

35 {
36  logger(DEBUG, "InfiniteWall::InfiniteWall ) finished");
37 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ DEBUG

References DEBUG, and logger.

Referenced by copy().

◆ InfiniteWall() [2/5]

InfiniteWall::InfiniteWall ( const InfiniteWall w)

Copy constructor, copy the given wall.

Parameters
[in]wInfiniteWall that has to be copied.

First copy the attributes of the BaseWall, then copy the ones that are specific for the InfiniteWall.

45  : BaseWall(w)
46 {
47  logger(DEBUG, "InfiniteWall::InfiniteWall(const InfiniteWall &p) finished");
48 }
BaseWall()
Default constructor.
Definition: BaseWall.cc:36

References DEBUG, and logger.

◆ InfiniteWall() [3/5]

InfiniteWall::InfiniteWall ( const ParticleSpecies species)
explicit

Constructor setting species.

51  : BaseWall()
52 {
53  setSpecies(s);
54  logger(DEBUG, "InfiniteWall::InfiniteWall(const ParticleSpecies* s) finished");
55 }
void setSpecies(const ParticleSpecies *species)
Defines the species of the current wall.
Definition: BaseWall.cc:169

References DEBUG, logger, and BaseWall::setSpecies().

◆ InfiniteWall() [4/5]

InfiniteWall::InfiniteWall ( Vec3D  normal,
Vec3D  point,
const ParticleSpecies species 
)

Constructor setting values.

58 {
59  setNormal(normal);
60  setPosition(point);
61  setSpecies(species);
62 }
virtual void setPosition(const Vec3D &position)
Sets the position of this BaseInteractable.
Definition: BaseInteractable.h:239
void setNormal(Vec3D normal)
Changes the normal of the InfiniteWall.
Definition: InfiniteWall.cc:127

References setNormal(), BaseInteractable::setPosition(), and BaseWall::setSpecies().

◆ InfiniteWall() [5/5]

InfiniteWall::InfiniteWall ( Vec3D  PointA,
Vec3D  PointB,
Vec3D  PointC,
const ParticleSpecies species 
)

Constructor setting values if 3 coordinates are given.

Parameters
PointAfirst coordinate
PointBsecond coordinate
PointCthird coordinate
species

Builds an infinite wall through 3 points, normal is defined with Right hand rule following the three input points.

73 {
74 //function returns infinite wall passing through 3 points
75 //subtract second and third from the first
76  Vec3D SubtB = PointA - PointB;
77  Vec3D SubtC = PointA - PointC;
78 
79  Vec3D WallNormal = Vec3D::cross(SubtB, SubtC);
80  //Check if walls coordinates inline, if true Do not build wall and give error message.
81 
82  if (WallNormal.getLengthSquared() == 0.0)
83  {
84  logger(ERROR,
85  "Error Building InfiniteWall out of 3 coordinates. Coordinates are in line, Wall not constructed.");
86  }
87  else
88  {
89  setNormal(WallNormal);
90  setPosition(PointA);
91  setSpecies(species);
92  }
93 
94 
95 }
@ ERROR
Definition: Vector.h:51
static Mdouble getLengthSquared(const Vec3D &a)
Calculates the squared length of a Vec3D: .
Definition: Vector.h:332
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:163

References Vec3D::cross(), ERROR, Vec3D::getLengthSquared(), logger, setNormal(), BaseInteractable::setPosition(), and BaseWall::setSpecies().

◆ ~InfiniteWall()

InfiniteWall::~InfiniteWall ( )
override

Default destructor.

98 {
99  logger(DEBUG, "InfiniteWall::~InfiniteWall finished");
100 }

References DEBUG, and logger.

Member Function Documentation

◆ copy()

InfiniteWall * InfiniteWall::copy ( ) const
overridevirtual

Wall copy method. It calls the copy constructor of this Wall, useful for polymorphism.

Wall copy method. It calls the copy constructor of this Wall, useful for polymorphism

Implements BaseWall.

106 {
107  return new InfiniteWall(*this);
108 }
InfiniteWall()
Default constructor, the normal is infinitely long.
Definition: InfiniteWall.cc:34

References InfiniteWall().

Referenced by RestrictedWall::set().

◆ createVTK() [1/2]

void InfiniteWall::createVTK ( std::vector< Vec3D > &  myPoints) const

Returns all intersection points of the infinite wall with the domain boundary.

Used to create an array of points for representing the wall in the VTK viewer. E.g., for a InfiniteWall through p=(0,0,0) with normal n=(0,0,-1) in a domain (-1,1)^3, createVTK returns {(0,-1,-1), (0,-1,1), (0,1,1), (0,1,-1)} Calling addToVTK will then create a triangle strip connecting these points with triangle faces.

225 {
226  Vec3D max = getHandler()->getDPMBase()->getMax();
227  Vec3D min = getHandler()->getDPMBase()->getMin();
228  createVTK(myPoints, min, max);
229 }
DPMBase * getDPMBase()
Gets the problem that is solved using this handler.
Definition: BaseHandler.h:725
WallHandler * getHandler() const
A function which returns the WallHandler that handles this BaseWall.
Definition: BaseWall.cc:134
Vec3D getMax() const
Definition: DPMBase.h:670
Vec3D getMin() const
Definition: DPMBase.h:664
void createVTK(std::vector< Vec3D > &myPoints) const
Definition: InfiniteWall.cc:224

References BaseHandler< T >::getDPMBase(), BaseWall::getHandler(), DPMBase::getMax(), and DPMBase::getMin().

Referenced by LevelSetWall::createVTK(), and writeVTK().

◆ createVTK() [2/2]

void InfiniteWall::createVTK ( std::vector< Vec3D > &  myPoints,
Vec3D  max,
Vec3D  min 
) const

Same as createVTK(), but with a self-defined domain size (useful for plotting AxisymmetricWall's).

232 {
233  const Vec3D& n = getOrientation().getAxis();
234  const Vec3D& p = getPosition();
235 
236  if (fabs(n.X) > 0.5)
237  {
238  // If the wall normal has a nonzero x-component,
239  // We first find four intersection points with the four domain edges pointing in x-direction.
240  // Because these points might be outside the domain in x-direction, we use the intersection function have points in the domain
241  myPoints.emplace_back(p.X - ((min.Y - p.Y) * n.Y + (min.Z - p.Z) * n.Z) / n.X, min.Y, min.Z);
242  myPoints.emplace_back(p.X - ((min.Y - p.Y) * n.Y + (max.Z - p.Z) * n.Z) / n.X, min.Y, max.Z);
243  myPoints.emplace_back(p.X - ((max.Y - p.Y) * n.Y + (max.Z - p.Z) * n.Z) / n.X, max.Y, max.Z);
244  myPoints.emplace_back(p.X - ((max.Y - p.Y) * n.Y + (min.Z - p.Z) * n.Z) / n.X, max.Y, min.Z);
245  intersectVTK(myPoints, Vec3D(+1, 0, 0), Vec3D(max.X, 0, 0));
246  intersectVTK(myPoints, Vec3D(-1, 0, 0), Vec3D(min.X, 0, 0));
247  }
248  else if (fabs(n.Y) > 0.5)
249  {
250  // Else, if the wall normal has a nonzero y-component ...
251  myPoints.emplace_back(min.X, p.Y - ((min.X - p.X) * n.X + (min.Z - p.Z) * n.Z) / n.Y, min.Z);
252  myPoints.emplace_back(min.X, p.Y - ((min.X - p.X) * n.X + (max.Z - p.Z) * n.Z) / n.Y, max.Z);
253  myPoints.emplace_back(max.X, p.Y - ((max.X - p.X) * n.X + (max.Z - p.Z) * n.Z) / n.Y, max.Z);
254  myPoints.emplace_back(max.X, p.Y - ((max.X - p.X) * n.X + (min.Z - p.Z) * n.Z) / n.Y, min.Z);
255  intersectVTK(myPoints, Vec3D(0, +1, 0), Vec3D(0, max.Y, 0));
256  intersectVTK(myPoints, Vec3D(0, -1, 0), Vec3D(0, min.Y, 0));
257  }
258  else
259  {
260  // Else, the wall normal has to have a a nonzero z-component ...
261  myPoints.emplace_back(min.X, min.Y, p.Z - ((min.Y - p.Y) * n.Y + (min.X - p.X) * n.X) / n.Z);
262  myPoints.emplace_back(min.X, max.Y, p.Z - ((max.Y - p.Y) * n.Y + (min.X - p.X) * n.X) / n.Z);
263  myPoints.emplace_back(max.X, max.Y, p.Z - ((max.Y - p.Y) * n.Y + (max.X - p.X) * n.X) / n.Z);
264  myPoints.emplace_back(max.X, min.Y, p.Z - ((min.Y - p.Y) * n.Y + (max.X - p.X) * n.X) / n.Z);
265  intersectVTK(myPoints, Vec3D(0, 0, +1), Vec3D(0, 0, max.Z));
266  intersectVTK(myPoints, Vec3D(0, 0, -1), Vec3D(0, 0, min.Z));
267  }
268 }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:32
const Quaternion & getOrientation() const
Returns the orientation of this BaseInteractable.
Definition: BaseInteractable.h:230
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:218
void intersectVTK(std::vector< Vec3D > &points, Vec3D normal, Vec3D position) const
Definition: BaseWall.cc:241
Vec3D getAxis() const
Converts the quaternions into a normal vector by rotating the vector x=(1,0,0); see See Wiki for deta...
Definition: Quaternion.cc:501
Mdouble Y
Definition: Vector.h:66
Mdouble Z
Definition: Vector.h:66
Mdouble X
the vector components
Definition: Vector.h:66

References Quaternion::getAxis(), BaseInteractable::getOrientation(), BaseInteractable::getPosition(), BaseWall::intersectVTK(), n, Vec3D::X, Vec3D::Y, and Vec3D::Z.

◆ getDistance()

Mdouble InfiniteWall::getDistance ( Vec3D  otherPosition) const

Returns the distance of the wall to the particle.

Parameters
[in]otherPositionThe position to which the distance must be computed to.
Returns
The distance of the wall to the particle.
152 {
153  return getOrientation().getDistance(otherPosition, getPosition());
154 }
static Mdouble getDistance(const Quaternion &a, const Quaternion &b)
Calculates the distance between two Quaternion: .
Definition: Quaternion.cc:188

References Quaternion::getDistance(), BaseInteractable::getOrientation(), and BaseInteractable::getPosition().

Referenced by getDistanceAndNormal(), RestrictedWall::getDistanceAndNormal(), getDistanceNormalOverlapSuperquadric(), RestrictedWall::getInteractionWith(), and RestrictedWall::writeVTK().

◆ getDistanceAndNormal()

bool InfiniteWall::getDistanceAndNormal ( const BaseParticle p,
Mdouble distance,
Vec3D normal_return 
) const
overridevirtual

Compute the distance from the wall for a given BaseParticle and return if there is a collision. If there is a collision, also return the normal vector.

Parameters
[in]pBaseParticle for which the distance to the wall must be computed.
[out]distanceDistance between the particle and the wall.
[out]normal_returnThe normal of this wall, will only be set if there is a collision.
Returns
A boolean value for whether or not there is a collision.

First the distance is checked. If there is no collision, this function will return false and give the distance. If there is a collision, the function will return true and give the distance and the normal vector of this wall. Since this function should be called before calculating any Particle-Wall interactions, it can also be used to set the normal vector in case of curved walls.

Implements BaseWall.

169 {
170  distance = getDistance(p.getPosition());
171  if (distance >= p.getWallInteractionRadius(this))
172  return false;
173  normal_return = getOrientation().getAxis();
174  return true;
175 }
Mdouble getWallInteractionRadius(const BaseWall *wall) const
returns the radius plus the interactionDistance
Definition: BaseParticle.h:386
Mdouble getDistance(Vec3D position) const
Returns the distance of the wall to the particle.
Definition: InfiniteWall.cc:151

References Quaternion::getAxis(), getDistance(), BaseInteractable::getOrientation(), BaseInteractable::getPosition(), and BaseParticle::getWallInteractionRadius().

Referenced by SimpleDrumSuperquadrics::getDistanceAndNormal().

◆ getDistanceNormalOverlapSuperquadric()

bool InfiniteWall::getDistanceNormalOverlapSuperquadric ( const SuperQuadricParticle p,
Mdouble distance,
Vec3D normal_return,
Mdouble overlap 
) const
overridevirtual

Compute the distance from the wall for a given BaseParticle and return if there is a collision. If there is a collision, also return the normal vector.

Reimplemented from BaseWall.

281 {
282  //first check: if the bounding sphere does not touch the wall, there is no contact.
284  {
285  return false;
286  }
287  Vec3D normalBodyFixed = getOrientation().getAxis();
288  p.getOrientation().rotateBack(normalBodyFixed);
289  Vec3D xWallBodyFixed = getPosition() - p.getPosition();
290  p.getOrientation().rotateBack(xWallBodyFixed);
291  Vec3D axes = p.getAxes();
292  Mdouble eps1 = p.getExponentEps1();
293  Mdouble eps2 = p.getExponentEps2();
294 
295  Vec3D furthestPoint = getFurthestPointSuperQuadric(normalBodyFixed, axes, eps1, eps2);
296  overlap = Vec3D::dot(xWallBodyFixed - furthestPoint, -normalBodyFixed);
297  if (overlap > 0)
298  {
299  Vec3D overlapBody = overlap * normalBodyFixed;
300  Vec3D contactPoint = furthestPoint - overlapBody / 2;
301  p.getOrientation().rotate(contactPoint);
302  contactPoint += p.getPosition();
303  distance = (contactPoint - overlapBody / 2 - p.getPosition()).getLength();
304  normal_return = getOrientation().getAxis();
305  return true;
306  }
307  return false;
308 }
double Mdouble
Definition: GeneralDefine.h:34
Vec3D getFurthestPointSuperQuadric(const Vec3D &normalBodyFixed, const Vec3D &axes, Mdouble eps1, Mdouble eps2) const override
Largely untested, use at your own risk for anything other than ellipsoids.
Definition: InfiniteWall.cc:312
void rotate(Vec3D &position) const
Definition: Quaternion.cc:563
void rotateBack(Vec3D &position) const
Definition: Quaternion.cc:610
Vec3D getAxes() const override
Get the axes-lengths of this superquadric. We use the super-ellipsoid definition stated in Chapter 2 ...
Definition: SuperQuadricParticle.cc:182
Mdouble getExponentEps2() const override
Get the second exponent of this superquadric. We use the super-ellipsoid definition stated in Chapter...
Definition: SuperQuadricParticle.cc:192
Mdouble getExponentEps1() const override
Get the first exponent of this superquadric. We use the super-ellipsoid definition stated in Chapter ...
Definition: SuperQuadricParticle.cc:187
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:76

References Vec3D::dot(), SuperQuadricParticle::getAxes(), Quaternion::getAxis(), getDistance(), SuperQuadricParticle::getExponentEps1(), SuperQuadricParticle::getExponentEps2(), getFurthestPointSuperQuadric(), BaseInteractable::getOrientation(), BaseInteractable::getPosition(), BaseParticle::getWallInteractionRadius(), Quaternion::rotate(), and Quaternion::rotateBack().

Referenced by SimpleDrumSuperquadrics::getDistanceNormalOverlapSuperquadric().

◆ getFurthestPointSuperQuadric()

Vec3D InfiniteWall::getFurthestPointSuperQuadric ( const Vec3D normalBodyFixed,
const Vec3D axes,
Mdouble  eps1,
Mdouble  eps2 
) const
overridevirtual

Largely untested, use at your own risk for anything other than ellipsoids.

Reimplemented from BaseWall.

314 {
315  Vec3D furthestPoint;
316  if (std::abs(normalBodyFixed.X) > 1e-10)
317  {
318  Mdouble alpha = std::abs(normalBodyFixed.Y * axes.Y / normalBodyFixed.X / axes.X);
319  Mdouble gamma = std::pow(1 + std::pow(alpha, 2 / eps2), eps2 / eps1 - 1);
320  Mdouble beta = std::pow(gamma * std::abs(normalBodyFixed.Z * axes.Z / normalBodyFixed.X / axes.X),
321  eps1 / (2 - eps1));
322  furthestPoint.X = axes.X * mathsFunc::sign(normalBodyFixed.X) /
323  (std::pow(std::pow(1 + std::pow(alpha, 2 / eps2), eps2 / eps1) + std::pow(beta, 2 / eps1),
324  eps1 / 2));
325  furthestPoint.Y = axes.Y * alpha / axes.X * std::abs(furthestPoint.X) * mathsFunc::sign(normalBodyFixed.Y);
326  furthestPoint.Z = axes.Z * beta / axes.X * std::abs(furthestPoint.X) * mathsFunc::sign(normalBodyFixed.Z);
327  }
328  else if (std::abs(normalBodyFixed.Y) > 1e-10)
329  {
330  Mdouble beta = std::pow(std::abs(normalBodyFixed.Z * axes.Z / normalBodyFixed.Y / axes.Y), eps1 / (2 - eps1));
331  furthestPoint.Y = axes.Y / std::pow((1 + std::pow(beta, 2/eps1)), eps1 / 2) * mathsFunc::sign(normalBodyFixed.Y);
332  furthestPoint.Z = axes.Z / axes.Y * std::abs(furthestPoint.Y) * beta * mathsFunc::sign(normalBodyFixed.Z);
333  }
334  else
335  {
336  furthestPoint.Z = axes.Z * mathsFunc::sign(normalBodyFixed.Z);
337  }
338  return furthestPoint;
339 }
Mdouble beta(Mdouble z, Mdouble w)
This is the beta function, returns the approximation based on cmath's implementation of ln(gamma)
Definition: ExtendedMath.cc:164
int sign(T val)
This is a sign function, it returns -1 for negative numbers, 1 for positive numbers and 0 for 0.
Definition: ExtendedMath.h:97
Mdouble gamma(Mdouble gamma_in)
This is the gamma function returns the true value for the half integer value.
Definition: ExtendedMath.cc:137

References mathsFunc::beta(), mathsFunc::gamma(), mathsFunc::sign(), Vec3D::X, Vec3D::Y, and Vec3D::Z.

Referenced by getDistanceNormalOverlapSuperquadric(), and SimpleDrumSuperquadrics::getFurthestPointSuperQuadric().

◆ getName()

std::string InfiniteWall::getName ( ) const
overridevirtual

Writes the InfiniteWall to an output stream, usually a restart file.

Returns the name of the object, in this case the string "InfiniteWall".

Returns
The string "InfiniteWall", which is the name of this class.

Implements BaseObject.

206 {
207  return "InfiniteWall";
208 }

◆ getNormal()

Vec3D InfiniteWall::getNormal ( ) const

Access function for normal.

Returns
The 3D vector that represents the normal to the wall.
214 {
215  return getOrientation().getAxis();
216 }

References Quaternion::getAxis(), and BaseInteractable::getOrientation().

Referenced by LawinenBox::actionsBeforeTimeStep(), CGHandler::computeContactPoints(), save(), FlowFrontChute::stretch(), LawinenBox::writeEneTimeStep(), and RestrictedWall::writeVTK().

◆ oldRead()

void InfiniteWall::oldRead ( std::istream &  is)

Reads InfiniteWall from an old-style restart file.

Parameters
[in]isThe input stream from which the InfiniteWall old style is read.
191 {
192  std::string dummy;
193  Vec3D velocity;
194  Vec3D position;
195  Vec3D normal;
196  is >> dummy >> normal >> dummy >> position >> dummy >> velocity;
197  setPosition(position);
198  setVelocity(velocity);
199  setOrientation(Quaternion(normal));
200 }
virtual void setOrientation(const Quaternion &orientation)
Sets the orientation of this BaseInteractable.
Definition: BaseInteractable.h:260
void setVelocity(const Vec3D &velocity)
set the velocity of the BaseInteractable.
Definition: BaseInteractable.cc:350
This class contains the 4 components of a quaternion and the standard operators and functions needed ...
Definition: Quaternion.h:63

References BaseInteractable::setOrientation(), BaseInteractable::setPosition(), and BaseInteractable::setVelocity().

◆ read()

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

Reads InfiniteWall from a restart file.

Parameters
[in]isThe input stream from which the InfiniteWall is read. Only needed for backward compatibility.

Reimplemented from BaseWall.

181 {
182  BaseWall::read(is);
183  Vec3D normal;
184  if (helpers::readOptionalVariable(is,"normal",normal)) setNormal(normal);
185 }
void read(std::istream &is) override
Function that reads a BaseWall from an input stream, usually a restart file.
Definition: BaseWall.cc:78
bool readOptionalVariable(std::istream &is, const std::string &name, T &variable)
Reads optional variables in the restart file.
Definition: FileIOHelpers.h:82

References BaseWall::read(), helpers::readOptionalVariable(), and setNormal().

Referenced by SimpleDrumSuperquadrics::read(), and RestrictedWall::read().

◆ set() [1/2]

void InfiniteWall::set ( Vec3D  normal,
Mdouble  positionInNormalDirection 
)

Defines a standard wall by computing normal*position = point and using the overloaded function set(Vec3D, vec3D).

Deprecated:
In Mercury 2, the user will have to use the new interface, namely set(Vec3D, Vec3D).

Defines a standard wall, given an normal vector pointing into the wall (i.e. out of the flow domain), to give a plane defined by normal*x=position

Parameters
[in]normalA Vec3D that represents the normal vector to the wall.
[in]positionInNormalDirectionThe position of the wall in the direction of the normal vector.
Deprecated:
InfiniteWall::set(Vec3D, Mdouble) is deprecated. Use set(Vec3D, Vec3D) instead.
142 {
143  logger(WARN, "InfiniteWall::set(Vec3D, Mdouble) is deprecated. Use set(Vec3D, Vec3D) instead.");
144  set(normal, positionInNormalDirection * normal);
145 }
@ WARN
void set(Vec3D normal, Vec3D point)
Defines a standard wall, given an outward normal vector s.t. normal*x=normal*point for all x of the w...
Definition: InfiniteWall.cc:118

References logger, set(), and WARN.

◆ set() [2/2]

void InfiniteWall::set ( Vec3D  normal,
Vec3D  point 
)

Defines a standard wall, given an outward normal vector s.t. normal*x=normal*point for all x of the wall.

119 {
120  setNormal(normal);
121  setPosition(point);
122 }

References setNormal(), and BaseInteractable::setPosition().

Referenced by IntersectionOfWalls::addObject(), NautaMixer::addTopWall(), ClosedCSCWalls::ClosedCSCWalls(), Chute::createBottom(), LevelSetWall::createVTK(), CSCWalls::CSCWalls(), SimpleDrumSuperquadrics::getDistanceAndNormal(), SimpleDrumSuperquadrics::getDistanceNormalOverlapSuperquadric(), GranuDrum::GranuDrum(), InitialConditions< SpeciesType >::InitialConditions(), loadingTest(), main(), MercuryCGSelfTest::MercuryCGSelfTest(), normalAndTangentialLoadingTest(), ParticleCreation::ParticleCreation(), WallHandler::readAndCreateOldObject(), DPMBase::readParAndIniFiles(), set(), Slide::set_Walls(), T_protectiveWall::setupInitialConditions(), LawinenBox::setupInitialConditions(), ClosedCSCWalls::setupInitialConditions(), MercuryLogo::setupInitialConditions(), SmoothChute::setupInitialConditions(), multiParticleT1::setupInitialConditions(), VerticalMixer::setupInitialConditions(), Binary::setupInitialConditions(), my_problem::setupInitialConditions(), FreeCooling2DinWallsDemo::setupInitialConditions(), FreeCooling3DinWallsDemo::setupInitialConditions(), HeaterBoundaryTest::setupInitialConditions(), HourGlass2D::setupInitialConditions(), HourGlass::setupInitialConditions(), Cstatic2d::setupInitialConditions(), AngleOfRepose::setupInitialConditions(), SilbertPeriodic::setupInitialConditions(), Drum::setupInitialConditions(), HertzSelfTest::setupInitialConditions(), MindlinSelfTest::setupInitialConditions(), Penetration::setupInitialConditions(), SegregationPeriodic::setupInitialConditions(), Chutebelt::setupInitialConditions(), ConstantMassFlowMaserBoundaryMixedSpeciesSelfTest::setupInitialConditions(), ConstantMassFlowMaserSelfTest::setupInitialConditions(), InsertionBoundarySelfTest::setupInitialConditions(), PolydisperseInsertionBoundarySelfTest::setupInitialConditions(), CGHandlerSelfTest::setupInitialConditions(), NewtonsCradleSelftest::setupInitialConditions(), NewtonsCradleSelfTest::setupInitialConditions(), SquarePacking::setupInitialConditions(), ChargedBondedParticleUnitTest::setupInitialConditions(), WallParticleCollision::setupInitialConditions(), my_problem_HGRID::setupInitialConditions(), FreeFallInteractionSelfTest::setupInitialConditions(), FreeFallSelfTest::setupInitialConditions(), LiquidMigrationSelfTest::setupInitialConditions(), DPM::setupInitialConditions(), ObliqueImpactSelfTest::setupInitialConditions(), TwoBondedParticleElasticCollision::setupInitialConditions(), TwoParticleElasticCollisionInteraction::setupInitialConditions(), TwoParticleElasticCollision::setupInitialConditions(), CoilSelfTest::setupInitialConditions(), UnionOfWalls::setupInitialConditions(), DrumRot::setupInitialConditions(), RotatingDrum::setupInitialConditions(), Wall::setupInitialConditions(), BouncingSuperQuadric::setupInitialConditions(), EllipsoidsBouncingOnWallDemo::setupInitialConditions(), GranularCollapse::setupInitialConditions(), EllipticalSuperQuadricCollision::setupInitialConditions(), protectiveWall::setupInitialConditions(), Tutorial11::setupInitialConditions(), Tutorial12::setupInitialConditions(), Tutorial3::setupInitialConditions(), Tutorial4::setupInitialConditions(), Tutorial7::setupInitialConditions(), Tutorial8::setupInitialConditions(), Tutorial9::setupInitialConditions(), ParticleWallInteraction::setupInitialConditions(), EvaporationAndHeatTest::setupInitialConditions(), ExtremeOverlapWithWallsUnitTest::setupInitialConditions(), FreeFallHertzMindlinUnitTest::setupInitialConditions(), FreeFall::setupInitialConditions(), FullRestartTest::setupInitialConditions(), InclinedPlane::setupInitialConditions(), MpiMaserChuteTest::setupInitialConditions(), MovingWalls::setupInitialConditions(), MultiParticlesInsertion::setupInitialConditions(), MpiPeriodicBoundaryUnitTest::setupInitialConditions(), WallSpecies::setupInitialConditions(), AreaVTK::setupInitialConditions(), ChuteBottom::setupInitialConditions(), CoupledProblem::setupMercury(), ContactDetectionWithWallTester::setupParticleAndWall(), Chute::setupSideWalls(), SimpleDrumSuperquadrics::SimpleDrumSuperquadrics(), and Slide::Slide().

◆ setNormal()

void InfiniteWall::setNormal ( Vec3D  normal)

Changes the normal of the InfiniteWall.

Parameters
[in]normalThe vector normal to the wall.
128 {
129  setOrientationViaNormal(normal);
130 }
void setOrientationViaNormal(Vec3D normal)
Sets the orientation of this BaseInteractable by defining the vector that results from the rotation o...
Definition: BaseInteractable.cc:199

References BaseInteractable::setOrientationViaNormal().

Referenced by InfiniteWall(), load(), main(), read(), set(), MaserRepeatedOutInMPI2Test::setupInitialConditions(), and MovingWall::setupInitialConditions().

◆ writeVTK()

void InfiniteWall::writeVTK ( VTKContainer vtk) const
overridevirtual

Adds the vtk wall representation to the VTK container

Reimplemented from BaseWall.

271 {
272  std::vector<Vec3D> points;
273  createVTK(points);
274  addToVTK(points, vtk);
275 }
static void addToVTK(const std::vector< Vec3D > &points, VTKContainer &vtk)
Takes the points provided and adds a triangle strip connecting these points to the vtk container.
Definition: BaseWall.cc:501

References BaseWall::addToVTK(), and createVTK().


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