WearableTriangleMeshWall Class Reference

#include <WearableTriangleMeshWall.h>

+ Inheritance diagram for WearableTriangleMeshWall:

Public Member Functions

 WearableTriangleMeshWall ()=default
 Default constructor. More...
 
 WearableTriangleMeshWall (const std::vector< Vec3D > &points, const std::vector< std::array< unsigned, 3 >> &cells, const ParticleSpecies *species=nullptr)
 Constructor creating triangles. More...
 
 WearableTriangleMeshWall (const Vec3D &P0, const Vec3D &P1, const Vec3D &P, Mdouble resolutionU, Mdouble resolutionV, const ParticleSpecies *species=nullptr, bool periodicInU=false, bool periodicInV=false)
 Constructor creating parallelogram shaped mesh with a given resolution in u- and v-direction. More...
 
 WearableTriangleMeshWall (const WearableTriangleMeshWall &other)
 Copy constructor. More...
 
 ~WearableTriangleMeshWall ()=default
 Destructor. More...
 
WearableTriangleMeshWalloperator= (const WearableTriangleMeshWall &other)
 Copy assignment operator. More...
 
WearableTriangleMeshWallcopy () const override
 Wall copy method. More...
 
void read (std::istream &is) override
 Reads a WearableTriangleMeshWall from an input stream, for example a restart file. More...
 
void write (std::ostream &os) const override
 Writes a WearableTriangleMeshWall to an output stream, for example a restart file. More...
 
std::string getName () const override
 Returns the name of the object. More...
 
void computeWear () override
 
void setWearCoefficient (Mdouble wearCoefficient)
 
void setHardness (Mdouble hardness)
 
void setWearAcceleration (Mdouble wearAcceleration)
 
- Public Member Functions inherited from TriangleMeshWall
 TriangleMeshWall ()=default
 Default constructor. More...
 
 TriangleMeshWall (const TriangleMeshWall &other)
 Copy constructor. More...
 
 TriangleMeshWall (const std::vector< Vec3D > &points, const std::vector< std::array< unsigned, 3 >> &cells, const ParticleSpecies *species=nullptr)
 Constructor creating triangles. More...
 
 TriangleMeshWall (const Vec3D &P0, const Vec3D &P1, const Vec3D &P, Mdouble resolutionU, Mdouble resolutionV, const ParticleSpecies *species=nullptr, bool periodicInU=false, bool periodicInV=false)
 Constructor creating parallelogram shaped mesh with a given resolution in u- and v-direction. More...
 
 ~TriangleMeshWall () override=default
 Destructor. More...
 
TriangleMeshWalloperator= (const TriangleMeshWall &other)
 Copy assignment operator. More...
 
void set (const std::vector< Vec3D > &points, const std::vector< std::array< unsigned, 3 >> &cells)
 Set function creating triangles. More...
 
bool getDistanceAndNormal (const BaseParticle &p, Mdouble &distance, Vec3D &normal_return) const override
 Checks if a collision is happening with any of the TriangleWalls. NOTE: this does not handle actual interactions, but can only be used to know IF one (of possible many) interaction occurs. More...
 
BaseInteractiongetInteractionWith (BaseParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler) override
 Gets the interaction between a given BaseParticle and all TriangleWalls in this mesh at a given time. In case of multiple contacts, chooses the one with greatest overlap. More...
 
void resetForceTorque (int numberOfOMPthreads) override
 
void writeVTK (VTKContainer &vtk) const override
 
void setSpecies (const ParticleSpecies *species)
 
void moveVertex (unsigned index, const Vec3D &dP)
 Updates a vertex by a given change in position. More...
 
void removeTriangle (unsigned index, bool removeFreeVertex=true)
 Removes a single triangle from the mesh. More...
 
void refineTriangle (unsigned index, unsigned numberOfTimes=1)
 Makes a really basic refinement by adding a vertex in the middle of the triangle to create 3 new ones. Note: Multiple refinements will definitely NOT create a nice mesh. More...
 
void setPeriodicCompanions (const std::vector< std::pair< unsigned, unsigned >> &periodicCompanions)
 Sets the vector with pairs of indices of vertices which lie on a periodic boundary and are each others 'companion'. E.g. when moving one vertex, the other is also moved the same amount. A vertex can have multiple companions (when lying on multiple periodic boundaries) and those should be in the vector as multiple pairs. More...
 
void moveVerticesToMatchVolume (std::vector< Vec3D > displacements, Mdouble targetVolume, int maxNumRecursiveCalls=15)
 Moves all vertices by given displacements and corrects to match the change in volume with the target volume. The final displacements will be proportional to the initial ones, but their actual values will have change. More...
 
Mdouble getVolumeTetrahedron (const Vec3D &a, const Vec3D &b, const Vec3D &c, const Vec3D &d)
 Calculates the volume of a tetrahedron formed by 4 vertices. More...
 
void addToMesh (TriangleMeshWall mesh)
 Adds a mesh to this mesh. When vertices share a position, only the one in the current mesh is kept and both meshes are then "connected" at that point. More...
 
void createParallelogramMesh (const Vec3D &P0, const Vec3D &P1, const Vec3D &P2, Mdouble resolutionU, Mdouble resolutionV, bool periodicInU=false, bool periodicInV=false)
 Creates a parallelogram shaped mesh with a given resolution in u- and v-direction. More...
 
void createFourPointMesh (const Vec3D &P0, const Vec3D &P1, const Vec3D &P2, const Vec3D &P3, int numSegmentsU, int numSegmentsV)
 Creates mesh consisting of four points, with linearly interpolated segments. The points therefore don't necessarily have to lie in plane. More...
 
bool isLocal (Vec3D &min, Vec3D &max) const override
 
void setPosition (const Vec3D &position) override
 Sets the position of this BaseInteractable. More...
 
void setOrientation (const Quaternion &orientation) override
 Sets the orientation of this BaseInteractable. More...
 
void move (const Vec3D &move) override
 Moves this BaseInteractable by adding an amount to the position. More...
 
void rotate (const Vec3D &angularVelocity) override
 Rotates this BaseInteractable. 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
 
virtual void setHandler (WallHandler *handler)
 A function which sets the WallHandler for this BaseWall. More...
 
WallHandlergetHandler () const
 A function which returns the WallHandler that handles this BaseWall. More...
 
void setIndSpecies (unsigned int indSpecies) override
 Define the species of this wall using the index of the species in the SpeciesHandler in this DPMBase. More...
 
void setSpecies (const ParticleSpecies *species)
 Defines the species of the current wall. More...
 
bool isFixed () const override
 
void setForceControl (Vec3D forceGoal, Vec3D gainFactor, Vec3D baseVelocity={0, 0, 0})
 Slowly adjusts the force on a wall towards a specified goal, by adjusting (prescribing) the velocity of the wall. More...
 
bool getLinePlaneIntersect (Vec3D &intersect, const Vec3D &p0, const Vec3D &p1, const Vec3D &n, const Vec3D &p)
 
bool isInsideWallVTK (const Vec3D &point, const Vec3D &normal, const Vec3D &position) const
 
void projectOntoWallVTK (Vec3D &point0, const Vec3D &point1, const Vec3D &normal, const Vec3D &position) const
 
void intersectVTK (std::vector< Vec3D > &points, Vec3D normal, Vec3D position) const
 
virtual BaseInteractiongetInteractionWithSuperQuad (SuperQuadricParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler)
 
void getVTK (std::vector< Vec3D > &points, std::vector< std::vector< double >> &triangleStrips)
 
const Vec3D getAxis () const
 
virtual void actionsOnRestart ()
 No implementation but can be overidden in its derived classes. More...
 
virtual void actionsAfterParticleGhostUpdate ()
 No implementation but can be overidden in its derived classes. More...
 
virtual void handleParticleAddition (unsigned int id, BaseParticle *p)
 Handles the addition of particles to the particleHandler. More...
 
virtual void handleParticleRemoval (unsigned int id)
 Handles the addition of particles to the particleHandler. More...
 
virtual void checkInteractions (InteractionHandler *interactionHandler, unsigned int timeStamp)
 Check if all interactions are valid. More...
 
bool getVTKVisibility () const
 
void setVTKVisibility (bool vtkVisibility)
 
void addRenderedWall (BaseWall *w)
 
BaseWallgetRenderedWall (size_t i) const
 
std::vector< BaseWall * > getRenderedWalls () const
 
void removeRenderedWalls ()
 
void renderWall (VTKContainer &vtk)
 
void addParticlesAtWall (unsigned numElements=50)
 
void setVelocityControl (Vec3D forceGoal, Vec3D gainFactor, Vec3D baseVelocity)
 
virtual void writeWallDetailsVTK (VTKData &data) const
 
- 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...
 
void sumForceTorqueOMP ()
 
const Vec3DgetPosition () const
 Returns the position of this BaseInteractable. More...
 
const QuaterniongetOrientation () const
 Returns the orientation of this BaseInteractable. More...
 
void setOrientationViaNormal (Vec3D normal)
 Sets the orientation of this BaseInteractable by defining the vector that results from the rotation of the (1,0,0) vector. More...
 
void setOrientationViaEuler (Vec3D eulerAngle)
 Sets the orientation of this BaseInteractable by defining the euler angles. More...
 
const std::vector< BaseInteraction * > & getInteractions () const
 Returns a list of interactions which belong to this interactable. More...
 
void addInteraction (BaseInteraction *I)
 Adds an interaction to this BaseInteractable. More...
 
bool removeInteraction (BaseInteraction *I)
 Removes an interaction from this BaseInteractable. More...
 
void copyInteractionsForPeriodicParticles (const BaseInteractable &p)
 Copies interactions to this BaseInteractable whenever a periodic copy made. More...
 
void setVelocity (const Vec3D &velocity)
 set the velocity of the BaseInteractable. More...
 
void setAngularVelocity (const Vec3D &angularVelocity)
 set the angular velocity of the BaseInteractble. More...
 
void addVelocity (const Vec3D &velocity)
 adds an increment to the velocity. More...
 
void addAngularVelocity (const Vec3D &angularVelocity)
 add an increment to the angular velocity. More...
 
virtual const Vec3DgetVelocity () const
 Returns the velocity of this interactable. More...
 
virtual const Vec3DgetAngularVelocity () const
 Returns the angular velocity of this interactable. More...
 
void setPrescribedPosition (const std::function< Vec3D(double)> &prescribedPosition)
 Allows the position of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedPosition (double time)
 Computes the position from the user defined prescribed position function. More...
 
void setPrescribedVelocity (const std::function< Vec3D(double)> &prescribedVelocity)
 Allows the velocity of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedVelocity (double time)
 Computes the velocity from the user defined prescribed velocity function. More...
 
void setPrescribedOrientation (const std::function< Quaternion(double)> &prescribedOrientation)
 Allows the orientation of the infinite mass interactbale to be prescribed. More...
 
void applyPrescribedOrientation (double time)
 Computes the orientation from the user defined prescribed orientation function. More...
 
void setPrescribedAngularVelocity (const std::function< Vec3D(double)> &prescribedAngularVelocity)
 Allows the angular velocity of the infinite mass interactable to be prescribed. More...
 
void applyPrescribedAngularVelocity (double time)
 Computes the angular velocity from the user defined prescribed angular velocity. More...
 
virtual const Vec3D getVelocityAtContact (const Vec3D &contact) const
 Returns the velocity at the contact point, use by many force laws. More...
 
void integrateBeforeForceComputation (double time, double timeStep)
 This is part of integrate routine for objects with infinite mass. More...
 
void integrateAfterForceComputation (double time, double timeStep)
 This is part of the integration routine for objects with infinite mass. More...
 
virtual Mdouble getInvMass () const
 
virtual Mdouble getCurvature (const Vec3D &labFixedCoordinates) const
 
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
 

Private Member Functions

void storeDebris (const Triangle &triangle, const Vec3D &position, const Vec3D &debris, std::vector< Vec3D > &debrisContainer)
 Proportionally assigns the debris located at a certain position on a triangle to the triangle vertices. More...
 

Private Attributes

Mdouble wearCoefficient_
 Dimensionless wear coefficient. More...
 
Mdouble hardness_
 Hardness. More...
 
Mdouble wearAcceleration_
 Accelerates the wear process, to reduce simulation time needed to get results. More...
 

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...
 
- Protected Attributes inherited from TriangleMeshWall
std::vector< Triangletriangles_
 
std::vector< Vertexvertices_
 

Constructor & Destructor Documentation

◆ WearableTriangleMeshWall() [1/4]

WearableTriangleMeshWall::WearableTriangleMeshWall ( )
default

Default constructor.

Referenced by copy().

◆ WearableTriangleMeshWall() [2/4]

WearableTriangleMeshWall::WearableTriangleMeshWall ( const std::vector< Vec3D > &  points,
const std::vector< std::array< unsigned, 3 >> &  cells,
const ParticleSpecies species = nullptr 
)

Constructor creating triangles.

Parameters
pointsVector of vertices
cellsVector of cells, each cell holds three indices to the points vector
species
33  : TriangleMeshWall(points, cells, species)
34 {
35 }
TriangleMeshWall()=default
Default constructor.

◆ WearableTriangleMeshWall() [3/4]

WearableTriangleMeshWall::WearableTriangleMeshWall ( const Vec3D P0,
const Vec3D P1,
const Vec3D P,
Mdouble  resolutionU,
Mdouble  resolutionV,
const ParticleSpecies species = nullptr,
bool  periodicInU = false,
bool  periodicInV = false 
)

Constructor creating parallelogram shaped mesh with a given resolution in u- and v-direction.

Parameters
P0,P1,P2The three points, in clockwise order, forming the parallelogram (fourth point is implied).
resolutionUMaximum length of a segment in u-direction.
resolutionVMaximum length of a segment in v-direction.
species
periodicInUWhether or not the first and last column of vertices lie in a periodic boundary and should be set as periodic companions.
periodicInVWhether or not the first and last row of vertices lie in a periodic boundary and should be set as periodic companions.
39  :
40  TriangleMeshWall(P0, P1, P2, resolutionU, resolutionV, species, periodicInU, periodicInV)
41 {
42 }

◆ WearableTriangleMeshWall() [4/4]

WearableTriangleMeshWall::WearableTriangleMeshWall ( const WearableTriangleMeshWall other)

Copy constructor.

Parameters
otherThe WearableTriangleMeshWall that must be copied.
44  : TriangleMeshWall(other)
45 {
47  hardness_ = other.hardness_;
49 }
Mdouble wearCoefficient_
Dimensionless wear coefficient.
Definition: WearableTriangleMeshWall.h:128
Mdouble hardness_
Hardness.
Definition: WearableTriangleMeshWall.h:130
Mdouble wearAcceleration_
Accelerates the wear process, to reduce simulation time needed to get results.
Definition: WearableTriangleMeshWall.h:132

References hardness_, wearAcceleration_, and wearCoefficient_.

◆ ~WearableTriangleMeshWall()

WearableTriangleMeshWall::~WearableTriangleMeshWall ( )
default

Destructor.

Member Function Documentation

◆ computeWear()

void WearableTriangleMeshWall::computeWear ( )
overridevirtual

Reimplemented from BaseWall.

86 {
87  // Archard wear model, from https://en.wikipedia.org/wiki/Archard_equation
88  // Q = KWL/H
89  // Q is the total volume of wear debris produced
90  // K is a dimensionless constant (typically mild wear 1e-8, severe wear 1e-2)
91  // W is the total normal load
92  // L is the sliding distance
93  // H is the hardness of the softest contacting surfaces
94 
95  // Initialize debris vector with all zero Vec3Ds.
96  std::vector<Vec3D> debrisContainer(vertices_.size(), Vec3D(0.0, 0.0, 0.0));
97  Mdouble totalDebris = 0.0;
98 
99  for (auto& t : triangles_)
100  {
101  for (BaseInteraction* interaction : t.wall.getInteractions())
102  {
103  // Ignore old interactions and interactions from periodic particles.
104  if (interaction->getTimeStamp() <= getHandler()->getDPMBase()->getNumberOfTimeSteps() ||
105  static_cast<BaseParticle*>(interaction->getP())->getPeriodicFromParticle() != nullptr)
106  continue;
107 
108  const Mdouble W = interaction->getAbsoluteNormalForce();
109  // Pythagoras to get tangential magnitude from the normal magnitude and the relative velocity vector.
110  const Mdouble tangentialRelativeVelocity = std::sqrt(interaction->getRelativeVelocity().getLengthSquared() -
111  interaction->getNormalRelativeVelocity() * interaction->getNormalRelativeVelocity());
112  const Mdouble L = tangentialRelativeVelocity * getHandler()->getDPMBase()->getTimeStep();
113  const Mdouble Q = wearAcceleration_ * wearCoefficient_ * W * L / hardness_;
114 
115  totalDebris += Q;
116  // The contact point in normal direction is halfway the overlap. However, to properly store the debris the
117  // position exactly on the wall should be used.
118  Vec3D contactPointOnWall = interaction->getContactPoint() + 0.5 * interaction->getOverlap() * interaction->getNormal();
119  // Debris is removed in the direction of the interaction normal (rotated back to lab frame).
120  storeDebris(t, getOrientation().rotateBack(contactPointOnWall - getPosition()), getOrientation().rotateBack(interaction->getNormal()) * -Q, debrisContainer);
121  }
122  }
123 
124  moveVerticesToMatchVolume(debrisContainer, totalDebris);
125 }
double Mdouble
Definition: GeneralDefine.h:34
DPMBase * getDPMBase()
Gets the problem that is solved using this handler.
Definition: BaseHandler.h:725
const Quaternion & getOrientation() const
Returns the orientation of this BaseInteractable.
Definition: BaseInteractable.h:230
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:218
Stores information about interactions between two interactable objects; often particles but could be ...
Definition: BaseInteraction.h:60
Definition: BaseParticle.h:54
BaseParticle * getPeriodicFromParticle() const
Returns the 'original' particle this one's a periodic copy of.
Definition: BaseParticle.h:341
WallHandler * getHandler() const
A function which returns the WallHandler that handles this BaseWall.
Definition: BaseWall.cc:134
Mdouble getTimeStep() const
Returns the simulation time step.
Definition: DPMBase.cc:1250
std::vector< Triangle > triangles_
Definition: TriangleMeshWall.h:258
void moveVerticesToMatchVolume(std::vector< Vec3D > displacements, Mdouble targetVolume, int maxNumRecursiveCalls=15)
Moves all vertices by given displacements and corrects to match the change in volume with the target ...
Definition: TriangleMeshWall.cc:553
std::vector< Vertex > vertices_
Definition: TriangleMeshWall.h:259
Definition: Vector.h:51
void storeDebris(const Triangle &triangle, const Vec3D &position, const Vec3D &debris, std::vector< Vec3D > &debrisContainer)
Proportionally assigns the debris located at a certain position on a triangle to the triangle vertice...
Definition: WearableTriangleMeshWall.cc:127

References BaseHandler< T >::getDPMBase(), BaseWall::getHandler(), BaseInteractable::getOrientation(), BaseParticle::getPeriodicFromParticle(), BaseInteractable::getPosition(), DPMBase::getTimeStep(), hardness_, TriangleMeshWall::moveVerticesToMatchVolume(), storeDebris(), TriangleMeshWall::triangles_, TriangleMeshWall::vertices_, wearAcceleration_, and wearCoefficient_.

◆ copy()

WearableTriangleMeshWall * WearableTriangleMeshWall::copy ( ) const
overridevirtual

Wall copy method.

Returns
Pointer to the copy.

Reimplemented from TriangleMeshWall.

61 {
62  return new WearableTriangleMeshWall(*this);
63 }
WearableTriangleMeshWall()=default
Default constructor.

References WearableTriangleMeshWall().

Referenced by operator=().

◆ getName()

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

Returns the name of the object.

Returns
"WearableTriangleMeshWall"

Reimplemented from TriangleMeshWall.

81 {
82  return "WearableTriangleMeshWall";
83 }

◆ operator=()

WearableTriangleMeshWall & WearableTriangleMeshWall::operator= ( const WearableTriangleMeshWall other)

Copy assignment operator.

Parameters
otherThe WearableTriangleMeshWall that must be copied.
Returns
52 {
53  if (this == &other)
54  {
55  return *this;
56  }
57  return *(other.copy());
58 }
WearableTriangleMeshWall * copy() const override
Wall copy method.
Definition: WearableTriangleMeshWall.cc:60

References copy().

◆ read()

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

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

Reimplemented from TriangleMeshWall.

66 {
68 
69  std::string dummy;
70  is >> dummy >> wearCoefficient_ >> dummy >> hardness_ >> dummy >> wearAcceleration_;
71 }
void read(std::istream &is) override
Reads a TriangleMeshWall from an input stream, for example a restart file.
Definition: TriangleMeshWall.cc:73

References hardness_, TriangleMeshWall::read(), wearAcceleration_, and wearCoefficient_.

◆ setHardness()

void WearableTriangleMeshWall::setHardness ( Mdouble  hardness)
Parameters
hardnessHardness of the wall.
178 {
179  hardness_ = hardness;
180 }

References hardness_.

◆ setWearAcceleration()

void WearableTriangleMeshWall::setWearAcceleration ( Mdouble  wearAcceleration)
Parameters
wearAccelerationValue to accelerate the wear process with, to reduce simulation time needed to get results.
183 {
184  wearAcceleration_ = wearAcceleration;
185 }

References wearAcceleration_.

◆ setWearCoefficient()

void WearableTriangleMeshWall::setWearCoefficient ( Mdouble  wearCoefficient)
Parameters
wearCoefficientNon-dimensionless wear constant.
173 {
174  wearCoefficient_ = wearCoefficient;
175 }

References wearCoefficient_.

◆ storeDebris()

void WearableTriangleMeshWall::storeDebris ( const Triangle triangle,
const Vec3D position,
const Vec3D debris,
std::vector< Vec3D > &  debrisContainer 
)
private

Proportionally assigns the debris located at a certain position on a triangle to the triangle vertices.

Parameters
triangleTriangle on which the point lies.
positionPoint of debris.
debrisThe debris, which has a direction, i.e. the vector length is be the actual debris volume.
debrisContainerVector storing the debris for all vertices (index wise related).
128 {
129  const unsigned index0 = triangle.vertexIndices[0];
130  const unsigned index1 = triangle.vertexIndices[1];
131  const unsigned index2 = triangle.vertexIndices[2];
132 
133 // // Assign 1/3 to each vertex
134 // debrisContainer[index0] += debris / 3.0;
135 // debrisContainer[index1] += debris / 3.0;
136 // debrisContainer[index2] += debris / 3.0;
137 
138 
139 // // Distance from contact point to each vertex
140 // const Mdouble distanceSquared0 = Vec3D::getLengthSquared(vertices_[index0].position - position);
141 // const Mdouble distanceSquared1 = Vec3D::getLengthSquared(vertices_[index1].position - position);
142 // const Mdouble distanceSquared2 = Vec3D::getLengthSquared(vertices_[index2].position - position);
143 //
144 // // Assign to closest vertex
145 // const unsigned indexFinal = distanceSquared0 < distanceSquared1 ?
146 // (distanceSquared0 < distanceSquared2 ? index0 : index2) :
147 // (distanceSquared1 < distanceSquared2 ? index1 : index2);
148 // debrisContainer[indexFinal] += debris;
149 
150 
151  // Barycentric interpolation
152  // Get vertices positions and shift by debris position.
153  const Vec3D v0 = vertices_[index0].position - position;
154  const Vec3D v1 = vertices_[index1].position - position;
155  const Vec3D v2 = vertices_[index2].position - position;
156 
157  // Area of each section, and total area.
158  // When the debris position does not lie perfectly in plane, the sum of the area of the sections differs from the
159  // area of the triangle itself. It is therefore better to use the sum, even though still a small error in the
160  // fractions remain, depending on how far off plane the position is.
161  Mdouble A0 = 0.5 * Vec3D::cross(v0, v1).getLength();
162  Mdouble A1 = 0.5 * Vec3D::cross(v1, v2).getLength();
163  Mdouble A2 = 0.5 * Vec3D::cross(v2, v0).getLength();
164  Mdouble A = A0 + A1 + A2;
165 
166  // Fraction of debris is equal to the area opposite to vertex divided by total area.
167  debrisContainer[index0] += debris * A1 / A;
168  debrisContainer[index1] += debris * A2 / A;
169  debrisContainer[index2] += debris * A0 / A;
170 }
@ A
Definition: StatisticsVector.h:42
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:163
static Mdouble getLength(const Vec3D &a)
Calculates the length of a Vec3D: .
Definition: Vector.cc:331

References A, Vec3D::cross(), Vec3D::getLength(), TriangleMeshWall::Triangle::vertexIndices, and TriangleMeshWall::vertices_.

Referenced by computeWear().

◆ write()

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

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

Reimplemented from TriangleMeshWall.

74 {
76 
77  os << " wearCoefficient " << wearCoefficient_ << " hardness " << hardness_ << " wearAcceleration " << wearAcceleration_;
78 }
void write(std::ostream &os) const override
Writes a TriangleMeshWall to an output stream, for example a restart file.
Definition: TriangleMeshWall.cc:113

References hardness_, wearAcceleration_, wearCoefficient_, and TriangleMeshWall::write().

Member Data Documentation

◆ hardness_

Mdouble WearableTriangleMeshWall::hardness_
private

◆ wearAcceleration_

Mdouble WearableTriangleMeshWall::wearAcceleration_
private

Accelerates the wear process, to reduce simulation time needed to get results.

Referenced by computeWear(), read(), setWearAcceleration(), WearableTriangleMeshWall(), and write().

◆ wearCoefficient_

Mdouble WearableTriangleMeshWall::wearCoefficient_
private

Dimensionless wear coefficient.

Referenced by computeWear(), read(), setWearCoefficient(), WearableTriangleMeshWall(), and write().


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