MeshTriangle Class Reference

MeshTriangle implements a triangle whose vertex positions are defined by three particles. More...

#include <MeshTriangle.h>

+ Inheritance diagram for MeshTriangle:

Public Member Functions

 MeshTriangle ()=default
 Default constructor. More...
 
 MeshTriangle (const MeshTriangle &other)=default
 Copy constructor. More...
 
 ~MeshTriangle () override=default
 Destructor. More...
 
MeshTrianglecopy () 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 "MeshTriangle". More...
 
void read (std::istream &is) override
 Reads an MeshTriangle from an input stream, for example a restart file. More...
 
void write (std::ostream &os) const override
 Writes an MeshTriangle 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...
 
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 setVertexVelocities (Vec3D A, Vec3D B, Vec3D C)
 Sets the velocity of the vertex points. More...
 
std::array< Vec3D, 3 > getVertices () const
 Returns an array of the vertex coordinates. More...
 
Vec3D getFaceNormal () const
 Returns the face normal. More...
 
Mdouble getArea () const
 Returns the area of the triangle. More...
 
void move (const Vec3D &move) override
 
void setVertexIds (unsigned int i, unsigned int j, unsigned int k)
 sets the ids of the vertex particles. Calls retrieveVertexParticles. More...
 
std::array< unsigned int, 3 > getVertexIds () const
 Returns an array containing the ids of the vertex particles. More...
 
void writeVTK (VTKContainer &vtk) const override
 
BaseInteractiongetInteractionWith (BaseParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler) override
 
void checkInteractions (InteractionHandler *interactionHandler, unsigned int timeStamp) override
 Checks, if the forces of all interctions should be applied. 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 getDistanceNormalOverlapType (const BaseParticle &p, Mdouble &distance, Vec3D &normal, Mdouble &overlap, unsigned int &type) const
 
const Vec3D getVelocityAtContact (const Vec3D &contact) const override
 Calculates the local velocity at a specified point. More...
 
const Vec3D getBaricentricWeight (const Vec3D &contact) const
 Calculates the barycentric weight of a specified point. More...
 
void rotate (const Vec3D &angularVelocity) override
 Rotates this BaseInteractable. More...
 
bool isLocal (Vec3D &min, Vec3D &max) const override
 Determines if the triangle is considered local. More...
 
bool isInsideTriangle (const Vec3D &point) const
 Determines if a given point is within the triangle. More...
 
Mdouble getInvMass () const override
 
void setMass (Mdouble mass)
 
void actionsOnRestart () override
 Actions executed on restart. More...
 
void actionsAfterParticleGhostUpdate () override
 actionsPerformed after the position update of (ghost-) particles. More...
 
void handleParticleRemoval (unsigned int id) override
 Handles the removal of particles to the particle Handler. More...
 
void handleParticleAddition (unsigned int id, BaseParticle *p) override
 Handles the addition of particles to the particle Handler. More...
 
void retrieveVertexParticles ()
 Tries to get pointers to all vertex particles from the handler. More...
 
void checkActive ()
 Check if the triangle is considered active. More...
 
bool getActive ()
 
void setHandler (WallHandler *handler) override
 Set the handler. More...
 
void applyPressure (Mdouble presure)
 Apply a force pointing in normal direction corresponding to the specified pressure. More...
 
void applyForce (Vec3D force)
 Apply the given force to the triangle. More...
 
- 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
 
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
 
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...
 
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...
 
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 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
 

Public Attributes

std::array< MeshTriangle *, 3 > neighbor = {{nullptr}}
 
std::vector< std::vector< unsigned int > > vertexNeighbors
 

Private Member Functions

void updateVertexAndNormal ()
 Update vertexMin_, vertexMax_ and faceNormal_ for an updated position. More...
 
void updateVerticesFromParticles ()
 Retrieve new positions from updated vertex particles. More...
 

Private Attributes

std::array< Vec3D, 3 > vertex_
 
std::array< Vec3D, 3 > vertexVelocity_
 
std::array< BaseParticle *, 3 > vertexParticle_ = {{nullptr}}
 
Vec3D vertexMin_
 
Vec3D vertexMax_
 
std::array< Vec3D, 3 > edgeNormal_
 
std::array< Vec3D, 3 > edge_
 
std::array< unsigned int, 3 > vertexIds_
 
std::array< double, 3 > edgeLength_
 
Vec3D faceNormal_
 
Mdouble area_
 
Mdouble invMass_ = 0.0
 
bool isActive = 0
 

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

MeshTriangle implements a triangle whose vertex positions are defined by three particles.

It tracks the movement of the specified particles and updates its own potition every timestep. If neighboring objects along the edges and vertices are known, their contacts will be taken into account to consider if particle contacts should produce a force. Calculated contact forces will be transferred to the certex particles.

A triangle may e.g. be constructed with the following code.

p0 = Vec3D(0, 0.1, 1);
p1 = Vec3D(0, 0.1, 0.9);
p2 = Vec3D(0, 0 , 1);
p.setSpecies(someSpecies);
p.setRadius(2e-2);
p.setVelocity(Vec3D(0.0, 0.0, 0.0));
p.setHandler(&particleHandler);
unsigned int Id1 = particleHandler.copyAndAddObject(p)->getId();
unsigned int Id2 = particleHandler.copyAndAddObject(p)->getId();
unsigned int Id3 = particleHandler.copyAndAddObject(p)->getId();
f.setVertices(p0, p1, p2);
// Create all faces with their initial positions
f.setVertexIds(0, 1, 2);
f.setSpecies(someSpecies);
f.setMass(mass);
// Retrieve the correct positions from the particles
void setVelocity(const Vec3D &velocity)
set the velocity of the BaseInteractable.
Definition: BaseInteractable.cc:350
virtual void setPosition(const Vec3D &position)
Sets the position of this BaseInteractable.
Definition: BaseInteractable.h:239
void setHandler(ParticleHandler *handler)
Sets the pointer to the particle's ParticleHandler.
Definition: BaseParticle.cc:663
virtual void setRadius(Mdouble radius)
Sets the particle's radius_ (and adjusts the mass_ accordingly, based on the particle's species)
Definition: BaseParticle.cc:553
void setSpecies(const ParticleSpecies *species)
Definition: BaseParticle.cc:818
void setSpecies(const ParticleSpecies *species)
Defines the species of the current wall.
Definition: BaseWall.cc:169
MeshTriangle implements a triangle whose vertex positions are defined by three particles.
Definition: MeshTriangle.h:75
void setMass(Mdouble mass)
Definition: MeshTriangle.cc:681
void actionsAfterParticleGhostUpdate() override
actionsPerformed after the position update of (ghost-) particles.
Definition: MeshTriangle.cc:599
void setVertices(Vec3D A, Vec3D B, Vec3D C)
Sets member variables such that the wall represents a triangle with vertices A, B,...
Definition: MeshTriangle.cc:456
void setVertexIds(unsigned int i, unsigned int j, unsigned int k)
sets the ids of the vertex particles. Calls retrieveVertexParticles.
Definition: MeshTriangle.cc:481
A spherical particle is the most simple particle used in MercuryDPM.
Definition: SphericalParticle.h:37
Definition: Vector.h:51

For a longer example, please have a look at the class Membrane.

Constructor & Destructor Documentation

◆ MeshTriangle() [1/2]

MeshTriangle::MeshTriangle ( )
default

Default constructor.

Referenced by copy().

◆ MeshTriangle() [2/2]

MeshTriangle::MeshTriangle ( const MeshTriangle other)
default

Copy constructor.

◆ ~MeshTriangle()

MeshTriangle::~MeshTriangle ( )
overridedefault

Destructor.

Member Function Documentation

◆ actionsAfterParticleGhostUpdate()

void MeshTriangle::actionsAfterParticleGhostUpdate ( )
overridevirtual

actionsPerformed after the position update of (ghost-) particles.

After the particles get new positions, these need to be retrived to update the triangle.

Reimplemented from BaseWall.

600 {
602 }
void updateVerticesFromParticles()
Retrieve new positions from updated vertex particles.
Definition: MeshTriangle.cc:510

References updateVerticesFromParticles().

Referenced by Membrane::buildMesh().

◆ actionsOnRestart()

void MeshTriangle::actionsOnRestart ( )
overridevirtual

Actions executed on restart.

On restart, try to get all vertex particle pointers.

Reimplemented from BaseWall.

590 {
591  // Note this can not be done in the read sequence, as the particles are not yet available there
593 }
void retrieveVertexParticles()
Tries to get pointers to all vertex particles from the handler.
Definition: MeshTriangle.cc:563

References retrieveVertexParticles().

◆ applyForce()

void MeshTriangle::applyForce ( Vec3D  force)

Apply the given force to the triangle.

Parameters
[in]force

Applies a given force to the triangle, by splitting it between the vertex particles.

714 {
715  if (isActive)
716  {
717  force /= 3.0;
718  for (unsigned int j=0; j<3; j++){
719  vertexParticle_[j]->addForce(force);
720  }
721  }
722 
723 }
bool isActive
Definition: MeshTriangle.h:309
std::array< BaseParticle *, 3 > vertexParticle_
Definition: MeshTriangle.h:282

References isActive, and vertexParticle_.

◆ applyPressure()

void MeshTriangle::applyPressure ( Mdouble  pressure)

Apply a force pointing in normal direction corresponding to the specified pressure.

Parameters
[in]pressureThe pressure value in units of 1 Pa

Calculates the force acting on the triangle using the pressure and the triangles surface area. Applies the force to the vertex particles.

695 {
696  if (isActive)
697  {
698  // Area in Normaldirection
699  // Calculate F = p*A/3.0
700  // The division by 3.0 is to split the force evenly between the points
701  Vec3D pressureForce = getArea()*pressure*getFaceNormal()/3.0;
702  for (unsigned int j=0; j<3; j++){
703  vertexParticle_[j]->addForce(pressureForce);
704  }
705  }
706 }
Mdouble getArea() const
Returns the area of the triangle.
Definition: MeshTriangle.h:147
Vec3D getFaceNormal() const
Returns the face normal.
Definition: MeshTriangle.h:142

References getArea(), getFaceNormal(), isActive, and vertexParticle_.

◆ checkActive()

void MeshTriangle::checkActive ( )

Check if the triangle is considered active.

Checks if the triangle is active. A triangle is considered active, if pointers to all references are known.

582 {
584 }

References isActive, and vertexParticle_.

Referenced by handleParticleAddition(), read(), retrieveVertexParticles(), and updateVerticesFromParticles().

◆ checkInteractions()

void MeshTriangle::checkInteractions ( InteractionHandler interactionHandler,
unsigned int  timeStamp 
)
overridevirtual

Checks, if the forces of all interctions should be applied.

Parameters
[in]interactionHandlerPointer to InteractionHandler that contains all the interactions.
[in]timeStampThe time at which we want to look at the interactions.

Determine, if a contact os valid, i.e. that the forces due to this contact should be applied to both the particle in the wall. The evaluation is done by looking at the interaction a specific particle has with the current wall as well as neighboring walls. If a particle has multiple contacts, the selection criteria noted in doi:10.1002/nme.4487 are applied. Note: At this time this leads to issues when the particles are much bigger than the triangles.

Reimplemented from BaseWall.

111 {
112  unsigned j, id;
113  // Iterate through the interactions and check if there is one with higher priority
114  // Note: The id should already be taken into account when creating the interactions
115  for (j = 0; j<getInteractions().size(); j++)
116  {
117  bool interactionValid = true;
118  auto i = getInteractions()[j];
119  if (i->getTimeStamp()!=timeStamp)
120  {
121  continue; // Old interactions are not needed
122  logger(WARN, "Taking old interaction into account");
123  }
124  if (i->getMultiContactIdentifier()>3) // Vertex contact
125  {
126  id = i->getMultiContactIdentifier() - 4;
127  for (auto triangleId : vertexNeighbors[id])
128  {
129  // logger(INFO, "id % %", id, getId());
130  auto k = getHandler()->getObjectById(triangleId);
131  auto interactionNeighbor = interactionHandler->getExistingInteraction(i->getP(), k);
132  if (interactionNeighbor && interactionNeighbor->getTimeStamp()==timeStamp)
133  {
134  if (interactionNeighbor->getMultiContactIdentifier() <= 3)
135  {
136  // logger(INFO, "Found invalid vertex contact");
137  interactionValid = false;
138  break;
139  }
140  }
141  }
142  }
143  else if (i->getMultiContactIdentifier()>0) // Edge contact
144  {
145  id = i->getMultiContactIdentifier() - 1;
146  if (neighbor[id])
147  {
148  auto interactionNeighbor = interactionHandler->getExistingInteraction(i->getP(), neighbor[id]);
149  if (interactionNeighbor && interactionNeighbor->getTimeStamp()==timeStamp)
150  {
151  if (interactionNeighbor->getMultiContactIdentifier() == 0)
152  {
153  // logger(INFO, "Found invalid edge contact");
154  interactionValid = false;
155  }
156  }
157  }
158  }
159 
160  // Undo the interaction
161  if (!interactionValid)
162  {
163  i->getP()->addForce(-i->getForce());
164  this->addForce(i->getForce());
165  if (getHandler()->getDPMBase()->getRotation()) {
166  i->getP()->addTorque(-i->getTorque() + Vec3D::cross(i->getP()->getPosition() - i->getContactPoint(), i->getForce()));
167  this->addTorque(i->getTorque() - Vec3D::cross(this->getPosition() - i->getContactPoint(), i->getForce()));
168  }
169  i->setForce(Vec3D(0,0,0));
170  }
171  else
172  {
173  Vec3D weight = getBaricentricWeight(i->getContactPoint());
174 
175  // Contact seems valid, distribute its force to the edges
176  Vec3D force = -i->getForce();
177 
178  for(unsigned k=0; k<3; k++)
179  {
180  vertexParticle_[k]->addForce(weight.getComponent(k)*force);
181  }
182  }
183  }
184 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ WARN
T * getObjectById(const unsigned int id)
Gets a pointer to the Object at the specified index in the BaseHandler.
Definition: BaseHandler.h:565
void addTorque(const Vec3D &addTorque)
Adds an amount to the torque on this BaseInteractable.
Definition: BaseInteractable.cc:132
const std::vector< BaseInteraction * > & getInteractions() const
Returns a list of interactions which belong to this interactable.
Definition: BaseInteractable.h:277
void addForce(const Vec3D &addForce)
Adds an amount to the force on this BaseInteractable.
Definition: BaseInteractable.cc:116
WallHandler * getHandler() const
A function which returns the WallHandler that handles this BaseWall.
Definition: BaseWall.cc:134
BaseInteraction * getExistingInteraction(const BaseInteractable *P, const BaseInteractable *I) const
Returns the Interaction between the BaseInteractable's P and I if it exists, otherwise returns a null...
Definition: InteractionHandler.cc:103
std::vector< std::vector< unsigned int > > vertexNeighbors
Definition: MeshTriangle.h:262
std::array< MeshTriangle *, 3 > neighbor
Definition: MeshTriangle.h:257
const Vec3D getBaricentricWeight(const Vec3D &contact) const
Calculates the barycentric weight of a specified point.
Definition: MeshTriangle.cc:326
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:163
Mdouble getComponent(int index) const
Returns the requested component of this Vec3D.
Definition: Vector.cc:194
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51

References BaseInteractable::addForce(), BaseInteractable::addTorque(), Vec3D::cross(), getBaricentricWeight(), Vec3D::getComponent(), InteractionHandler::getExistingInteraction(), BaseWall::getHandler(), BaseInteractable::getInteractions(), BaseHandler< T >::getObjectById(), constants::i, logger, neighbor, vertexNeighbors, vertexParticle_, and WARN.

◆ copy()

MeshTriangle* MeshTriangle::copy ( ) const
inlineoverridevirtual

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

Implements BaseWall.

98  { return new MeshTriangle(*this); }
MeshTriangle()=default
Default constructor.

References MeshTriangle().

◆ getActive()

bool MeshTriangle::getActive ( )
inline
236 { return isActive; }

References isActive.

◆ getArea()

Mdouble MeshTriangle::getArea ( ) const
inline

Returns the area of the triangle.

147 {return area_;}
Mdouble area_
Definition: MeshTriangle.h:304

References area_.

Referenced by applyPressure(), Membrane::buildMesh(), and Membrane::initializeEdgeBendingQuantities().

◆ getBaricentricWeight()

const Vec3D MeshTriangle::getBaricentricWeight ( const Vec3D contact) const

Calculates the barycentric weight of a specified point.

Parameters
[in]contactThe point, at which the weights should be calculated.
Returns
A Vec3D giving the weights. (x is the weight of corner 0,...)

Calculates baricentric weight of a given point in the triangle

327 {
328  Vec3D m;
329  // The area of a triangle is
330  // Taken from https://gamedev.stackexchange.com/questions/23743/whats-the-most-efficient-way-to-find-barycentric-coordinates
331  // Mdouble areaABC = Vec3D::dot( getFaceNormal(), Vec3D::cross(vertex_[1]-vertex_[0], vertex_[2]-vertex_[0] ) ) ;
332  // logger.assert_always(std::fabs(areaABC-area_*2)<1e-8, "OOPS, triangle area calculared wrong % %", areaABC, area_*2);
333  Mdouble areaPBC = 0.5*Vec3D::dot( getFaceNormal(), Vec3D::cross(vertex_[1]-contact, vertex_[2]-contact ) ) ;
334  Mdouble areaPCA = 0.5*Vec3D::dot( getFaceNormal(), Vec3D::cross(vertex_[2]-contact, vertex_[0]-contact ) ) ;
335 
336  m.X = areaPBC / area_ ; // alpha
337  m.Y = areaPCA / area_ ; // beta
338  m.Z = 1.0f - m.X - m.Y ; // gamma
339  return m;
340 }
double Mdouble
Definition: GeneralDefine.h:34
std::array< Vec3D, 3 > vertex_
Definition: MeshTriangle.h:280
Mdouble Y
Definition: Vector.h:66
Mdouble Z
Definition: Vector.h:66
Mdouble X
the vector components
Definition: Vector.h:66
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:76

References area_, Vec3D::cross(), Vec3D::dot(), getFaceNormal(), vertex_, Vec3D::X, Vec3D::Y, and Vec3D::Z.

Referenced by checkInteractions(), getVelocityAtContact(), and isInsideTriangle().

◆ getDistanceAndNormal()

bool MeshTriangle::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 MeshTriangle. If there is a collision, this function also computes the distance between the BaseParticle and MeshTriangle and the normal of the MeshTriangle at the intersection point. It does this by calling MeshTriangle::getDistanceNormalOverlapType. 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.

200 {
201  if (!isActive) return false;
202  Mdouble overlap;
203  unsigned int type;
204  return getDistanceNormalOverlapType(p, distance, normal, overlap, type);
205 }
bool getDistanceNormalOverlapType(const BaseParticle &p, Mdouble &distance, Vec3D &normal, Mdouble &overlap, unsigned int &type) const
Definition: MeshTriangle.cc:225

References getDistanceNormalOverlapType(), and isActive.

◆ getDistanceNormalOverlapType()

bool MeshTriangle::getDistanceNormalOverlapType ( const BaseParticle p,
Mdouble distance,
Vec3D normal,
Mdouble overlap,
unsigned int &  type 
) const
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.
[out]overlap_returnIf there was a collision, the overlap will be placed here.
[out]type_returnIf there was a collision, the contact type 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 MeshTriangle. If there is a collision, this function also computes the distance between the BaseParticle and MeshTriangle and the normal of the MeshTriangle at the intersection point as well as the contact overlap and type. The type is set to the following: 0: face contact 1, 2, 3: contact with type-1th edge 4, 5, 6; contact with the type-4th vertex

226 {
227  if (!isActive) return false;
228 
229  // TODO: Note that this if may lead to contacts beeing ignored in md->checkParticleForInteraction because unadded
230  // particles have the Id 0, which might indeed be the Id of a node particle
231  if ( std::find(vertexIds_.begin(), vertexIds_.end(), p.getId()) != vertexIds_.end() )
232  {
233  // Do not detect any contact between particles that correspond to the walls nodes
234  return false;
235  }
236 
237  const Vec3D position = p.getPosition(); // note: no transfer to lab coordinates
238  const Mdouble distanceMax = p.getWallInteractionRadius(this);
239 
240  // compute distance from face
241  // get distance from particle position to the face
242  const Mdouble signedDistance = Vec3D::dot(position-vertex_[0], faceNormal_);
243  distance = std::abs(signedDistance);
244 
245  // check if any contact is possible
246  if (distance >= distanceMax) return false;
247 
248  // compute distance from edges
249  const std::array<Vec3D,3> edgeBranch {position - vertex_[0], position - vertex_[1], position - vertex_[2]};
250  const std::array<double,3> edgeDistance {Vec3D::dot(edgeBranch[0], edgeNormal_[0]), Vec3D::dot(edgeBranch[1], edgeNormal_[1]), Vec3D::dot(edgeBranch[2], edgeNormal_[2])};
251 
252  // find edge with largest distance (this will be the edge if there is a edge contact)
253  const size_t id = (edgeDistance[1] > edgeDistance[2]) ?
254  (edgeDistance[0] > edgeDistance[1] ? 0 : 1) : (edgeDistance[0] > edgeDistance[2] ? 0 : 2);
255 
256  // check if there will be contact
257  if (edgeDistance[id] > distanceMax) return false;
258 
259  // determine if it is a face contact
260  const Vec3D posProjected = position - signedDistance * faceNormal_;
261  if (edgeDistance[id] <= 0 && isInsideTriangle(posProjected)){
262  normal = (signedDistance >= 0) ? -faceNormal_ : faceNormal_;
263  overlap = p.getRadius() - distance;
264  type=0; // Face contact
265  return true;
266  }
267 
268  // Then the neighbor will handle this interaction
269  if (neighbor[id] && neighbor[id]->getId()<this->getId()) return false;
270  // determine if it is an edge or vertex contact
271  const double positionAlongEdge = Vec3D::dot(edgeBranch[id], edge_[id]);
272  if (positionAlongEdge <= 0) {
273  //possible contact with left vertex
274  distance = edgeBranch[id].getLength();
275  if (distance > distanceMax) return false;
276  // check vertex ids
277  for (auto triangleId: vertexNeighbors[id])
278  {
279  if (triangleId<this->getId()) return false;
280  }
281  normal = edgeBranch[id] / -distance;
282  type = 4 + id; // Vertex contact
283  } else if (positionAlongEdge >= edgeLength_[id]) {
284  //contact with right vertex
285  const size_t idRight = (id + 1) % 3;
286  distance = edgeBranch[idRight].getLength();
287  if (distance > distanceMax) return false;
288  // check vertex ids
289  for (auto triangleId: vertexNeighbors[idRight])
290  {
291  if (triangleId<this->getId()) return false;
292  }
293  type = 4 + idRight; // Vertex contact
294  normal = edgeBranch[idRight] / -distance;
295 
296  } else {
297  // edge contact
298  normal = edge_[id] * positionAlongEdge - edgeBranch[id];
299  distance = normal.getLength();
300  if (distance > distanceMax) return false;
301  normal /= distance;
302  type = 1 + id;
303  }
304 
305  overlap = p.getRadius() - distance;
306  return true;
307 }
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:218
unsigned int getId() const
Returns the unique identifier of any particular object.
Definition: BaseObject.h:125
Mdouble getRadius() const
Returns the particle's radius.
Definition: BaseParticle.h:348
Mdouble getWallInteractionRadius(const BaseWall *wall) const
returns the radius plus the interactionDistance
Definition: BaseParticle.h:386
Vec3D faceNormal_
Definition: MeshTriangle.h:303
std::array< unsigned int, 3 > vertexIds_
Definition: MeshTriangle.h:296
std::array< double, 3 > edgeLength_
Definition: MeshTriangle.h:297
bool isInsideTriangle(const Vec3D &point) const
Determines if a given point is within the triangle.
Definition: MeshTriangle.cc:660
std::array< Vec3D, 3 > edgeNormal_
Definition: MeshTriangle.h:293
std::array< Vec3D, 3 > edge_
Definition: MeshTriangle.h:294
static Mdouble getLength(const Vec3D &a)
Calculates the length of a Vec3D: .
Definition: Vector.cc:331

References Vec3D::dot(), edge_, edgeLength_, edgeNormal_, faceNormal_, BaseObject::getId(), Vec3D::getLength(), BaseInteractable::getPosition(), BaseParticle::getRadius(), BaseParticle::getWallInteractionRadius(), isActive, isInsideTriangle(), neighbor, vertex_, vertexIds_, and vertexNeighbors.

Referenced by getDistanceAndNormal(), and getInteractionWith().

◆ getFaceNormal()

Vec3D MeshTriangle::getFaceNormal ( ) const
inline

Returns the face normal.

142 {return faceNormal_;}

References faceNormal_.

Referenced by applyPressure(), and getBaricentricWeight().

◆ getInteractionWith()

BaseInteraction * MeshTriangle::getInteractionWith ( BaseParticle p,
unsigned  timeStamp,
InteractionHandler interactionHandler 
)
overridevirtual
Parameters
[in]pPointer to the BaseParticle which we want to check the interaction for.
[in]timeStampThe time at which we want to look at the interaction.
[in]interactionHandlerA pointer to the InteractionHandler in which the interaction can be found.
Returns
A pointer to the BaseInteraction that happened between this BaseWall and the BaseParticle at the timeStamp.

It is determined, if there is any contact between the particle and the triangular wall. If a contact exists, all neccesary quantities are determined and set. This includes the contact type is determined (corner, edge or face contact, see laos doi:10.1002/nme.4487) and set with setMultiContactIdentifier.

Todo:

Reimplemented from BaseWall.

49 {
50  if (!isActive) return nullptr;
51 
52  // TODO: Note that this if may lead to contacts beeing ignored in md->checkParticleForInteraction because unadded
53  // particles have the Id 0, which might indeed be the Id of an node particle
54  if ( std::find(vertexIds_.begin(), vertexIds_.end(), p->getId()) != vertexIds_.end() )
55  {
56  // Do not detect any contact between particles that correspond to the walls nodes
57  return nullptr;
58  }
59 
60  Mdouble distance;
61  Vec3D normal;
62  Mdouble overlap;
63  unsigned int type;
64 
65  if (!(p->isSphericalParticle()))
66  {
67  logger(ERROR, "MeshTriangle::getInteractionWith not implemented for particles of type %", p->getName());
68  }
69 
70  if (getDistanceNormalOverlapType(*p, distance, normal, overlap, type))
71  {
72  // look for an existing interaction, or create a new one
73  BaseInteraction *c = nullptr;
74  if (getGroupId() > 0 && p->getInteractions().size() > 0)
75  {
76  // Do not care for the result in that case.
77  return BaseWall::getInteractionWith(p, timeStamp, interactionHandler);
78  }
79 
80  // look for an existing interaction, or create a new one
81  c = interactionHandler->getInteraction(p, this, timeStamp);
82 
83  c->setNormal(-normal);
84  c->setDistance(distance);
85  c->setOverlap(overlap);
87  c->setWallInteraction(1);
88 
90  c->setContactPoint(p->getPosition() - (p->getRadius() - 0.5 * c->getOverlap()) * c->getNormal());
91 
92  logger(DEBUG, "Particle contact with wall at %", c->getContactPoint());
93  return c;
94  }
95  return nullptr;
96 }
@ ERROR
@ DEBUG
Stores information about interactions between two interactable objects; often particles but could be ...
Definition: BaseInteraction.h:60
void setDistance(Mdouble distance)
Sets the interaction distance between the two interacting objects.
Definition: BaseInteraction.cc:221
const Vec3D & getContactPoint() const
Gets constant reference to contact point (vector).
Definition: BaseInteraction.h:234
void setContactPoint(Vec3D contactPoint)
Set the location of the contact point between the two interacting objects.
Definition: BaseInteraction.cc:240
void setMultiContactIdentifier(unsigned int multiContactIdentifier_)
Definition: BaseInteraction.cc:768
const Vec3D & getNormal() const
Gets the normal vector between the two interacting objects.
Definition: BaseInteraction.h:226
void setWallInteraction(bool flag)
Definition: BaseInteraction.cc:950
void setNormal(Vec3D normal)
Sets the normal vector between the two interacting objects.
Definition: BaseInteraction.cc:212
Mdouble getOverlap() const
Returns a Mdouble with the current overlap between the two interacting objects.
Definition: BaseInteraction.h:240
void setOverlap(Mdouble overlap)
Set the overlap between the two interacting object.
Definition: BaseInteraction.cc:231
unsigned getGroupId() const
Definition: BaseObject.h:137
virtual bool isSphericalParticle() const =0
std::string getName() const override
Returns the name of the object.
Definition: BaseParticle.cc:348
BaseInteraction * getInteractionWith(BaseParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler) override
Returns the interaction between this wall and a given particle, nullptr if there is no interaction.
Definition: BaseWall.cc:367
BaseInteraction * getInteraction(BaseInteractable *P, BaseInteractable *I, unsigned timeStamp)
Returns the Interaction between the BaseInteractable's P and I.
Definition: InteractionHandler.cc:146

References DEBUG, ERROR, BaseInteraction::getContactPoint(), getDistanceNormalOverlapType(), BaseObject::getGroupId(), BaseObject::getId(), InteractionHandler::getInteraction(), BaseInteractable::getInteractions(), BaseWall::getInteractionWith(), BaseParticle::getName(), BaseInteraction::getNormal(), BaseInteraction::getOverlap(), BaseInteractable::getPosition(), BaseParticle::getRadius(), isActive, BaseParticle::isSphericalParticle(), logger, BaseInteraction::setContactPoint(), BaseInteraction::setDistance(), BaseInteraction::setMultiContactIdentifier(), BaseInteraction::setNormal(), BaseInteraction::setOverlap(), BaseInteraction::setWallInteraction(), and vertexIds_.

◆ getInvMass()

Mdouble MeshTriangle::getInvMass ( ) const
overridevirtual

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

Reimplemented from BaseInteractable.

672 {
673  return invMass_;
674 }
Mdouble invMass_
Definition: MeshTriangle.h:307

References invMass_.

◆ getName()

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

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

Implements BaseObject.

104  { return "MeshTriangle"; }

◆ getVelocityAtContact()

const Vec3D MeshTriangle::getVelocityAtContact ( const Vec3D contact) const
overridevirtual

Calculates the local velocity at a specified point.

Parameters
[in]contactThe point, at which the velocity should be determined.
Returns
A Vec3D giving the calculated velocity.

Calculates the velocity at the contact position by interpolating the velocity of the triangle nodes using barycentric coordinates.

Reimplemented from BaseInteractable.

316 {
317  Vec3D m = getBaricentricWeight(contact);
318  return m.x()*vertexVelocity_[0] + m.y()*vertexVelocity_[1] + m.z()*vertexVelocity_[2];
319 };
std::array< Vec3D, 3 > vertexVelocity_
Definition: MeshTriangle.h:281
Mdouble & y()
RW reference to Y.
Definition: Vector.h:372
Mdouble & z()
RW reference to Z.
Definition: Vector.h:384
Mdouble & x()
RW reference to X.
Definition: Vector.h:360

References getBaricentricWeight(), vertexVelocity_, Vec3D::x(), Vec3D::y(), and Vec3D::z().

◆ getVertexIds()

std::array<unsigned int, 3> MeshTriangle::getVertexIds ( ) const
inline

Returns an array containing the ids of the vertex particles.

160 { return vertexIds_; }

References vertexIds_.

Referenced by Membrane::initializeEdgeBendingQuantities(), and Membrane::updateFaceNeighbors().

◆ getVertices()

std::array<Vec3D,3> MeshTriangle::getVertices ( ) const
inline

Returns an array of the vertex coordinates.

137 {return vertex_;}

References vertex_.

◆ handleParticleAddition()

void MeshTriangle::handleParticleAddition ( unsigned int  id,
BaseParticle p 
)
overridevirtual

Handles the addition of particles to the particle Handler.

Parameters
[in]idThe id of the added particle
[in]pPointer to the particle

If thie given id is equal to one of the vertex Particles, the reference is added. If this added particle makes the triangle active, the vertex positions are uodated using updateVerticesFromParticles().

Reimplemented from BaseWall.

546 {
547  unsigned int i;
548  for (i=0; i<3; i++)
549  {
550  if (vertexIds_[i] == id)
551  {
552  vertexParticle_[i] = p;
553  checkActive();
555  }
556  }
557 }
void checkActive()
Check if the triangle is considered active.
Definition: MeshTriangle.cc:581

References checkActive(), constants::i, updateVerticesFromParticles(), vertexIds_, and vertexParticle_.

◆ handleParticleRemoval()

void MeshTriangle::handleParticleRemoval ( unsigned int  id)
overridevirtual

Handles the removal of particles to the particle Handler.

Parameters
[in]idThe id of the removed particle

If the given id is equal to one of the vertexParticles, the reference to that particle is removed and the triangle is set inactive

Reimplemented from BaseWall.

526 {
527  unsigned int i;
528  for (i=0; i<3; i++)
529  {
530  if (vertexIds_[i] == id)
531  {
532  vertexParticle_[i] = nullptr;
533  isActive = false;
534  }
535  }
536 }

References constants::i, isActive, vertexIds_, and vertexParticle_.

◆ isInsideTriangle()

bool MeshTriangle::isInsideTriangle ( const Vec3D point) const

Determines if a given point is within the triangle.

Parameters
[in]pointPosition to check
Returns
boolean specifying if point is within the triangle
661 {
662  Vec3D weights = getBaricentricWeight(point);
663 
664  Mdouble eps = 1e-12;
665  return ((1-eps) > weights.X > eps && (1-eps) > weights.Y > eps && (1-eps) > weights.Z > eps) && ((1-eps) > weights.X+weights.Y > eps && (1-eps) > weights.Y+weights.Z > eps && (1-eps) > weights.Z+weights.X > eps);
666  // return ((1-eps) > s0 > eps && (1-eps) > s1 > eps && (1-eps) > s2 > eps) && ((1-eps) > s0+s1 > eps && (1-eps) > s1+s2 > eps && (1-eps) > s2+s0 > eps);
667 
668 }

References getBaricentricWeight(), Vec3D::X, Vec3D::Y, and Vec3D::Z.

Referenced by getDistanceNormalOverlapType().

◆ isLocal()

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

Determines if the triangle is considered local.

Parameters
[out]minthe minimum of the triangle in all directions.
[out]maxthe maximum of the triangle in all directions
Returns
true The triangle is local

Reimplemented from BaseWall.

650 {
651  min = vertexMin_;
652  max = vertexMax_;
653  return true;
654 }
Vec3D vertexMin_
Definition: MeshTriangle.h:287
Vec3D vertexMax_
Definition: MeshTriangle.h:288

References vertexMax_, and vertexMin_.

◆ move()

void MeshTriangle::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.

611 {
614 }
virtual void move(const Vec3D &move)
Moves this BaseInteractable by adding an amount to the position.
Definition: BaseInteractable.cc:215
void updateVertexAndNormal()
Update vertexMin_, vertexMax_ and faceNormal_ for an updated position.
Definition: MeshTriangle.cc:622
void move(const Vec3D &move) override
Definition: MeshTriangle.cc:610

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

◆ read()

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

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

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

Reimplemented from BaseWall.

355 {
356  BaseWall::read(is);
357 
358  unsigned int i, j, s1, s2, Id;
359  std::string dummy;
360 
361  is >> dummy >> s1;
362  vertexNeighbors.reserve(s1);
363  for (i=0; i<s1; i++)
364  {
365  is >> dummy >> s2;
366  std::vector<unsigned int> vN;
367  vN.reserve(s2);
368 
369  for (j=0; j<s2; j++)
370  {
371  is >> Id;
372  vN.push_back(Id);
373  }
374  vertexNeighbors.push_back(vN);
375 
376  }
377 
378 
379 
380  is >> dummy;
381  for (int i = 0; i < 3; i++)
382  {
383  is >> vertex_[i];
384  }
385  is >> dummy;
386  for (int i = 0; i < 3; i++)
387  {
388  is >> vertexIds_[i];
389  }
390 
391  is >> dummy >> invMass_;
392  is >> dummy >> isActive;
393 
394  checkActive();
395  // updateVerticesFromParticles();
397 }
void read(std::istream &is) override
Function that reads a BaseWall from an input stream, usually a restart file.
Definition: BaseWall.cc:78

References checkActive(), constants::i, invMass_, isActive, BaseWall::read(), updateVertexAndNormal(), vertex_, vertexIds_, and vertexNeighbors.

◆ retrieveVertexParticles()

void MeshTriangle::retrieveVertexParticles ( )

Tries to get pointers to all vertex particles from the handler.

Tries to get the pointer to the vertex particles from the particleHandler. Afterwards it checks if the triangle is active and updates the vertex positions.

564 {
565  if (getHandler())
566  {
567  unsigned int i;
568  for (i=0; i<3; i++)
569  {
571  }
572  checkActive();
574  }
575 }
DPMBase * getDPMBase()
Gets the problem that is solved using this handler.
Definition: BaseHandler.h:725
ParticleHandler particleHandler
An object of the class ParticleHandler, contains the pointers to all the particles created.
Definition: DPMBase.h:1437

References checkActive(), BaseHandler< T >::getDPMBase(), BaseWall::getHandler(), BaseHandler< T >::getObjectById(), constants::i, DPMBase::particleHandler, updateVerticesFromParticles(), vertexIds_, and vertexParticle_.

Referenced by actionsOnRestart(), setHandler(), and setVertexIds().

◆ rotate()

void MeshTriangle::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.

343 {
344  if (!angularVelocityDt.isZero())
345  {
346  BaseInteractable::rotate(angularVelocityDt);
348  }
349 }
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().

◆ setHandler()

void MeshTriangle::setHandler ( WallHandler handler)
overridevirtual

Set the handler.

Parameters
[in]handlerPointer to the wallHandler.

Sets the wall handler and calls retrieveVertexParticles

Reimplemented from BaseWall.

494 {
495  BaseWall::setHandler(handler);
496 
497  if (handler->getDPMBase()->particleHandler.getNumberOfObjects()==0)
498  {
499  // Without this, restart does not work
500  return;
501  }
502 
504 }
virtual void setHandler(WallHandler *handler)
A function which sets the WallHandler for this BaseWall.
Definition: BaseWall.cc:125
unsigned int getNumberOfObjects() const override
Returns the number of objects in the container. In parallel code this practice is forbidden to avoid ...
Definition: ParticleHandler.cc:1325

References BaseHandler< T >::getDPMBase(), ParticleHandler::getNumberOfObjects(), DPMBase::particleHandler, retrieveVertexParticles(), and BaseWall::setHandler().

◆ setMass()

void MeshTriangle::setMass ( Mdouble  mass)
Parameters
[in]massValue of the mass assigned to the triangle.

The mass is neccesary for contact force determination. If not given, infinite mass is assumed.

682 {
683  logger.assert_always(mass > 0.0,
684  "Error in MeshTriangle::setMass, the given mass to be set must be positive.");
685  invMass_ = 1.0 / mass;
686 
687 }

References invMass_, and logger.

Referenced by Membrane::buildMesh().

◆ setVertexIds()

void MeshTriangle::setVertexIds ( unsigned int  i,
unsigned int  j,
unsigned int  k 
)

sets the ids of the vertex particles. Calls retrieveVertexParticles.

Parameters
[in]i,j,kThe ids of the particles representing the corners.
482 {
483  vertexIds_[0] = i;
484  vertexIds_[1] = j;
485  vertexIds_[2] = k;
487 }

References constants::i, retrieveVertexParticles(), and vertexIds_.

Referenced by Membrane::buildMesh().

◆ setVertexVelocities()

void MeshTriangle::setVertexVelocities ( Vec3D  A,
Vec3D  B,
Vec3D  C 
)

Sets the velocity of the vertex points.

472 {
473  vertexVelocity_[0] = A;
474  vertexVelocity_[1] = B;
475  vertexVelocity_[2] = C;
476 }
@ A
Definition: StatisticsVector.h:42

References A, and vertexVelocity_.

Referenced by updateVerticesFromParticles().

◆ setVertices() [1/2]

void MeshTriangle::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
  • updateVertexAndNormal is called to set the remaining variables
457 {
458  setVertices(A, B, C, (A + B + C) / 3);
459 }

References A.

Referenced by Membrane::buildMesh(), and updateVerticesFromParticles().

◆ setVertices() [2/2]

void MeshTriangle::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.

462 {
463  setPosition(position);
464  setOrientation({1, 0, 0, 0});
465  vertex_[0] = A;
466  vertex_[1] = B;
467  vertex_[2] = C;
469 }
virtual void setOrientation(const Quaternion &orientation)
Sets the orientation of this BaseInteractable.
Definition: BaseInteractable.h:260

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

◆ updateVertexAndNormal()

void MeshTriangle::updateVertexAndNormal ( )
private

Update vertexMin_, vertexMax_ and faceNormal_ for an updated position.

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)
623 {
626 
627  edge_ = {vertex_[1] - vertex_[0], vertex_[2] - vertex_[1], vertex_[0] - vertex_[2]};
629  area_ = 0.5*faceNormal_.getLength();
630  logger.assert_always(0.5*sqrt(Vec3D::cross(vertex_[1] - vertex_[0], vertex_[2] - vertex_[1]).getLengthSquared())==area_, "OOPS, face area wrong");
631 
633 
634  for (int i = 0; i < 3; i++)
635  {
637  edgeLength_[i] = edge_[i].getLength();
638  edge_[i] /= edgeLength_[i];
639  edgeNormal_[i].normalise();
640  }
641  //logger(INFO,"vertex %,%,% edge %,%,% face %",vertex_[0],vertex_[1],vertex_[2],edgeNormal_[0],edgeNormal_[1],edgeNormal_[2],faceNormal_);
642 }
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 area_, Vec3D::cross(), edge_, edgeLength_, edgeNormal_, faceNormal_, Vec3D::getLength(), constants::i, logger, Vec3D::max(), Vec3D::min(), Vec3D::normalise(), vertex_, vertexMax_, and vertexMin_.

Referenced by move(), read(), rotate(), and setVertices().

◆ updateVerticesFromParticles()

void MeshTriangle::updateVerticesFromParticles ( )
private

Retrieve new positions from updated vertex particles.

Update the triangle position and velocity based on the vertex particles, if the triangle is active.

511 {
512  // Need to get references to the particles
513  checkActive();
514  if (!isActive) return;
515 
518 }
virtual const Vec3D & getVelocity() const
Returns the velocity of this interactable.
Definition: BaseInteractable.cc:329
void setVertexVelocities(Vec3D A, Vec3D B, Vec3D C)
Sets the velocity of the vertex points.
Definition: MeshTriangle.cc:471

References checkActive(), BaseInteractable::getPosition(), BaseInteractable::getVelocity(), isActive, setVertexVelocities(), setVertices(), and vertexParticle_.

Referenced by actionsAfterParticleGhostUpdate(), handleParticleAddition(), and retrieveVertexParticles().

◆ write()

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

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

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

Reimplemented from BaseWall.

404 {
405  BaseWall::write(os);
406  unsigned int i, j, s1, s2;
407  s1 = vertexNeighbors.size();
408  os << " vertexNeighbors " << s1;
409  for (i=0; i<s1; i++)
410  {
411  s2 = vertexNeighbors[i].size();
412  os << " vertexNeighborsElements " << s2;
413  for (j=0; j<s2; j++)
414  {
415  os << " " << vertexNeighbors[i][j];
416  }
417  }
418 
419  os << " vertex ";
420  for (int i = 0; i < 3; i++)
421  {
422  os << ' ' << vertex_[i];
423  }
424  os << " edgeParticleIds ";
425  for (int i = 0; i < 3; i++)
426  {
427  os << ' ' << vertexIds_[i];
428  }
429 
430  os << " invMass " << invMass_;
431 
432  os << " isActive " << isActive;
433 
434  // Note: need not to write faceNormal_ and area_, as these are recalculated
435  // on read.
436 }
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, invMass_, isActive, vertex_, vertexIds_, vertexNeighbors, and BaseWall::write().

◆ writeVTK()

void MeshTriangle::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.

439 {
440  if (!isActive) return;
441 
442  const unsigned long s = vtk.points.size();
443  for (auto v : vertex_)
444  {
445  vtk.points.push_back(v);
446  }
447  std::vector<double> cell;
448  cell.reserve(3);
449  cell.push_back(s);
450  cell.push_back(s + 1);
451  cell.push_back(s + 2);
452 
453  vtk.triangleStrips.push_back(cell);
454 }
std::vector< std::vector< double > > triangleStrips
Definition: BaseWall.h:40
std::vector< Vec3D > points
Definition: BaseWall.h:39

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

Member Data Documentation

◆ area_

Mdouble MeshTriangle::area_
private

◆ edge_

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

◆ edgeLength_

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

◆ edgeNormal_

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

stores the wall normal n in n.x=p

Referenced by getDistanceNormalOverlapType(), and updateVertexAndNormal().

◆ faceNormal_

Vec3D MeshTriangle::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 getDistanceNormalOverlapType(), getFaceNormal(), and updateVertexAndNormal().

◆ invMass_

Mdouble MeshTriangle::invMass_ = 0.0
private

Referenced by getInvMass(), read(), setMass(), and write().

◆ isActive

◆ neighbor

std::array<MeshTriangle*, 3> MeshTriangle::neighbor = {{nullptr}}

stores references to the neighbors along the edges.

Referenced by checkInteractions(), getDistanceNormalOverlapType(), and Membrane::updateFaceNeighbors().

◆ vertex_

std::array<Vec3D, 3> MeshTriangle::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 getBaricentricWeight(), getDistanceNormalOverlapType(), getVertices(), read(), setVertices(), updateVertexAndNormal(), write(), and writeVTK().

◆ vertexIds_

◆ vertexMax_

Vec3D MeshTriangle::vertexMax_
private

Referenced by isLocal(), and updateVertexAndNormal().

◆ vertexMin_

Vec3D MeshTriangle::vertexMin_
private

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

Referenced by isLocal(), and updateVertexAndNormal().

◆ vertexNeighbors

std::vector<std::vector<unsigned int> > MeshTriangle::vertexNeighbors

stores references to the neighbors on corners.

Referenced by checkInteractions(), getDistanceNormalOverlapType(), read(), Membrane::updateFaceNeighbors(), and write().

◆ vertexParticle_

◆ vertexVelocity_

std::array<Vec3D, 3> MeshTriangle::vertexVelocity_
private

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