MercuryDPM  Alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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
 Struct used to store the properties of a face needed for contact detection. More...
 

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...
 
void readVTK (std::string filename)
 
void writeVTK (VTKContainer &vtk) const override
 
virtual ~TriangulatedWall ()
 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...
 
std::vector< BaseInteraction * > getInteractionWith (BaseParticle *p, Mdouble 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. It makes an empty BaseWall. More...
 
 BaseWall (const BaseWall &w)
 Copy constructor. More...
 
virtual ~BaseWall ()
 Default destructor. More...
 
virtual MERCURY_DEPRECATED void clear ()
 A function that removes all data from this BaseWall, so sets handler_ to nullptr. More...
 
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)
 Define the species of this wall using the index of the species in the SpeciesHandler in this DPMBase. More...
 
void setSpecies (const ParticleSpecies *species)
 Define the species of this 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, const Vec3D normal, const Vec3D position) const
 
- Public Member Functions inherited from BaseInteractable
 BaseInteractable ()
 Default BaseInteractable constructor, it simply creates an empty BaseInteractable. More...
 
 BaseInteractable (const BaseInteractable &p)
 Copy constructor. It copies the BaseInteractable and all objects it contains. More...
 
virtual ~BaseInteractable ()
 Destructor, it simply destructs the BaseInteractable and all the objects it contains. More...
 
unsigned int getIndSpecies () const
 Returns the index of the Species of this BaseInteractable. 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...
 
const Vec3DgetPosition () const
 Returns the position of this BaseInteractable. More...
 
const Vec3DgetOrientation () const
 Returns the orientation of this BaseInteractable. More...
 
void setPosition (const Vec3D &position)
 Sets the position of this BaseInteractable. More...
 
void setOrientation (const Vec3D &orientation)
 Sets the orientation of this BaseInteractable. More...
 
void rotate (const Vec3D &rotate)
 Rotates this BaseInteractable. More...
 
const std::list
< BaseInteraction * > & 
getInteractions () const
 Returns a reference to the list of interactions in this BaseInteractable. 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< Vec3D(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...
 
- Public Member Functions inherited from BaseObject
 BaseObject ()
 Default constructor. More...
 
 BaseObject (const BaseObject &p)
 Copy constructor, copies all the objects BaseObject contains. More...
 
virtual ~BaseObject ()
 virtual destructor More...
 
virtual void moveInHandler (const unsigned int index)
 Except that it is virtual, it does the same thing as setIndex() does. More...
 
void setIndex (const unsigned int index)
 Allows one to assign an index to an object in the handler/container. More...
 
void setId (const unsigned int 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...
 

Private Attributes

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

Additional Inherited Members

- Static Public Member Functions inherited from BaseWall
static void addToVTK (const std::vector< Vec3D > &points, VTKContainer &vtk)
 
- Public Attributes inherited from BaseWall
WallHandlerhandler_
 

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());

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).

Definition at line 51 of file TriangulatedWall.h.

Constructor & Destructor Documentation

TriangulatedWall::TriangulatedWall ( )

Default constructor.

Definition at line 33 of file TriangulatedWall.cc.

References DEBUG, and logger.

Referenced by copy().

34 {
35  logger(DEBUG, "TriangulatedWall() constructed.");
36 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
TriangulatedWall::TriangulatedWall ( const TriangulatedWall other)

Copy constructor.

Parameters
[in]otherThe TriangulatedWall that must be copied.

Definition at line 41 of file TriangulatedWall.cc.

References DEBUG, face_, logger, and vertex_.

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  face_[f].neighbor[i] = &face_[other.face_[f].neighbor[i] - &other.face_[0]];
53  } //else nullptr
54  face_[f].vertex[i] = &vertex_[other.face_[f].vertex[i] - &other.vertex_[0]];
55  }
56  }
57  logger(DEBUG, "TriangulatedWall(TriangulatedWall&) constructed.");
58 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
std::vector< Vec3D > vertex_
BaseWall()
Default constructor. It makes an empty BaseWall.
Definition: BaseWall.cc:31
std::vector< Face > face_
TriangulatedWall::TriangulatedWall ( std::string  filename,
const ParticleSpecies species 
)

Constructor setting values.

Definition at line 60 of file TriangulatedWall.cc.

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

61 {
62  setSpecies(species);
63  readVTK(filename);
64 }
void readVTK(std::string filename)
void setSpecies(const ParticleSpecies *species)
Define the species of this wall.
Definition: BaseWall.cc:113
TriangulatedWall::~TriangulatedWall ( )
virtual

Destructor.

Definition at line 167 of file TriangulatedWall.cc.

References DEBUG, and logger.

168 {
169  logger(DEBUG, "~TriangulatedWall() has been called.");
170 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")

Member Function Documentation

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.

Definition at line 188 of file TriangulatedWall.cc.

References TriangulatedWall().

Referenced by operator=().

189 {
190  return new TriangulatedWall(*this);
191 }
TriangulatedWall()
Default constructor.
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.

Definition at line 209 of file TriangulatedWall.cc.

References face_.

210 {
211  //it's important to use a reference here, as the faces use pointers
212  for (const auto& face : face_)
213  {
214  if (face.getDistanceAndNormal(p, distance, normal_return))
215  return true;
216  }
217  return false;
218 }
std::vector< Face > face_
std::vector< BaseInteraction * > TriangulatedWall::getInteractionWith ( BaseParticle p,
Mdouble  timeStamp,
InteractionHandler interactionHandler 
)
overridevirtual

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

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

Reimplemented from BaseWall.

Definition at line 292 of file TriangulatedWall.cc.

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

294 {
295  Mdouble distance;
296  Vec3D normal;
297  Face* neighbor;
298  //write(std::cout);
299  std::vector<BaseInteraction*> interactions;
300  for (const auto& face : face_)
301  {
302  if (face.getDistanceAndNormal(*p, distance, normal))
303  {
304  normal = -normal;
305  //logger(INFO,"t=% n=% f=%(%) c=%", timeStamp,normal,&face,&face-&face_[0],p->getPosition() - distance * normal);
306  BaseInteraction* const c = interactionHandler->getInteraction(p, this, timeStamp, normal);
307  //logger.assert(c,"BaseInteraction not set");
308  c->setNormal(normal);
309  c->setDistance(distance);
310  c->setOverlap(p->getRadius() - distance);
311  c->setContactPoint(p->getPosition() - distance * normal);
312  interactions.push_back(c);
313  }
314  }
315  return interactions;
316 }
void setNormal(Vec3D normal)
Sets the normal vector between the two interacting objects.
void setOverlap(Mdouble overlap)
Set the overlap between the two interacting object.
BaseInteraction * getInteraction(BaseInteractable *P, BaseInteractable *I, Mdouble timeStamp)
Returns the Interaction between the BaseInteractable's P and I.
double Mdouble
void setContactPoint(Vec3D contactPoint)
Set the location of the contact point between the two interacting objects.
void setDistance(Mdouble distance)
Sets the interaction distance between the two interacting objects.
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Stores information about interactions between two interactable objects; often particles but could be ...
Mdouble getRadius() const
Returns the particle's radius_.
std::vector< Face > face_
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:45
std::string TriangulatedWall::getName ( ) const
overridevirtual

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

Returns
The string "TriangulatedWall".

Implements BaseObject.

Definition at line 280 of file TriangulatedWall.cc.

281 {
282  return "TriangulatedWall";
283 }
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.

Definition at line 231 of file TriangulatedWall.cc.

References BaseInteractable::move(), and vertex_.

232 {
234  for (auto& v : vertex_)
235  {
236  v += move;
237  }
238 }
std::vector< Vec3D > vertex_
void move(const Vec3D &move) override
Move the TriangulatedWall to a new position, which is a Vec3D from the old position.
virtual void move(const Vec3D &move)
Moves this BaseInteractable by adding an amount to the position.
TriangulatedWall & TriangulatedWall::operator= ( const TriangulatedWall other)

Copy assignment operator.

Parameters
[in]otherThe TriangulatedWall that must be copied.

Definition at line 175 of file TriangulatedWall.cc.

References copy(), DEBUG, and logger.

176 {
177  logger(DEBUG, "TriangulatedWall::operator= called.");
178  if (this == &other)
179  {
180  return *this;
181  }
182  return *(other.copy());
183 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
TriangulatedWall * copy() const override
Wall copy method. It calls the copy constructor of this Wall, useful for polymorphism.
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.

Definition at line 243 of file TriangulatedWall.cc.

244 {
246 }
void TriangulatedWall::readVTK ( std::string  filename)
Parameters
[in]filenamename of vtk input file, e.g. TriangulatedWallSelfTest.vtk

Definition at line 69 of file TriangulatedWall.cc.

References Vec3D::cross(), face_, Vec3D::getUnitVector(), logger, TriangulatedWall::Face::vertex, vertex_, Vec3D::X, Vec3D::Y, and Vec3D::Z.

Referenced by TriangulatedWall().

70 {
71  std::fstream file;
72  file.open(filename.c_str(), std::ios::in);
73  logger.assert_always(file.is_open(), "File opening failed: %",filename);
74 
75  std::string dummy;
76  getline(file, dummy);
77  getline(file, dummy);
78  getline(file, dummy);
79  getline(file, dummy);
80  //read points
81  unsigned num;
82  file >> dummy >> num >> dummy;
83  vertex_.reserve(num);
84  Vec3D v;
85  for (unsigned i = 0; i < num; i++)
86  {
87  file >> v.X >> v.Y >> v.Z;
88  vertex_.push_back(v);
89  }
90  //read faces
91  file >> dummy >> num >> dummy;
92  face_.reserve(num);
93  Face f;
94  unsigned id0, id1, id2;
95  for (unsigned i = 0; i < num; i++)
96  {
97  file >> dummy >> id0 >> id1 >> id2;
98  f.vertex[0] = &vertex_[id0];
99  f.vertex[1] = &vertex_[id1];
100  f.vertex[2] = &vertex_[id2];
101  face_.push_back(f);
102  }
103  file >> dummy;
104  //set normals and positions
105  for (auto& face: face_)
106  {
107  face.normal = Vec3D::getUnitVector(Vec3D::cross(*face.vertex[1]-*face.vertex[0],*face.vertex[2]-*face.vertex[0]));
108  }
109  //set neighbours
110  for (auto face0=face_.begin(); face0+1!=face_.end(); face0++)
111  {
112  for (auto face1=face0+1; face1!=face_.end(); face1++)
113  {
114  if (face0->vertex[0]==face1->vertex[0]) {
115  if (face0->vertex[1]==face1->vertex[2]) { //edge 0=2
116  face0->neighbor[0]=&*face1;
117  face1->neighbor[2]=&*face0;
118  } else if (face0->vertex[2]==face1->vertex[1]) { //edge 2=0
119  face0->neighbor[2]=&*face1;
120  face1->neighbor[0]=&*face0;
121  }
122  } else if (face0->vertex[0]==face1->vertex[1]) {
123  if (face0->vertex[1]==face1->vertex[0]) { //edge 0=0
124  face0->neighbor[0]=&*face1;
125  face1->neighbor[0]=&*face0;
126  } else if (face0->vertex[2]==face1->vertex[2]) { //edge 2=1
127  face0->neighbor[2]=&*face1;
128  face1->neighbor[1]=&*face0;
129  }
130  } else if (face0->vertex[0]==face1->vertex[2]) {
131  if (face0->vertex[1]==face1->vertex[1]) { //edge 0=1
132  face0->neighbor[0]=&*face1;
133  face1->neighbor[1]=&*face0;
134  } else if (face0->vertex[2]==face1->vertex[0]) { //edge 2=2
135  face0->neighbor[2]=&*face1;
136  face1->neighbor[2]=&*face0;
137  }
138  } else if (face0->vertex[1]==face1->vertex[0]) {
139  if (face0->vertex[2]==face1->vertex[2]) { //edge 1=2
140  face0->neighbor[1]=&*face1;
141  face1->neighbor[2]=&*face0;
142  }
143  } else if (face0->vertex[1]==face1->vertex[1]) {
144  if (face0->vertex[2]==face1->vertex[0]) { //edge 1=0
145  face0->neighbor[1]=&*face1;
146  face1->neighbor[0]=&*face0;
147  }
148  } else if (face0->vertex[1]==face1->vertex[2]) {
149  if (face0->vertex[2]==face1->vertex[1]) { //edge 1=1
150  face0->neighbor[1]=&*face1;
151  face1->neighbor[1]=&*face0;
152  }
153  }
154 
155  }
156  }
157  //set edge normals (inwards facing)
158  for (auto& face: face_)
159  {
160  face.edgeNormal[0] = Vec3D::getUnitVector(Vec3D::cross(face.normal,*face.vertex[1]-*face.vertex[0]));
161  face.edgeNormal[1] = Vec3D::getUnitVector(Vec3D::cross(face.normal,*face.vertex[2]-*face.vertex[1]));
162  face.edgeNormal[2] = Vec3D::getUnitVector(Vec3D::cross(face.normal,*face.vertex[0]-*face.vertex[2]));
163  }
164  file.close();
165 }
Mdouble X
the vector components
Definition: Vector.h:52
static Vec3D getUnitVector(const Vec3D &a)
Returns a unit Vec3D based on a.
Definition: Vector.cc:428
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
std::vector< Vec3D > vertex_
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:255
Mdouble Y
Definition: Vector.h:52
std::vector< Face > face_
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:45
Mdouble Z
Definition: Vector.h:52
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.

Definition at line 252 of file TriangulatedWall.cc.

References face_, and vertex_.

253 {
254  os << "Vertices " << vertex_.size();
255  for (const auto& vertex: vertex_)
256  os << " " << vertex;
257  os << std::endl;
258  //os << "Faces " << face_.size() << std::endl;
259  unsigned counter = 0;
260  for (const auto& face: face_)
261  {
262  os << "Face " << counter++
263  << " vertex " << face.vertex[0] - &vertex_[0]
264  << " " << face.vertex[1] - &vertex_[0]
265  << " " << face.vertex[2] - &vertex_[0]
266  << " neighbor " << (face.neighbor[0]?(face.neighbor[0] - &face_[0]):-1)
267  << " " << (face.neighbor[1]?(face.neighbor[1] - &face_[0]):-1)
268  << " " << (face.neighbor[2]?(face.neighbor[2] - &face_[0]):-1)
269  << " normal " << face.normal
270  << " edgeNormal " << face.edgeNormal[0]
271  << " " << face.edgeNormal[1]
272  << " " << face.edgeNormal[2]
273  << std::endl;
274  }
275 }
std::vector< Vec3D > vertex_
std::vector< Face > face_
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.

Definition at line 405 of file TriangulatedWall.cc.

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

406 {
407  const int s = vtk.points.size();
408  vtk.points.reserve(s+vertex_.size());
409  for (auto v : vertex_) {
410  vtk.points.push_back(v);
411  }
412  for (auto f : face_) {
413  std::vector<double> cell;
414  cell.reserve(3);
415  cell.push_back(s+f.vertex[0]- &vertex_[0]);
416  cell.push_back(s+f.vertex[1]- &vertex_[0]);
417  cell.push_back(s+f.vertex[2]- &vertex_[0]);
418  vtk.triangleStrips.push_back(cell);
419  }
420 }
std::vector< Vec3D > vertex_
std::vector< std::vector< double > > triangleStrips
Definition: BaseWall.h:36
std::vector< Vec3D > points
Definition: BaseWall.h:35
std::vector< Face > face_

Member Data Documentation

std::vector<Face> TriangulatedWall::face_
private

stores the face properties

Definition at line 149 of file TriangulatedWall.h.

Referenced by getDistanceAndNormal(), getInteractionWith(), readVTK(), TriangulatedWall(), write(), and writeVTK().

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

stores the vertex coordinates

Definition at line 144 of file TriangulatedWall.h.

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


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