TriangleWall Class Reference

A TriangleWall is convex polygon defined as an intersection of InfiniteWall's. More...

#include <TriangleWall.h>

+ Inheritance diagram for TriangleWall:

Public Member Functions

 TriangleWall ()=default
 Default constructor. More...
 
 TriangleWall (const TriangleWall &other)=default
 Copy constructor. More...
 
 ~TriangleWall () override=default
 Destructor. More...
 
TriangleWallcopy () const override
 Wall copy method. It calls the copy constructor of this Wall, useful for polymorphism. More...
 
std::string getName () const override
 Returns the name of the object, here the string "TriangleWall". More...
 
void read (std::istream &is) override
 Reads an TriangleWall from an input stream, for example a restart file. More...
 
void write (std::ostream &os) const override
 Writes an TriangleWall to an output stream, for example a restart file. More...
 
void setVertices (Vec3D A, Vec3D B, Vec3D C)
 Sets member variables such that the wall represents a triangle with vertices A, B, C. More...
 
std::array< Vec3D, 3 > getVertices () const
 
void move (const Vec3D &move) override
 
void setPosition (const Vec3D &position) override
 Sets the position of this BaseInteractable. More...
 
void setOrientation (const Quaternion &orientation) override
 Sets the orientation of this BaseInteractable. More...
 
void setVertices (Vec3D A, Vec3D B, Vec3D C, Vec3D position)
 Same as setVertices(A,B,C), but sets the position explicitly. The position is important when you rotate the wall, as the wall will be rotated around this position. More...
 
void writeVTK (VTKContainer &vtk) const override
 
const Vec3DgetVertex (unsigned i) const
 Returns the position of a vertex. 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...
 
void rotate (const Vec3D &angularVelocity) override
 Rotates this BaseInteractable. More...
 
bool isLocal (Vec3D &min, Vec3D &max) const override
 
bool isInsideTriangle (const Vec3D &point) const
 
void moveVertex (unsigned index, const Vec3D &dP)
 Updates the indexed vertex by a given change in position. More...
 
void moveVertices (const std::array< Vec3D, 3 > &dPs)
 Updates all vertices by the given changes in position. More...
 
bool isFaceContact (const Vec3D &normal) const override
 
- Public Member Functions inherited from BaseWall
 BaseWall ()
 Default constructor. More...
 
 BaseWall (const BaseWall &w)
 Copy constructor. More...
 
 ~BaseWall () override
 Default destructor. More...
 
virtual bool getDistanceNormalOverlap (const BaseParticle &P, Mdouble &distance, Vec3D &normal_return, Mdouble &overlap) const
 
virtual bool getDistanceNormalOverlapSuperquadric (const SuperQuadricParticle &p, Mdouble &distance, Vec3D &normal_return, Mdouble &overlap) const
 
virtual Vec3D getFurthestPointSuperQuadric (const Vec3D &normalBodyFixed, const Vec3D &axes, Mdouble eps1, Mdouble eps2) 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...
 
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...
 
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...
 
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
 
- 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
 

Private Member Functions

void updateVertexAndNormal ()
 

Private Attributes

std::array< Vec3D, 3 > vertexInLabFrame_
 
std::array< Vec3D, 3 > vertex_
 
Vec3D vertexMin_
 
Vec3D vertexMax_
 
std::array< Vec3D, 3 > edgeNormal_
 
std::array< Vec3D, 3 > edge_
 
std::array< double, 3 > edgeLength_
 
Vec3D faceNormal_
 

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 TriangleWall is convex polygon defined as an intersection of InfiniteWall's.

It can be defined as the intersection of a set of InfiniteWall's, defined by the normal vector into the wall and a point on the wall:

//for each wall, specify a normal and position vector
w.addObject(Vec3D(-1.0,0.0,0.0),Vec3D(1.0,0.0,0.0));
w.addObject(Vec3D(1.0,0.0,0.0),Vec3D(0.0,0.0,0.0));
w.addObject(Vec3D(0.0,-1.0,0.0),Vec3D(0.0,1.0,0.0));
w.addObject(Vec3D(0.0,1.0,0.0),Vec3D(0.0,0.0,0.0));
wallHandler.copyAndAddObject(w);
A TriangleWall is convex polygon defined as an intersection of InfiniteWall's.
Definition: TriangleWall.h:57
Definition: Vector.h:51

A particle of radius r and position x touches an InfiniteWall with normal n and position p if \(p-n\cdot x\leq r\) (note 'touching particles' also includes particles that are completely enclosed inside the wall). A particle touches an TriangleWall if it touches all InfiniteWall objects (shown in the image below).

For a demonstration on how to use this class, see T8: Motion of a particle in a box with an obstacle and Flow through a 3D hourglass/silo (shown in the image below).

Constructor & Destructor Documentation

◆ TriangleWall() [1/2]

TriangleWall::TriangleWall ( )
default

Default constructor.

Referenced by copy().

◆ TriangleWall() [2/2]

TriangleWall::TriangleWall ( const TriangleWall other)
default

Copy constructor.

◆ ~TriangleWall()

TriangleWall::~TriangleWall ( )
overridedefault

Destructor.

Member Function Documentation

◆ copy()

TriangleWall* TriangleWall::copy ( ) const
inlineoverridevirtual

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

Implements BaseWall.

79  { return new TriangleWall(*this); }
TriangleWall()=default
Default constructor.

References TriangleWall().

◆ getDistanceAndNormal()

bool TriangleWall::getDistanceAndNormal ( const BaseParticle p,
Mdouble distance,
Vec3D normal 
) 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 we want to calculate the distance and whether it collided of.
[out]distanceThe distance of the BaseParticle to this wall.
[out]normal_returnIf there was a collision, the normal vector to this wall will be placed here.
Returns
A boolean which says whether or not there was a collision.

This function computes whether or not there is a collision between a given BaseParticle and this TriangleWall. If there is a collision, this function also computes the distance between the BaseParticle and TriangleWall and the normal of the TriangleWall at the intersection point. It does this by calling TriangleWall::getDistanceAndNormal(const Vec3D& , Mdouble , Mdouble&, Vec3D&) const. 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.

58 {
59  const Vec3D position = p.getPosition(); // note: no transfer to lab coordinates
60  const Mdouble distanceMax = p.getWallInteractionRadius(this);
61 
62  // compute distance from face
63  //get distance from particle position to the face
64  const Mdouble signedDistance = Vec3D::dot(position-vertex_[0], faceNormal_);
65  distance = fabs(signedDistance);
66 
67  // check if any contact is possible
68  if (distance >= distanceMax) return false;
69 
70  // compute distance from edges
71  const std::array<Vec3D,3> edgeBranch {position - vertex_[0], position - vertex_[1], position - vertex_[2]};
72  const std::array<double,3> edgeDistance {Vec3D::dot(edgeBranch[0], edgeNormal_[0]), Vec3D::dot(edgeBranch[1], edgeNormal_[1]), Vec3D::dot(edgeBranch[2], edgeNormal_[2])};
73 
74  // find edge with largest distance (this will be the edge if there is a edge contact)
75  const size_t id = (edgeDistance[1] > edgeDistance[2]) ?
76  (edgeDistance[0] > edgeDistance[1] ? 0 : 1) : (edgeDistance[0] > edgeDistance[2] ? 0 : 2);
77 
78  // check if there will be contact
79  if (edgeDistance[id] > distanceMax) return false;
80 
81  // determine if it is a face contact
82  const Vec3D posProjected = position - signedDistance * faceNormal_;
83  if (edgeDistance[id] <= 0 && isInsideTriangle(posProjected)){
84  normal = (signedDistance >= 0) ? -faceNormal_ : faceNormal_;
85  return true;
86  }
87 
88  // determine if it is an edge or vertex contact
89  const double positionAlongEdge = Vec3D::dot(edgeBranch[id], edge_[id]);
90  if (positionAlongEdge < 0) {
91  //possible contact with left vertex
92  distance = edgeBranch[id].getLength();
93  if (distance > distanceMax) return false;
94  normal = edgeBranch[id] / -distance;
95  } else if (positionAlongEdge > edgeLength_[id]) {
96  //contact with right vertex
97  const size_t idRight = (id + 1) % 3;
98  distance = edgeBranch[idRight].getLength();
99  if (distance > distanceMax) return false;
100  normal = edgeBranch[idRight] / -distance;
101  } else {
102  // edge contact
103  normal = edge_[id] * positionAlongEdge - edgeBranch[id];
104  distance = normal.getLength();
105  if (distance > distanceMax) return false;
106  normal /= distance;
107  }
108  return true;
109 }
double Mdouble
Definition: GeneralDefine.h:34
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:218
Mdouble getWallInteractionRadius(const BaseWall *wall) const
returns the radius plus the interactionDistance
Definition: BaseParticle.h:386
Vec3D faceNormal_
Definition: TriangleWall.h:185
std::array< Vec3D, 3 > edgeNormal_
Definition: TriangleWall.h:178
std::array< Vec3D, 3 > vertex_
Definition: TriangleWall.h:167
std::array< double, 3 > edgeLength_
Definition: TriangleWall.h:180
std::array< Vec3D, 3 > edge_
Definition: TriangleWall.h:179
bool isInsideTriangle(const Vec3D &point) const
Definition: TriangleWall.cc:250
static Mdouble getLength(const Vec3D &a)
Calculates the length of a Vec3D: .
Definition: Vector.cc:331
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:76

References Vec3D::dot(), edge_, edgeLength_, edgeNormal_, faceNormal_, Vec3D::getLength(), BaseInteractable::getPosition(), BaseParticle::getWallInteractionRadius(), isInsideTriangle(), and vertex_.

◆ getName()

std::string TriangleWall::getName ( ) const
inlineoverridevirtual

Returns the name of the object, here the string "TriangleWall".

Implements BaseObject.

85  { return "TriangleWall"; }

◆ getVertex()

const Vec3D& TriangleWall::getVertex ( unsigned  i) const
inline

Returns the position of a vertex.

124 {return vertex_[i];}
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51

References constants::i, and vertex_.

◆ getVertices()

std::array<Vec3D,3> TriangleWall::getVertices ( ) const
inline
105 {return vertex_;}

References vertex_.

Referenced by SCoupling< M, O >::updateTriangleWall().

◆ isFaceContact()

bool TriangleWall::isFaceContact ( const Vec3D normal) const
inlineoverridevirtual

Reimplemented from BaseInteractable.

151  {
152  return (normal==faceNormal_ or normal == -faceNormal_) ? true : false;
153  }

References faceNormal_.

◆ isInsideTriangle()

bool TriangleWall::isInsideTriangle ( const Vec3D point) const
251 {
252  const Vec3D branch0 = point - vertex_[0];
253  const Vec3D branch1 = point - vertex_[1];
254  const Vec3D branch2 = point - vertex_[2];
255 
256  //compute total area
257  Mdouble s = sqrt(Vec3D::cross(vertex_[1] - vertex_[0], vertex_[2] - vertex_[1]).getLengthSquared());
258  //compute barycentric coordinates
259  Mdouble s0 = sqrt(Vec3D::cross(branch0, branch1).getLengthSquared())/s;
260  Mdouble s1 = sqrt(Vec3D::cross(branch1, branch2).getLengthSquared())/s;
261  Mdouble s2 = sqrt(Vec3D::cross(branch2, branch0).getLengthSquared())/s;
262  return (1 > s0 > 0 && 1 > s1 > 0 && 1 > s2 > 0);
263 }
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:163

References Vec3D::cross(), and vertex_.

Referenced by getDistanceAndNormal().

◆ isLocal()

bool TriangleWall::isLocal ( Vec3D min,
Vec3D max 
) const
overridevirtual

if isLocal returns true and the DPM class derived from MercuryBase, the hGrid will be used to find wall-particle contacts, using min/max.

Reimplemented from BaseWall.

244 {
245  min = vertexMin_;
246  max = vertexMax_;
247  return true;
248 }
Vec3D vertexMax_
Definition: TriangleWall.h:173
Vec3D vertexMin_
Definition: TriangleWall.h:172

References vertexMax_, and vertexMin_.

◆ move()

void TriangleWall::move ( const Vec3D move)
overridevirtual

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 from BaseInteractable.

185 {
188 }
virtual void move(const Vec3D &move)
Moves this BaseInteractable by adding an amount to the position.
Definition: BaseInteractable.cc:215
void updateVertexAndNormal()
Definition: TriangleWall.cc:218
void move(const Vec3D &move) override
Definition: TriangleWall.cc:184

References BaseInteractable::move(), and updateVertexAndNormal().

◆ moveVertex()

void TriangleWall::moveVertex ( unsigned  index,
const Vec3D dP 
)

Updates the indexed vertex by a given change in position.

Parameters
indexVertex index, i.e. 0, 1 or 2
dPChange in position
266 {
267  vertexInLabFrame_[index] += dP;
269 }
std::array< Vec3D, 3 > vertexInLabFrame_
Definition: TriangleWall.h:162

References updateVertexAndNormal(), and vertexInLabFrame_.

◆ moveVertices()

void TriangleWall::moveVertices ( const std::array< Vec3D, 3 > &  dPs)

Updates all vertices by the given changes in position.

Parameters
dPsChange in position for each vertex
272 {
273  vertexInLabFrame_[0] += dPs[0];
274  vertexInLabFrame_[1] += dPs[1];
275  vertexInLabFrame_[2] += dPs[2];
277 }

References updateVertexAndNormal(), and vertexInLabFrame_.

◆ read()

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

Reads an TriangleWall from an input stream, for example a restart file.

Parameters
[in]isThe input stream from which the TriangleWall is read, usually a restart file.

Reimplemented from BaseWall.

124 {
125  BaseWall::read(is);
126  std::string dummy;
127  is >> dummy;
128  for (int i = 0; i < 3; i++)
129  {
130  is >> vertexInLabFrame_[i];
131  }
133 }
void read(std::istream &is) override
Function that reads a BaseWall from an input stream, usually a restart file.
Definition: BaseWall.cc:78

References constants::i, BaseWall::read(), updateVertexAndNormal(), and vertexInLabFrame_.

◆ rotate()

void TriangleWall::rotate ( const Vec3D angularVelocityDt)
overridevirtual

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 from BaseInteractable.

112 {
113  if (!angularVelocityDt.isZero())
114  {
115  BaseWall::rotate(angularVelocityDt);
117  }
118 }
virtual void rotate(const Vec3D &angularVelocityDt)
Rotates this BaseInteractable.
Definition: BaseInteractable.cc:230
bool isZero() const
Checks if ALL elements are zero.
Definition: Vector.h:113

References Vec3D::isZero(), BaseInteractable::rotate(), and updateVertexAndNormal().

◆ setOrientation()

void TriangleWall::setOrientation ( const Quaternion orientation)
overridevirtual

Sets the orientation of this BaseInteractable.

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

Parameters
[in]orientationReference to Vec3D storing the orientation of the particle.

Reimplemented from BaseInteractable.

197 {
198  BaseWall::setOrientation(orientation);
200 }
virtual void setOrientation(const Quaternion &orientation)
Sets the orientation of this BaseInteractable.
Definition: BaseInteractable.h:260

References BaseInteractable::setOrientation(), and updateVertexAndNormal().

◆ setPosition()

void TriangleWall::setPosition ( const Vec3D position)
overridevirtual

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 from BaseInteractable.

191 {
192  BaseWall::setPosition(position);
194 }
virtual void setPosition(const Vec3D &position)
Sets the position of this BaseInteractable.
Definition: BaseInteractable.h:239

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

◆ setVertices() [1/2]

void TriangleWall::setVertices ( Vec3D  A,
Vec3D  B,
Vec3D  C 
)

Sets member variables such that the wall represents a triangle with vertices A, B, C.

  • position_ is set to the center of mass of the wall
  • vertexInLabFrame_ is set relative to the position
  • updateVertexAndNormal is called to set the remaining variables
166 {
167  BaseWall::setPosition((A + B + C) / 3);
168  BaseWall::setOrientation({1, 0, 0, 0});
170  vertexInLabFrame_[1] = B - getPosition();
171  vertexInLabFrame_[2] = C - getPosition();
172  //vertexInLabFrame_[0] = A;
173  //vertexInLabFrame_[1] = B;
174  //vertexInLabFrame_[2] = C;
176 }
@ A
Definition: StatisticsVector.h:42

References A, BaseInteractable::getPosition(), BaseInteractable::setOrientation(), BaseInteractable::setPosition(), updateVertexAndNormal(), and vertexInLabFrame_.

Referenced by MercuryProblem::createTriangleWall(), SCoupling< M, O >::createTriangleWall(), main(), WallHandler::readTriangleWall(), TriangleMeshWall::refineTriangle(), TriangleMeshWall::set(), GetDistanceAndNormalForTriangleWall::setupInitialConditions(), Drum::setupInitialConditions(), Penetration::setupInitialConditions(), Silo::setupInitialConditions(), GetDistanceAndNormalForTriangleWalls::setupInitialConditions(), RollingOverTriangleWalls::setupInitialConditions(), and SCoupling< M, O >::updateTriangleWall().

◆ setVertices() [2/2]

void TriangleWall::setVertices ( Vec3D  A,
Vec3D  B,
Vec3D  C,
Vec3D  position 
)

Same as setVertices(A,B,C), but sets the position explicitly. The position is important when you rotate the wall, as the wall will be rotated around this position.

203 {
204  BaseWall::setPosition(position);
205  BaseWall::setOrientation({1, 0, 0, 0});
207  vertexInLabFrame_[1] = B - getPosition();
208  vertexInLabFrame_[2] = C - getPosition();
210 }

References A, BaseInteractable::getPosition(), BaseInteractable::setOrientation(), BaseInteractable::setPosition(), updateVertexAndNormal(), and vertexInLabFrame_.

◆ updateVertexAndNormal()

void TriangleWall::updateVertexAndNormal ( )
private

This function should be called after setting either position_ or vertexInLabFrame_.

  • vertex is set to the vertex position in the real coordinate system (rotated and shifted)
  • vertexMin_/vertexMax_ is set to define a bounding box around the wall (for contact detection)
  • edge_, edgeNormal_ and faceNormal_ is updated (stored for quick computation of contact point)
219 {
220  for (int i = 0; i < 3; i++)
221  {
224  vertex_[i] += getPosition();
225  }
228 
229  edge_ = {vertex_[1] - vertex_[0], vertex_[2] - vertex_[1], vertex_[0] - vertex_[2]};
232 
233  for (int i = 0; i < 3; i++)
234  {
236  edgeLength_[i] = edge_[i].getLength();
237  edge_[i] /= edgeLength_[i];
238  edgeNormal_[i].normalise();
239  }
240  //logger(INFO,"vertex %,%,% edge %,%,% face %",vertex_[0],vertex_[1],vertex_[2],edgeNormal_[0],edgeNormal_[1],edgeNormal_[2],faceNormal_);
241 }
const Quaternion & getOrientation() const
Returns the orientation of this BaseInteractable.
Definition: BaseInteractable.h:230
void rotate(Vec3D &position) const
Definition: Quaternion.cc:563
static Vec3D max(const Vec3D &a, const Vec3D &b)
Calculates the pointwise maximum of two Vec3D.
Definition: Vector.cc:89
static Vec3D min(const Vec3D &a, const Vec3D &b)
Calculates the pointwise minimum of two Vec3D.
Definition: Vector.cc:102
void normalise()
Makes this Vec3D unit length.
Definition: Vector.cc:123

References Vec3D::cross(), edge_, edgeLength_, edgeNormal_, faceNormal_, BaseInteractable::getOrientation(), BaseInteractable::getPosition(), constants::i, Vec3D::max(), Vec3D::min(), Vec3D::normalise(), Quaternion::rotate(), vertex_, vertexInLabFrame_, vertexMax_, and vertexMin_.

Referenced by move(), moveVertex(), moveVertices(), read(), rotate(), setOrientation(), setPosition(), and setVertices().

◆ write()

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

Writes an TriangleWall to an output stream, for example a restart file.

Parameters
[in]osThe output stream where the TriangleWall must be written to, usually a restart file.

Reimplemented from BaseWall.

140 {
141  BaseWall::write(os);
142  os << " vertexInLabFrame ";
143  for (int i = 0; i < 3; i++)
144  {
145  os << ' ' << vertexInLabFrame_[i];
146  }
147 }
void write(std::ostream &os) const override
Function that writes a BaseWall to an output stream, usually a restart file.
Definition: BaseWall.cc:100

References constants::i, vertexInLabFrame_, and BaseWall::write().

◆ writeVTK()

void TriangleWall::writeVTK ( VTKContainer vtk) const
overridevirtual

adds extra information to the points and triangleStrips vectors needed to plot the wall in vtk format

Parameters
pointsCoordinates of the vertices of the triangulated surfaces (in the VTK file this is called POINTS)
triangleStripsIndices of three vertices forming one triangulated surface (in the VTK file this is called CELL)

Reimplemented from BaseWall.

150 {
151  const unsigned long s = vtk.points.size();
152  for (auto v : vertex_)
153  {
154  vtk.points.push_back(v);
155  }
156  std::vector<double> cell;
157  cell.reserve(3);
158  cell.push_back(s);
159  cell.push_back(s + 1);
160  cell.push_back(s + 2);
161  //cell.push_back(s);
162  vtk.triangleStrips.push_back(cell);
163 }
std::vector< std::vector< double > > triangleStrips
Definition: BaseWall.h:40
std::vector< Vec3D > points
Definition: BaseWall.h:39

References VTKContainer::points, VTKContainer::triangleStrips, and vertex_.

Member Data Documentation

◆ edge_

std::array<Vec3D, 3> TriangleWall::edge_
private

◆ edgeLength_

std::array<double, 3> TriangleWall::edgeLength_
private

◆ edgeNormal_

std::array<Vec3D, 3> TriangleWall::edgeNormal_
private

stores the wall normal n in n.x=p

Referenced by getDistanceAndNormal(), and updateVertexAndNormal().

◆ faceNormal_

Vec3D TriangleWall::faceNormal_
private

stores the face normal, not rotated into the lab frame; thus, if the wall rotates, this normal has to be rotated as well

Referenced by getDistanceAndNormal(), isFaceContact(), and updateVertexAndNormal().

◆ vertex_

std::array<Vec3D, 3> TriangleWall::vertex_
private

stores the position of the vertices relative to the position of the wall but not rotated into the lab frame; thus, if the wall rotates, these vertices have to be rotated as well

Referenced by getDistanceAndNormal(), getVertex(), getVertices(), isInsideTriangle(), updateVertexAndNormal(), and writeVTK().

◆ vertexInLabFrame_

std::array<Vec3D, 3> TriangleWall::vertexInLabFrame_
private

stores the position of the vertices relative to the position of the wall and rotated into the lab frame;

Referenced by moveVertex(), moveVertices(), read(), setVertices(), updateVertexAndNormal(), and write().

◆ vertexMax_

Vec3D TriangleWall::vertexMax_
private

Referenced by isLocal(), and updateVertexAndNormal().

◆ vertexMin_

Vec3D TriangleWall::vertexMin_
private

stores the min and max coordinate values of the vertices (needed for hGrid)

Referenced by isLocal(), and updateVertexAndNormal().


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