TriangulatedWall Class Reference

A TriangulatedWall is a triangulation created from a set of vertices and a n-by-3 connectivity matrix defining n faces. More...

#include <TriangulatedWall.h>

+ Inheritance diagram for TriangulatedWall:

Classes

struct  Face
 

Public Member Functions

 TriangulatedWall ()
 Default constructor. More...
 
 TriangulatedWall (const TriangulatedWall &other)
 Copy constructor. More...
 
 TriangulatedWall (std::string filename, const ParticleSpecies *species)
 Constructor setting values. More...
 
 TriangulatedWall (const std::vector< Vec3D > &points, const std::vector< std::vector< unsigned >> &cells, const ParticleSpecies *species)
 Constructor setting values. More...
 
void set (const std::vector< Vec3D > &points, const std::vector< std::vector< unsigned >> &cells)
 
void readVTK (std::string filename)
 
void writeVTK (VTKContainer &vtk) const override
 
 ~TriangulatedWall () override
 Destructor. More...
 
TriangulatedWalloperator= (const TriangulatedWall &other)
 
TriangulatedWallcopy () const override
 Wall copy method. It calls the copy constructor of this Wall, useful for polymorphism. More...
 
bool getDistanceAndNormal (const BaseParticle &p, Mdouble &distance, Vec3D &normal_return) const override
 Compute the distance from the wall for a given BaseParticle and return if there is a collision. If there is a collision, also return the normal vector. More...
 
void move (const Vec3D &move) override
 Move the TriangulatedWall to a new position, which is a Vec3D from the old position. More...
 
void read (std::istream &is) override
 Reads an TriangulatedWall from an input stream, for example a restart file. More...
 
void write (std::ostream &os) const override
 Writes an TriangulatedWall to an output stream, for example a restart file. More...
 
std::string getName () const override
 Returns the name of the object, here the string "TriangulatedWall". More...
 
BaseInteractiongetInteractionWith (BaseParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler) override
 Get the interaction between this TriangulatedWall and given BaseParticle at a given time. 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...
 
virtual bool isLocal (Vec3D &min, Vec3D &max) const
 
bool getLinePlaneIntersect (Vec3D &intersect, const Vec3D &p0, const Vec3D &p1, const Vec3D &n, const Vec3D &p)
 
bool isInsideWallVTK (const Vec3D &point, const Vec3D &normal, const Vec3D &position) const
 
void projectOntoWallVTK (Vec3D &point0, const Vec3D &point1, const Vec3D &normal, const Vec3D &position) const
 
void intersectVTK (std::vector< Vec3D > &points, Vec3D normal, Vec3D position) const
 
virtual BaseInteractiongetInteractionWithSuperQuad (SuperQuadricParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler)
 
void getVTK (std::vector< Vec3D > &points, std::vector< std::vector< double >> &triangleStrips)
 
const Vec3D getAxis () const
 
virtual void actionsOnRestart ()
 No implementation but can be overidden in its derived classes. More...
 
virtual void actionsAfterParticleGhostUpdate ()
 No implementation but can be overidden in its derived classes. More...
 
virtual void handleParticleAddition (unsigned int id, BaseParticle *p)
 Handles the addition of particles to the particleHandler. More...
 
virtual void handleParticleRemoval (unsigned int id)
 Handles the addition of particles to the particleHandler. More...
 
virtual void checkInteractions (InteractionHandler *interactionHandler, unsigned int timeStamp)
 Check if all interactions are valid. More...
 
bool getVTKVisibility () const
 
void setVTKVisibility (bool vtkVisibility)
 
void addRenderedWall (BaseWall *w)
 
BaseWallgetRenderedWall (size_t i) const
 
std::vector< BaseWall * > getRenderedWalls () const
 
void removeRenderedWalls ()
 
void renderWall (VTKContainer &vtk)
 
void addParticlesAtWall (unsigned numElements=50)
 
void setVelocityControl (Vec3D forceGoal, Vec3D gainFactor, Vec3D baseVelocity)
 
virtual void writeWallDetailsVTK (VTKData &data) const
 
virtual void computeWear ()
 
- Public Member Functions inherited from BaseInteractable
 BaseInteractable ()
 Default BaseInteractable constructor. More...
 
 BaseInteractable (const BaseInteractable &p)
 Copy constructor. More...
 
 ~BaseInteractable () override
 Destructor, it simply destructs the BaseInteractable and all the objects it contains. More...
 
unsigned int getIndSpecies () const
 Returns the index of the species associated with the interactable object. More...
 
const ParticleSpeciesgetSpecies () const
 Returns a pointer to the species of this BaseInteractable. More...
 
void setSpecies (const ParticleSpecies *species)
 Sets the species of this BaseInteractable. More...
 
const Vec3DgetForce () const
 Returns the force on this BaseInteractable. More...
 
const Vec3DgetTorque () const
 Returns the torque on this BaseInteractable. More...
 
void setForce (const Vec3D &force)
 Sets the force on this BaseInteractable. More...
 
void setTorque (const Vec3D &torque)
 Sets the torque on this BaseInteractable. More...
 
void addForce (const Vec3D &addForce)
 Adds an amount to the force on this BaseInteractable. More...
 
void addTorque (const Vec3D &addTorque)
 Adds an amount to the torque on this BaseInteractable. More...
 
virtual void resetForceTorque (int numberOfOMPthreads)
 
void sumForceTorqueOMP ()
 
const Vec3DgetPosition () const
 Returns the position of this BaseInteractable. More...
 
const QuaterniongetOrientation () const
 Returns the orientation of this BaseInteractable. More...
 
virtual void setPosition (const Vec3D &position)
 Sets the position of this BaseInteractable. More...
 
void setOrientationViaNormal (Vec3D normal)
 Sets the orientation of this BaseInteractable by defining the vector that results from the rotation of the (1,0,0) vector. More...
 
void setOrientationViaEuler (Vec3D eulerAngle)
 Sets the orientation of this BaseInteractable by defining the euler angles. More...
 
virtual void setOrientation (const Quaternion &orientation)
 Sets the orientation of this BaseInteractable. More...
 
virtual void rotate (const Vec3D &angularVelocityDt)
 Rotates this BaseInteractable. More...
 
const std::vector< BaseInteraction * > & getInteractions () const
 Returns a list of interactions which belong to this interactable. More...
 
void addInteraction (BaseInteraction *I)
 Adds an interaction to this BaseInteractable. More...
 
bool removeInteraction (BaseInteraction *I)
 Removes an interaction from this BaseInteractable. More...
 
void copyInteractionsForPeriodicParticles (const BaseInteractable &p)
 Copies interactions to this BaseInteractable whenever a periodic copy made. More...
 
void setVelocity (const Vec3D &velocity)
 set the velocity of the BaseInteractable. More...
 
void setAngularVelocity (const Vec3D &angularVelocity)
 set the angular velocity of the BaseInteractble. More...
 
void addVelocity (const Vec3D &velocity)
 adds an increment to the velocity. More...
 
void addAngularVelocity (const Vec3D &angularVelocity)
 add an increment to the angular velocity. More...
 
virtual const Vec3DgetVelocity () const
 Returns the velocity of this interactable. More...
 
virtual const Vec3DgetAngularVelocity () const
 Returns the angular velocity of this interactable. More...
 
void setPrescribedPosition (const std::function< Vec3D(double)> &prescribedPosition)
 Allows the position of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedPosition (double time)
 Computes the position from the user defined prescribed position function. More...
 
void setPrescribedVelocity (const std::function< Vec3D(double)> &prescribedVelocity)
 Allows the velocity of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedVelocity (double time)
 Computes the velocity from the user defined prescribed velocity function. More...
 
void setPrescribedOrientation (const std::function< Quaternion(double)> &prescribedOrientation)
 Allows the orientation of the infinite mass interactbale to be prescribed. More...
 
void applyPrescribedOrientation (double time)
 Computes the orientation from the user defined prescribed orientation function. More...
 
void setPrescribedAngularVelocity (const std::function< Vec3D(double)> &prescribedAngularVelocity)
 Allows the angular velocity of the infinite mass interactable to be prescribed. More...
 
void applyPrescribedAngularVelocity (double time)
 Computes the angular velocity from the user defined prescribed angular velocity. More...
 
virtual const Vec3D getVelocityAtContact (const Vec3D &contact) const
 Returns the velocity at the contact point, use by many force laws. More...
 
void integrateBeforeForceComputation (double time, double timeStep)
 This is part of integrate routine for objects with infinite mass. More...
 
void integrateAfterForceComputation (double time, double timeStep)
 This is part of the integration routine for objects with infinite mass. More...
 
virtual Mdouble getInvMass () const
 
virtual Mdouble getCurvature (const Vec3D &labFixedCoordinates) const
 
virtual bool isFaceContact (const Vec3D &normal) const
 
- Public Member Functions inherited from BaseObject
 BaseObject ()=default
 Default constructor. More...
 
 BaseObject (const BaseObject &p)=default
 Copy constructor, copies all the objects BaseObject contains. More...
 
virtual ~BaseObject ()=default
 virtual destructor More...
 
virtual void moveInHandler (unsigned int index)
 Except that it is virtual, it does the same thing as setIndex() does. More...
 
void setIndex (unsigned int index)
 Allows one to assign an index to an object in the handler/container. More...
 
void setId (unsigned long id)
 Assigns a unique identifier to each object in the handler (container) which remains constant even after the object is deleted from the container/handler. More...
 
unsigned int getIndex () const
 Returns the index of the object in the handler. More...
 
unsigned int getId () const
 Returns the unique identifier of any particular object. More...
 
void setGroupId (unsigned groupId)
 
unsigned getGroupId () const
 

Protected Attributes

std::vector< Vec3Dvertex_
 
std::vector< Faceface_
 

Private Member Functions

void setNormalsAndNeighbours ()
 

Additional Inherited Members

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

Detailed Description

A TriangulatedWall is a triangulation created from a set of vertices and a n-by-3 connectivity matrix defining n faces.

It is initialised by a unstructured grid vtk file:

TriangulatedWall w("TriangulatedWallSimple.vtk",speciesHandler.getLastObject());
A TriangulatedWall is a triangulation created from a set of vertices and a n-by-3 connectivity matrix...
Definition: TriangulatedWall.h:52

The file consists of a set of vertices and a n-by-3 connectivity matrix defining n faces. Three vertices form a face; the face normal is oriented such that the vertices are ordered in anticlockwise direction around the normal.

Particles interact with a TriangulatedWall when they contact a face (from either side), a edge, or a vertex.

For a demonstration on how to use this class, see TriangulatedWallSelfTest (shown in the image below).

Constructor & Destructor Documentation

◆ TriangulatedWall() [1/4]

TriangulatedWall::TriangulatedWall ( )

Default constructor.

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

References DEBUG, and logger.

Referenced by copy().

◆ TriangulatedWall() [2/4]

TriangulatedWall::TriangulatedWall ( const TriangulatedWall other)

Copy constructor.

Parameters
[in]otherThe TriangulatedWall that must be copied.
42  : BaseWall(other)
43 {
44  face_ = other.face_;
45  vertex_ = other.vertex_;
46  //now reset the pointers
47  for (unsigned f = 0; f < face_.size(); f++)
48  {
49  for (unsigned i = 0; i < 3; i++)
50  {
51  if (other.face_[f].neighbor[i])
52  {
53  face_[f].neighbor[i] = &face_[other.face_[f].neighbor[i] - &other.face_[0]];
54  } //else nullptr
55  face_[f].vertex[i] = &vertex_[other.face_[f].vertex[i] - &other.vertex_[0]];
56  }
57  }
58  logger(DEBUG, "TriangulatedWall(TriangulatedWall&) constructed.");
59 }
BaseWall()
Default constructor.
Definition: BaseWall.cc:36
std::vector< Face > face_
Definition: TriangulatedWall.h:158
std::vector< Vec3D > vertex_
Definition: TriangulatedWall.h:153
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51

References DEBUG, face_, constants::i, logger, and vertex_.

◆ TriangulatedWall() [3/4]

TriangulatedWall::TriangulatedWall ( std::string  filename,
const ParticleSpecies species 
)

Constructor setting values.

62 {
63  setSpecies(species);
64  readVTK(filename);
65 }
void setSpecies(const ParticleSpecies *species)
Defines the species of the current wall.
Definition: BaseWall.cc:169
void readVTK(std::string filename)
Definition: TriangulatedWall.cc:70

References readVTK(), and BaseWall::setSpecies().

◆ TriangulatedWall() [4/4]

TriangulatedWall::TriangulatedWall ( const std::vector< Vec3D > &  points,
const std::vector< std::vector< unsigned >> &  cells,
const ParticleSpecies species 
)

Constructor setting values.

◆ ~TriangulatedWall()

TriangulatedWall::~TriangulatedWall ( )
override

Destructor.

227 {
228  logger(DEBUG, "~TriangulatedWall() has been called.");
229 }

References DEBUG, and logger.

Member Function Documentation

◆ copy()

TriangulatedWall * TriangulatedWall::copy ( ) const
overridevirtual

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

Returns
pointer to a TriangulatedWall object allocated using new.

Implements BaseWall.

Reimplemented in WearableTriangulatedWall.

248 {
249  return new TriangulatedWall(*this);
250 }
TriangulatedWall()
Default constructor.
Definition: TriangulatedWall.cc:33

References TriangulatedWall().

Referenced by operator=().

◆ getDistanceAndNormal()

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

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

Parameters
[in]pBaseParticle 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 TriangulatedWall. If there is a collision, this function also computes the distance between the BaseParticle and TriangulatedWall and the normal of the TriangulatedWall at the intersection point. It does this by calling TriangulatedWall::getDistanceAndNormal(const Vec3D& , Mdouble , Mdouble&, Vec3D&) const. Since this function should be called before calculating any Particle-Wall interactions, it can also be used to set the normal vector in case of curved walls.

NOTE: THIS ONLY RETURNS ONE OF POSSIBLY MANY INTERACTIONS; it's only used for finding out if interactions exist, so should be fine

Implements BaseWall.

269 {
270  Mdouble interactionRadius = p.getWallInteractionRadius(this);
271  //it's important to use a reference here, as the faces use pointers
272  for (const auto& face : face_)
273  {
274  if (face.getDistanceAndNormal(p, distance, normal_return, interactionRadius))
275  return true;
276  }
277  return false;
278 }
double Mdouble
Definition: GeneralDefine.h:34
Mdouble getWallInteractionRadius(const BaseWall *wall) const
returns the radius plus the interactionDistance
Definition: BaseParticle.h:386

References face_, and BaseParticle::getWallInteractionRadius().

Referenced by TriangulatedStepWallSelfTest::setupInitialConditions().

◆ getInteractionWith()

BaseInteraction * TriangulatedWall::getInteractionWith ( BaseParticle p,
unsigned  timeStamp,
InteractionHandler interactionHandler 
)
overridevirtual

Get the interaction between this TriangulatedWall and given BaseParticle at a given time.

Gets the interaction with this wall. In case of multiple contact points return the one with the shortest distance (i.e. largest overlap)

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 wall and the BaseParticle at the timeStamp.

Reimplemented from BaseWall.

357 {
358  Mdouble distance, minDistance = constants::inf;
359  Vec3D normal, minNormal;
360  Mdouble interactionRadius = p->getWallInteractionRadius(this);
361  for (const auto& face : face_)
362  {
363  if (face.getDistanceAndNormal(*p, distance, normal,interactionRadius) && distance < minDistance)
364  {
365  minDistance = distance;
366  minNormal = normal;
367  }
368  }
369 
370  if (minDistance != constants::inf)
371  {
372  distance = minDistance;
373  normal = -minNormal;
374  BaseInteraction* const c = interactionHandler->getInteraction(p, this, timeStamp, normal);
375  c->setNormal(normal);
376  c->setDistance(distance);
377  c->setOverlap(p->getRadius() - distance);
378  c->setContactPoint(p->getPosition() - distance * normal);
379  return c;
380  }
381 
382  return nullptr;
383 }
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
void setDistance(Mdouble distance)
Sets the interaction distance between the two interacting objects.
Definition: BaseInteraction.cc:221
void setContactPoint(Vec3D contactPoint)
Set the location of the contact point between the two interacting objects.
Definition: BaseInteraction.cc:240
void setNormal(Vec3D normal)
Sets the normal vector between the two interacting objects.
Definition: BaseInteraction.cc:212
void setOverlap(Mdouble overlap)
Set the overlap between the two interacting object.
Definition: BaseInteraction.cc:231
Mdouble getRadius() const
Returns the particle's radius.
Definition: BaseParticle.h:348
BaseInteraction * getInteraction(BaseInteractable *P, BaseInteractable *I, unsigned timeStamp)
Returns the Interaction between the BaseInteractable's P and I.
Definition: InteractionHandler.cc:146
Definition: Vector.h:51
const Mdouble inf
Definition: GeneralDefine.h:44

References face_, InteractionHandler::getInteraction(), BaseInteractable::getPosition(), BaseParticle::getRadius(), BaseParticle::getWallInteractionRadius(), constants::inf, BaseInteraction::setContactPoint(), BaseInteraction::setDistance(), BaseInteraction::setNormal(), and BaseInteraction::setOverlap().

◆ getName()

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

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

Returns
The string "TriangulatedWall".

Implements BaseObject.

Reimplemented in WearableTriangulatedWall.

343 {
344  return "TriangulatedWall";
345 }

◆ move()

void TriangulatedWall::move ( const Vec3D move)
overridevirtual

Move the TriangulatedWall to a new position, which is a Vec3D from the old position.

Parameters
[in]moveA reference to a Vec3D that denotes the direction and length it should be moved with.

A function that moves the TriangulatedWall in a certain direction by both moving the walls and all intersections. Note that the directions of the intersections are not moved since they don't change when moving the TriangulatedWall as a whole.

Todo:
We should use the position_ and orientation_ of the TriangulatedWall; that way, TriangulatedWall can be moved with the standard BaseInteractable::move function, getting rid of an anomaly in the code and removing the virtual from the move function.
Author
weinhartt

Reimplemented from BaseInteractable.

292 {
294  for (auto& v : vertex_)
295  {
296  v += move;
297  }
298 }
virtual void move(const Vec3D &move)
Moves this BaseInteractable by adding an amount to the position.
Definition: BaseInteractable.cc:215
void move(const Vec3D &move) override
Move the TriangulatedWall to a new position, which is a Vec3D from the old position.
Definition: TriangulatedWall.cc:291

References BaseInteractable::move(), and vertex_.

◆ operator=()

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

Copy assignment operator.

Parameters
[in]otherThe TriangulatedWall that must be copied.
235 {
236  logger(DEBUG, "TriangulatedWall::operator= called.");
237  if (this == &other)
238  {
239  return *this;
240  }
241  return *(other.copy());
242 }
TriangulatedWall * copy() const override
Wall copy method. It calls the copy constructor of this Wall, useful for polymorphism.
Definition: TriangulatedWall.cc:247

References copy(), DEBUG, and logger.

◆ read()

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

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

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

Reimplemented from BaseWall.

Reimplemented in WearableTriangulatedWall.

304 {
306 }

Referenced by WearableTriangulatedWall::read().

◆ readVTK()

void TriangulatedWall::readVTK ( std::string  filename)
Parameters
[in]filenamename of vtk input file, e.g. TriangulatedWallSelfTest.vtk
71 {
72  std::fstream file;
73  file.open(filename.c_str(), std::ios::in);
74  logger.assert_always(file.is_open(), "File opening failed: %", filename);
75 
76  std::string dummy;
77  getline(file, dummy);
78  getline(file, dummy);
79  getline(file, dummy);
80  getline(file, dummy);
81  //read points
82  unsigned num;
83  file >> dummy >> num >> dummy;
84  vertex_.reserve(num);
85  Vec3D v;
86  for (unsigned i = 0; i < num; i++)
87  {
88  file >> v.X >> v.Y >> v.Z;
89  vertex_.push_back(v);
90  }
91  //read faces
92  file >> dummy >> num >> dummy;
93  face_.reserve(num);
94  Face f;
95  unsigned id0, id1, id2;
96  for (unsigned i = 0; i < num; i++)
97  {
98  file >> dummy >> id0 >> id1 >> id2;
99  f.vertex[0] = &vertex_[id0];
100  f.vertex[1] = &vertex_[id1];
101  f.vertex[2] = &vertex_[id2];
102  face_.push_back(f);
103  }
104 
105  file.close();
106 
108 }
void setNormalsAndNeighbours()
Definition: TriangulatedWall.cc:138
Mdouble Y
Definition: Vector.h:66
Mdouble Z
Definition: Vector.h:66
Mdouble X
the vector components
Definition: Vector.h:66

References face_, constants::i, logger, setNormalsAndNeighbours(), TriangulatedWall::Face::vertex, vertex_, Vec3D::X, Vec3D::Y, and Vec3D::Z.

Referenced by TriangulatedWall().

◆ set()

void TriangulatedWall::set ( const std::vector< Vec3D > &  points,
const std::vector< std::vector< unsigned >> &  cells 
)
118 {
119  //points
120  vertex_ = points;
121 
122  //cells
123  unsigned num = cells.size();
124  face_.reserve(num);
125  Face f;
126  for (unsigned i = 0; i < num; i++)
127  {
128  f.vertex[0] = &vertex_[cells[i][0]];
129  f.vertex[1] = &vertex_[cells[i][1]];
130  f.vertex[2] = &vertex_[cells[i][2]];
131  face_.push_back(f);
132  }
133 
134  //normals, positions, neighbours
136 }

References face_, constants::i, setNormalsAndNeighbours(), TriangulatedWall::Face::vertex, and vertex_.

Referenced by WearableTriangulatedWall::WearableTriangulatedWall().

◆ setNormalsAndNeighbours()

void TriangulatedWall::setNormalsAndNeighbours ( )
private
139 {
140  //set normals and positions
141  for (auto& face: face_)
142  {
143  face.normal = Vec3D::getUnitVector(
144  Vec3D::cross(*face.vertex[1] - *face.vertex[0], *face.vertex[2] - *face.vertex[0]));
145  }
146  //set neighbours
147  for (auto face0 = face_.begin(); face0 + 1 != face_.end(); face0++)
148  {
149  for (auto face1 = face0 + 1; face1 != face_.end(); face1++)
150  {
151  if (face0->vertex[0] == face1->vertex[0])
152  {
153  if (face0->vertex[1] == face1->vertex[2])
154  { //edge 0=2
155  face0->neighbor[0] = &*face1;
156  face1->neighbor[2] = &*face0;
157  }
158  else if (face0->vertex[2] == face1->vertex[1])
159  { //edge 2=0
160  face0->neighbor[2] = &*face1;
161  face1->neighbor[0] = &*face0;
162  }
163  }
164  else if (face0->vertex[0] == face1->vertex[1])
165  {
166  if (face0->vertex[1] == face1->vertex[0])
167  { //edge 0=0
168  face0->neighbor[0] = &*face1;
169  face1->neighbor[0] = &*face0;
170  }
171  else if (face0->vertex[2] == face1->vertex[2])
172  { //edge 2=1
173  face0->neighbor[2] = &*face1;
174  face1->neighbor[1] = &*face0;
175  }
176  }
177  else if (face0->vertex[0] == face1->vertex[2])
178  {
179  if (face0->vertex[1] == face1->vertex[1])
180  { //edge 0=1
181  face0->neighbor[0] = &*face1;
182  face1->neighbor[1] = &*face0;
183  }
184  else if (face0->vertex[2] == face1->vertex[0])
185  { //edge 2=2
186  face0->neighbor[2] = &*face1;
187  face1->neighbor[2] = &*face0;
188  }
189  }
190  else if (face0->vertex[1] == face1->vertex[0])
191  {
192  if (face0->vertex[2] == face1->vertex[2])
193  { //edge 1=2
194  face0->neighbor[1] = &*face1;
195  face1->neighbor[2] = &*face0;
196  }
197  }
198  else if (face0->vertex[1] == face1->vertex[1])
199  {
200  if (face0->vertex[2] == face1->vertex[0])
201  { //edge 1=0
202  face0->neighbor[1] = &*face1;
203  face1->neighbor[0] = &*face0;
204  }
205  }
206  else if (face0->vertex[1] == face1->vertex[2])
207  {
208  if (face0->vertex[2] == face1->vertex[1])
209  { //edge 1=1
210  face0->neighbor[1] = &*face1;
211  face1->neighbor[1] = &*face0;
212  }
213  }
214 
215  }
216  }
217  //set edge normals (inwards facing)
218  for (auto& face: face_)
219  {
220  face.edgeNormal[0] = Vec3D::getUnitVector(Vec3D::cross(face.normal, *face.vertex[1] - *face.vertex[0]));
221  face.edgeNormal[1] = Vec3D::getUnitVector(Vec3D::cross(face.normal, *face.vertex[2] - *face.vertex[1]));
222  face.edgeNormal[2] = Vec3D::getUnitVector(Vec3D::cross(face.normal, *face.vertex[0] - *face.vertex[2]));
223  }
224 }
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:163
static Vec3D getUnitVector(const Vec3D &a)
Returns a unit Vec3D based on a.
Definition: Vector.cc:345

References Vec3D::cross(), face_, and Vec3D::getUnitVector().

Referenced by readVTK(), and set().

◆ write()

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

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

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

Reimplemented from BaseWall.

Reimplemented in WearableTriangulatedWall.

313 {
314  os << "Vertices " << vertex_.size();
315  for (const auto& vertex: vertex_)
316  {
317  os << " " << vertex;
318  }
319  os << std::endl;
320  //os << "Faces " << face_.size() << std::endl;
321  unsigned counter = 0;
322  for (const auto& face: face_)
323  {
324  os << "Face " << counter++
325  << " vertex " << face.vertex[0] - &vertex_[0]
326  << " " << face.vertex[1] - &vertex_[0]
327  << " " << face.vertex[2] - &vertex_[0]
328  << " neighbor " << (face.neighbor[0] ? (face.neighbor[0] - &face_[0]) : -1)
329  << " " << (face.neighbor[1] ? (face.neighbor[1] - &face_[0]) : -1)
330  << " " << (face.neighbor[2] ? (face.neighbor[2] - &face_[0]) : -1)
331  << " normal " << face.normal
332  << " edgeNormal " << face.edgeNormal[0]
333  << " " << face.edgeNormal[1]
334  << " " << face.edgeNormal[2]
335  << std::endl;
336  }
337 }

References face_, and vertex_.

Referenced by TriangulatedStepSelfTest::setupInitialConditions(), TriangulatedStepWallSelfTest::setupInitialConditions(), TriangulatedWallSelfTest::setupInitialConditions(), and WearableTriangulatedWall::write().

◆ writeVTK()

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

493 {
494  const int s = vtk.points.size();
495  for (auto v : vertex_)
496  {
497  vtk.points.push_back(v);
498  }
499  for (auto f : face_)
500  {
501  std::vector<double> cell;
502  cell.reserve(3);
503  cell.push_back(s + f.vertex[0] - &vertex_[0]);
504  cell.push_back(s + f.vertex[1] - &vertex_[0]);
505  cell.push_back(s + f.vertex[2] - &vertex_[0]);
506  vtk.triangleStrips.push_back(cell);
507  }
508 }
std::vector< std::vector< double > > triangleStrips
Definition: BaseWall.h:40
std::vector< Vec3D > points
Definition: BaseWall.h:39

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

Member Data Documentation

◆ face_

◆ vertex_

std::vector<Vec3D> TriangulatedWall::vertex_
protected

stores the vertex coordinates

Referenced by move(), readVTK(), set(), TriangulatedWall(), write(), and writeVTK().


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