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);
p.setPosition(p0);
unsigned int Id1 = particleHandler.copyAndAddObject(p)->getId();
p.setPosition(p1);
unsigned int Id2 = particleHandler.copyAndAddObject(p)->getId();
p.setPosition(p2);
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
f.actionsAfterParticleGhostUpdate();
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Vector3f p0
Definition: MatrixBase_all.cpp:2
Vector3f p1
Definition: MatrixBase_all.cpp:2
float * p
Definition: Tutorial_Map_using.cpp:9
MeshTriangle implements a triangle whose vertex positions are defined by three particles.
Definition: MeshTriangle.h:54
A spherical particle is the most simple particle used in MercuryDPM.
Definition: SphericalParticle.h:16
Definition: Kernel/Math/Vector.h:30
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237

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.

587 {
589 }
void updateVerticesFromParticles()
Retrieve new positions from updated vertex particles.
Definition: MeshTriangle.cc:497

References updateVerticesFromParticles().

◆ actionsOnRestart()

void MeshTriangle::actionsOnRestart ( )
overridevirtual

Actions executed on restart.

On restart, try to get all vertex particle pointers.

Reimplemented from BaseWall.

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

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.

701 {
702  if (isActive)
703  {
704  force /= 3.0;
705  for (unsigned int j=0; j<3; j++){
706  vertexParticle_[j]->addForce(force);
707  }
708  }
709 
710 }
bool isActive
Definition: MeshTriangle.h:288
std::array< BaseParticle *, 3 > vertexParticle_
Definition: MeshTriangle.h:261
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References isActive, j, 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.

682 {
683  if (isActive)
684  {
685  // Area in Normaldirection
686  // Calculate F = p*A/3.0
687  // The division by 3.0 is to split the force evenly between the points
688  Vec3D pressureForce = getArea()*pressure*getFaceNormal()/3.0;
689  for (unsigned int j=0; j<3; j++){
690  vertexParticle_[j]->addForce(pressureForce);
691  }
692  }
693 }
Mdouble getArea() const
Returns the area of the triangle.
Definition: MeshTriangle.h:126
Vec3D getFaceNormal() const
Returns the face normal.
Definition: MeshTriangle.h:121

References getArea(), getFaceNormal(), isActive, j, 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.

569 {
571 }

References isActive, and vertexParticle_.

Referenced by handleParticleAddition(), 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.

90 {
91  unsigned j, id;
92  // Iterate through the interactions and check if there is one with higher priority
93  // Note: The id should already be taken into account when creating the interactions
94  for (j = 0; j<getInteractions().size(); j++)
95  {
96  bool interactionValid = true;
97  auto i = getInteractions()[j];
98  if (i->getTimeStamp()!=timeStamp)
99  {
100  continue; // Old interactions are not needed
101  logger(WARN, "Taking old interaction into account");
102  }
103  if (i->getMultiContactIdentifier()>3) // Vertex contact
104  {
105  id = i->getMultiContactIdentifier() - 4;
106 
107  for (auto triangleId : vertexNeighbors[id])
108  {
109  // If a neighboring particle has a contact and a lower Id with
110  // the same particle, disregard this contact
111  if (triangleId<this->getId())
112  {
113  interactionValid = false;
114  break;
115  }
116  // logger(INFO, "id % %", id, getId());
117  auto k = getHandler()->getObjectById(triangleId);
118  auto interactionNeighbor = interactionHandler->getExistingInteraction(i->getP(), k);
119  if (interactionNeighbor && interactionNeighbor->getTimeStamp()==timeStamp)
120  {
121  if (interactionNeighbor->getMultiContactIdentifier() <= 3)
122  {
123  // logger(INFO, "Found invalid vertex contact");
124  interactionValid = false;
125  break;
126  }
127  }
128  }
129  }
130  else if (i->getMultiContactIdentifier()>0) // Edge contact
131  {
132  id = i->getMultiContactIdentifier() - 1;
133  if (neighbor[id])
134  {
135  // If a neighboring particle has at least an edge contact with
136  // the same particle and a lower Id, disregard this contact
137  if (neighbor[id]->getId()<this->getId())
138  {
139  // logger(INFO, "Found invalid Edge contact due to Id");
140  interactionValid = false;
141  }
142  else
143  {
144  auto interactionNeighbor = interactionHandler->getExistingInteraction(i->getP(), neighbor[id]);
145  if (interactionNeighbor && interactionNeighbor->getTimeStamp()==timeStamp)
146  {
147  if (interactionNeighbor->getMultiContactIdentifier() == 0)
148  {
149  // logger(INFO, "Found invalid edge contact");
150  interactionValid = false;
151  }
152  }
153  }
154  }
155  }
156 
157  // Undo the interaction
158  if (!interactionValid)
159  {
160  i->getP()->addForce(-i->getForce());
161  this->addForce(i->getForce());
162  if (getHandler()->getDPMBase()->getRotation())
163  {
164  i->getP()->addTorque(-i->getTorque() + Vec3D::cross(i->getP()->getPosition() - i->getContactPoint(), i->getForce()));
165  this->addTorque(i->getTorque() - Vec3D::cross(this->getPosition() - i->getContactPoint(), i->getForce()));
166  }
167  i->setForce(Vec3D(0,0,0));
168  }
169  else
170  {
171  Vec3D weight = getBaricentricWeight(i->getContactPoint());
172 
173  // Contact seems valid, distribute its force to the edges
174  Vec3D force = -i->getForce();
175 
176  for(unsigned k=0; k<3; k++)
177  {
178  vertexParticle_[k]->addForce(weight.getComponent(k)*force);
179  }
180  }
181  }
182 }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
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:573
void addTorque(const Vec3D &addTorque)
Adds an amount to the torque on this BaseInteractable.
Definition: BaseInteractable.cc:110
const std::vector< BaseInteraction * > & getInteractions() const
Returns a list of interactions which belong to this interactable.
Definition: BaseInteractable.h:256
void addForce(const Vec3D &addForce)
Adds an amount to the force on this BaseInteractable.
Definition: BaseInteractable.cc:94
unsigned int getId() const
Returns the unique identifier of any particular object.
Definition: BaseObject.h:104
WallHandler * getHandler() const
A function which returns the WallHandler that handles this BaseWall.
Definition: BaseWall.cc:113
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:83
std::vector< std::vector< unsigned int > > vertexNeighbors
Definition: MeshTriangle.h:241
std::array< MeshTriangle *, 3 > neighbor
Definition: MeshTriangle.h:236
const Vec3D getBaricentricWeight(const Vec3D &contact) const
Calculates the barycentric weight of a specified point.
Definition: MeshTriangle.cc:315
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:143
Mdouble getComponent(int index) const
Returns the requested component of this Vec3D.
Definition: Vector.cc:174
char char char int int * k
Definition: level2_impl.h:374

References BaseInteractable::addForce(), BaseInteractable::addTorque(), Vec3D::cross(), getBaricentricWeight(), Vec3D::getComponent(), InteractionHandler::getExistingInteraction(), BaseWall::getHandler(), BaseObject::getId(), BaseInteractable::getInteractions(), BaseHandler< T >::getObjectById(), i, j, k, 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.

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

References MeshTriangle().

◆ getActive()

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

References isActive.

◆ getArea()

Mdouble MeshTriangle::getArea ( ) const
inline

Returns the area of the triangle.

126 {return area_;}
Mdouble area_
Definition: MeshTriangle.h:283

References area_.

Referenced by applyPressure().

◆ 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

316 {
317  Vec3D m;
318  // The area of a triangle is
319  // Taken from https://gamedev.stackexchange.com/questions/23743/whats-the-most-efficient-way-to-find-barycentric-coordinates
320  // Mdouble areaABC = Vec3D::dot( getFaceNormal(), Vec3D::cross(vertex_[1]-vertex_[0], vertex_[2]-vertex_[0] ) ) ;
321  // logger.assert_always(std::fabs(areaABC-area_*2)<1e-8, "OOPS, triangle area calculared wrong % %", areaABC, area_*2);
322  Mdouble areaPBC = 0.5*Vec3D::dot( getFaceNormal(), Vec3D::cross(vertex_[1]-contact, vertex_[2]-contact ) ) ;
323  Mdouble areaPCA = 0.5*Vec3D::dot( getFaceNormal(), Vec3D::cross(vertex_[2]-contact, vertex_[0]-contact ) ) ;
324 
325  m.X = areaPBC / area_ ; // alpha
326  m.Y = areaPCA / area_ ; // beta
327  m.Z = 1.0f - m.X - m.Y ; // gamma
328  return m;
329 }
std::array< Vec3D, 3 > vertex_
Definition: MeshTriangle.h:259
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:56
int * m
Definition: level2_cplx_impl.h:294

References area_, Vec3D::cross(), Vec3D::dot(), getFaceNormal(), m, and vertex_.

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.

198 {
199  if (!isActive) return false;
200  Mdouble overlap;
201  unsigned int type;
202  return getDistanceNormalOverlapType(p, distance, normal, overlap, type);
203 }
bool getDistanceNormalOverlapType(const BaseParticle &p, Mdouble &distance, Vec3D &normal, Mdouble &overlap, unsigned int &type) const
Definition: MeshTriangle.cc:223
void normal(const Vector< double > &x, Vector< double > &normal)
Definition: free_surface_rotation.cc:65
type
Definition: compute_granudrum_aor.py:141

References getDistanceNormalOverlapType(), isActive, WallFunction::normal(), p, and compute_granudrum_aor::type.

◆ 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

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

References abs(), Vec3D::dot(), edge_, edgeLength_, edgeNormal_, faceNormal_, isActive, isInsideTriangle(), WallFunction::normal(), p, compute_granudrum_aor::type, vertex_, and vertexIds_.

Referenced by getDistanceAndNormal(), and getInteractionWith().

◆ getFaceNormal()

Vec3D MeshTriangle::getFaceNormal ( ) const
inline

Returns the face normal.

121 {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.

28 {
29  if (!isActive) return nullptr;
30 
31  // TODO: Note that this if may lead to contacts beeing ignored in md->checkParticleForInteraction because unadded
32  // particles have the Id 0, which might indeed be the Id of an node particle
33  if ( std::find(vertexIds_.begin(), vertexIds_.end(), p->getId()) != vertexIds_.end() )
34  {
35  // Do not detect any contact between particles that correspond to the walls nodes
36  return nullptr;
37  }
38 
39  Mdouble distance;
40  Vec3D normal;
41  Mdouble overlap;
42  unsigned int type;
43 
44  if (!(p->isSphericalParticle()))
45  {
46  logger(ERROR, "MeshTriangle::getInteractionWith not implemented for particles of type %", p->getName());
47  }
48 
49  if (getDistanceNormalOverlapType(*p, distance, normal, overlap, type))
50  {
51  // look for an existing interaction, or create a new one
52  BaseInteraction *c = nullptr;
53  if (getGroupId() > 0 && p->getInteractions().size() > 0)
54  {
55  // Do not care for the result in that case.
56  return BaseWall::getInteractionWith(p, timeStamp, interactionHandler);
57  }
58 
59  // look for an existing interaction, or create a new one
60  c = interactionHandler->getInteraction(p, this, timeStamp);
61 
62  c->setNormal(-normal);
63  c->setDistance(distance);
64  c->setOverlap(overlap);
65  c->setMultiContactIdentifier(type);
66  c->setWallInteraction(1);
67 
69  c->setContactPoint(p->getPosition() - (p->getRadius() - 0.5 * c->getOverlap()) * c->getNormal());
70 
71  logger(DEBUG, "Particle contact with wall at %", c->getContactPoint());
72  return c;
73  }
74  return nullptr;
75 }
@ ERROR
@ DEBUG
Stores information about interactions between two interactable objects; often particles but could be ...
Definition: BaseInteraction.h:39
unsigned getGroupId() const
Definition: BaseObject.h:116
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:346
BaseInteraction * getInteraction(BaseInteractable *P, BaseInteractable *I, unsigned timeStamp)
Returns the Interaction between the BaseInteractable's P and I.
Definition: InteractionHandler.cc:126
int c
Definition: calibrate.py:100

References calibrate::c, DEBUG, ERROR, getDistanceNormalOverlapType(), BaseObject::getGroupId(), InteractionHandler::getInteraction(), BaseWall::getInteractionWith(), isActive, logger, WallFunction::normal(), p, compute_granudrum_aor::type, 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.

659 {
660  return invMass_;
661 }
Mdouble invMass_
Definition: MeshTriangle.h:286

References invMass_.

◆ getName()

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

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

Implements BaseObject.

83  { 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.

305 {
306  Vec3D m = getBaricentricWeight(contact);
307  return m.x()*vertexVelocity_[0] + m.y()*vertexVelocity_[1] + m.z()*vertexVelocity_[2];
308 };
std::array< Vec3D, 3 > vertexVelocity_
Definition: MeshTriangle.h:260

References getBaricentricWeight(), m, and vertexVelocity_.

◆ getVertexIds()

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

Returns an array containing the ids of the vertex particles.

139 { return vertexIds_; }

References vertexIds_.

◆ getVertices()

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

Returns an array of the vertex coordinates.

116 {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.

533 {
534  unsigned int i;
535  for (i=0; i<3; i++)
536  {
537  if (vertexIds_[i] == id)
538  {
539  vertexParticle_[i] = p;
540  checkActive();
542  }
543  }
544 }
void checkActive()
Check if the triangle is considered active.
Definition: MeshTriangle.cc:568

References checkActive(), i, p, 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.

513 {
514  unsigned int i;
515  for (i=0; i<3; i++)
516  {
517  if (vertexIds_[i] == id)
518  {
519  vertexParticle_[i] = nullptr;
520  isActive = false;
521  }
522  }
523 }

References 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
648 {
650 
651  Mdouble eps = 1e-12;
652  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);
653  // 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);
654 
655 }
double eps
Definition: crbond_bessel.cc:24
list weights
Definition: calibrate.py:94

References e(), CRBond_Bessel::eps, getBaricentricWeight(), and calibrate::weights.

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.

637 {
638  min = vertexMin_;
639  max = vertexMax_;
640  return true;
641 }
Vec3D vertexMin_
Definition: MeshTriangle.h:266
Vec3D vertexMax_
Definition: MeshTriangle.h:267
#define min(a, b)
Definition: datatypes.h:22
#define max(a, b)
Definition: datatypes.h:23

References max, min, 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.

598 {
601 }
virtual void move(const Vec3D &move)
Moves this BaseInteractable by adding an amount to the position.
Definition: BaseInteractable.cc:193
void updateVertexAndNormal()
Update vertexMin_, vertexMax_ and faceNormal_ for an updated position.
Definition: MeshTriangle.cc:609
void move(const Vec3D &move) override
Definition: MeshTriangle.cc:597

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.

344 {
345  BaseWall::read(is);
346 
347  unsigned int i, j, s1, s2, Id;
348  std::string dummy;
349 
350  is >> dummy >> s1;
351  vertexNeighbors.reserve(s1);
352  for (i=0; i<s1; i++)
353  {
354  is >> dummy >> s2;
355  std::vector<unsigned int> vN;
356  vN.reserve(s2);
357 
358  for (j=0; j<s2; j++)
359  {
360  is >> Id;
361  vN.push_back(Id);
362  }
363  vertexNeighbors.push_back(vN);
364 
365  }
366 
367 
368 
369  is >> dummy;
370  for (int i = 0; i < 3; i++)
371  {
372  is >> vertex_[i];
373  }
374  is >> dummy;
375  for (int i = 0; i < 3; i++)
376  {
377  is >> vertexIds_[i];
378  }
379 
380  is >> dummy >> invMass_;
381  is >> dummy >> isActive;
382 
383  // updateVerticesFromParticles();
384 }
void read(std::istream &is) override
Function that reads a BaseWall from an input stream, usually a restart file.
Definition: BaseWall.cc:57
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286

References i, invMass_, isActive, j, BaseWall::read(), oomph::Global_string_for_annotation::string(), 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.

551 {
552  if (getHandler())
553  {
554  unsigned int i;
555  for (i=0; i<3; i++)
556  {
558  }
559  checkActive();
561  }
562 }
DPMBase * getDPMBase()
Gets the problem that is solved using this handler.
Definition: BaseHandler.h:733
ParticleHandler particleHandler
An object of the class ParticleHandler, contains the pointers to all the particles created.
Definition: DPMBase.h:1443

References checkActive(), BaseHandler< T >::getDPMBase(), BaseWall::getHandler(), BaseHandler< T >::getObjectById(), 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.

332 {
333  if (!angularVelocityDt.isZero())
334  {
335  BaseInteractable::rotate(angularVelocityDt);
337  }
338 }
virtual void rotate(const Vec3D &angularVelocityDt)
Rotates this BaseInteractable.
Definition: BaseInteractable.cc:208
bool isZero() const
Checks if ALL elements are zero.
Definition: Kernel/Math/Vector.h:94

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.

481 {
482  BaseWall::setHandler(handler);
483 
484  if (handler->getDPMBase()->particleHandler.getNumberOfObjects()==0)
485  {
486  // Without this, restart does not work
487  return;
488  }
489 
491 }
virtual void setHandler(WallHandler *handler)
A function which sets the WallHandler for this BaseWall.
Definition: BaseWall.cc:104
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:1323

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.

669 {
670  logger.assert_always(mass > 0.0,
671  "Error in MeshTriangle::setMass, the given mass to be set must be positive.");
672  invMass_ = 1.0 / mass;
673 
674 }

References invMass_, and logger.

◆ 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.
469 {
470  vertexIds_[0] = i;
471  vertexIds_[1] = j;
472  vertexIds_[2] = k;
474 }

References i, j, k, retrieveVertexParticles(), and vertexIds_.

◆ setVertexVelocities()

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

Sets the velocity of the vertex points.

459 {
460  vertexVelocity_[0] = A;
461  vertexVelocity_[1] = B;
462  vertexVelocity_[2] = C;
463 }
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
Matrix< Scalar, Dynamic, Dynamic > C
Definition: bench_gemm.cpp:49
Matrix< SCALARB, Dynamic, Dynamic, opt_B > B
Definition: bench_gemm.cpp:48

References 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
444 {
445  setVertices(A, B, C, (A + B + C) / 3);
446 }
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
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:443
Definition: matrices.h:74

Referenced by 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.

449 {
450  setPosition(position);
451  setOrientation({1, 0, 0, 0});
452  vertex_[0] = A;
453  vertex_[1] = B;
454  vertex_[2] = C;
456 }
virtual void setOrientation(const Quaternion &orientation)
Sets the orientation of this BaseInteractable.
Definition: BaseInteractable.h:239
virtual void setPosition(const Vec3D &position)
Sets the position of this BaseInteractable.
Definition: BaseInteractable.h:218

References 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)
610 {
613 
614  edge_ = {vertex_[1] - vertex_[0], vertex_[2] - vertex_[1], vertex_[0] - vertex_[2]};
616  area_ = 0.5*faceNormal_.getLength();
617  logger.assert_always(0.5*sqrt(Vec3D::cross(vertex_[1] - vertex_[0], vertex_[2] - vertex_[1]).getLengthSquared())==area_, "OOPS, face area wrong");
618 
620 
621  for (int i = 0; i < 3; i++)
622  {
624  edgeLength_[i] = edge_[i].getLength();
625  edge_[i] /= edgeLength_[i];
626  edgeNormal_[i].normalise();
627  }
628  //logger(INFO,"vertex %,%,% edge %,%,% face %",vertex_[0],vertex_[1],vertex_[2],edgeNormal_[0],edgeNormal_[1],edgeNormal_[2],faceNormal_);
629 }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
static Vec3D max(const Vec3D &a, const Vec3D &b)
Calculates the pointwise maximum of two Vec3D.
Definition: Vector.cc:69
static Vec3D min(const Vec3D &a, const Vec3D &b)
Calculates the pointwise minimum of two Vec3D.
Definition: Vector.cc:82
static Mdouble getLength(const Vec3D &a)
Calculates the length of a Vec3D: .
Definition: Vector.cc:350
void normalise()
Makes this Vec3D unit length.
Definition: Vector.cc:103

References area_, Vec3D::cross(), edge_, edgeLength_, edgeNormal_, faceNormal_, Vec3D::getLength(), i, logger, Vec3D::max(), Vec3D::min(), Vec3D::normalise(), sqrt(), vertex_, vertexMax_, and vertexMin_.

Referenced by move(), 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.

498 {
499  // Need to get references to the particles
500  checkActive();
501  if (!isActive) return;
502 
505 }
virtual const Vec3D & getVelocity() const
Returns the velocity of this interactable.
Definition: BaseInteractable.cc:307
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:197
void setVertexVelocities(Vec3D A, Vec3D B, Vec3D C)
Sets the velocity of the vertex points.
Definition: MeshTriangle.cc:458

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.

391 {
392  BaseWall::write(os);
393  unsigned int i, j, s1, s2;
394  s1 = vertexNeighbors.size();
395  os << " vertexNeighbors " << s1;
396  for (i=0; i<s1; i++)
397  {
398  s2 = vertexNeighbors[i].size();
399  os << " vertexNeighborsElements " << s2;
400  for (j=0; j<s2; j++)
401  {
402  os << " " << vertexNeighbors[i][j];
403  }
404  }
405 
406  os << " vertex ";
407  for (int i = 0; i < 3; i++)
408  {
409  os << ' ' << vertex_[i];
410  }
411  os << " edgeParticleIds ";
412  for (int i = 0; i < 3; i++)
413  {
414  os << ' ' << vertexIds_[i];
415  }
416 
417  os << " invMass " << invMass_;
418 
419  os << " isActive " << isActive;
420 
421  // Note: need not to write faceNormal_ and area_, as these are recalculated
422  // on read.
423 }
void write(std::ostream &os) const override
Function that writes a BaseWall to an output stream, usually a restart file.
Definition: BaseWall.cc:79

References i, invMass_, isActive, j, 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.

426 {
427  if (!isActive) return;
428 
429  const unsigned long s = vtk.points.size();
430  for (auto v : vertex_)
431  {
432  vtk.points.push_back(v);
433  }
434  std::vector<double> cell;
435  cell.reserve(3);
436  cell.push_back(s);
437  cell.push_back(s + 1);
438  cell.push_back(s + 2);
439 
440  vtk.triangleStrips.push_back(cell);
441 }
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
RealScalar s
Definition: level1_cplx_impl.h:130
std::vector< std::vector< double > > triangleStrips
Definition: BaseWall.h:19
std::vector< Vec3D > points
Definition: BaseWall.h:18

References isActive, VTKContainer::points, s, VTKContainer::triangleStrips, v, 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().

◆ 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(), read(), and write().

◆ vertexParticle_

◆ vertexVelocity_

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

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