MercuryDPM  Beta
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
IntersectionOfWalls Class Reference

A IntersectionOfWalls is convex polygon defined as an intersection of InfiniteWall's. More...

#include <IntersectionOfWalls.h>

+ Inheritance diagram for IntersectionOfWalls:

Public Member Functions

 IntersectionOfWalls ()
 Default constructor. More...
 
 IntersectionOfWalls (const IntersectionOfWalls &other)
 Copy constructor. More...
 
virtual ~IntersectionOfWalls ()
 Destructor. More...
 
IntersectionOfWallsoperator= (const IntersectionOfWalls &other)
 
IntersectionOfWallscopy () const override
 Wall copy method. It calls the copy constructor of this Wall, useful for polymorphism. More...
 
MERCURY_DEPRECATED void clear ()
 Removes all parts of the walls. More...
 
void addObject (Vec3D normal, Vec3D point)
 Adds a wall to the set of infinite walls, given an outward normal vector s.t. normal*x=normal*point. More...
 
MERCURY_DEPRECATED void addObject (Vec3D normal, Mdouble position)
 Adds a wall to the set of finite walls, given an outward normal vector s. t. normal*x=position. More...
 
void createOpenPrism (std::vector< Vec3D > points, Vec3D prismAxis)
 Creates an open prism which is a polygon between the points, except the first and last point, and extends infinitely in the PrismAxis direction. More...
 
void createPrism (std::vector< Vec3D > points, Vec3D prismAxis)
 Creates an open prism which is a polygon between the points and extends infinitely in the PrismAxis direction. More...
 
void createOpenPrism (std::vector< Vec3D > points)
 Creates an open prism which is a polygon between the points, except the first and last point, and extends infinitely in the direction perpendicular to the first and second wall. More...
 
void createPrism (std::vector< Vec3D > points)
 Creates an open prism which is a polygon between the points and extends infinitely in the direction perpendicular to the first and second wall. 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 getDistanceAndNormal (const Vec3D &postition, Mdouble wallInteractionRadius, Mdouble &distance, Vec3D &normal_return) const
 Compute the distance from the wall for a given BaseParticle and return if there is an interaction. If there is an interaction, also return the normal vector. More...
 
void move (const Vec3D &move)
 Move the IntersectionOfWalls to a new position, which is a Vec3D from the old position. More...
 
void read (std::istream &is)
 Reads an IntersectionOfWalls from an input stream, for example a restart file. More...
 
void write (std::ostream &os) const
 Writes an IntersectionOfWalls to an output stream, for example a restart file. More...
 
std::string getName () const override
 Returns the name of the object, here the string "IntersectionOfWalls". More...
 
BaseInteractiongetInteractionWith (BaseParticle *p, Mdouble timeStamp, InteractionHandler *interactionHandler)
 Get the interaction between this IntersectionOfWalls 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...
 
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...
 
- 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 (Vec3D force)
 Sets the force on this BaseInteractable. More...
 
void setTorque (Vec3D torque)
 Sets the torque on this BaseInteractable. More...
 
void addForce (Vec3D addForce)
 Adds an amount to the force on this BaseInteractable. More...
 
void addTorque (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 (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 (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 (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 (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< InfiniteWallwallObjects_
 The wall "segments"/directions that together make up the finite wall. More...
 
std::vector< Vec3DA_
 A vector that stores a point for each intersecting line between two different InfiniteWall. More...
 
std::vector< Vec3DAB_
 A vector that stores the direction of the intersecting lines between two different InfiniteWall. More...
 
std::vector< Vec3DC_
 A vector that stores the intersection point of three different InfiniteWall. More...
 

Detailed Description

A IntersectionOfWalls is convex polygon defined as an intersection of InfiniteWall's.

A IntersectionOfWalls can be defined as the intersection of a set of InfiniteWalls, defined by the normal vector into the wall and a point on the wall. For an example see Tutorial T8. A particle touches a finite wall if position_i-normal_i*x<=radius for all InfiniteWall's i.

finiteWall.jpg

Definition at line 44 of file IntersectionOfWalls.h.

Constructor & Destructor Documentation

IntersectionOfWalls::IntersectionOfWalls ( )

Default constructor.

Definition at line 30 of file IntersectionOfWalls.cc.

References DEBUG, and logger.

Referenced by copy().

31 {
32  logger(DEBUG, "IntersectionOfWalls() constructed.");
33 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
IntersectionOfWalls::IntersectionOfWalls ( const IntersectionOfWalls other)

Copy constructor.

Parameters
[in]otherThe IntersectionOfWalls that must be copied.

Definition at line 38 of file IntersectionOfWalls.cc.

References A_, AB_, C_, DEBUG, logger, and wallObjects_.

39  : BaseWall(other)
40 {
41  wallObjects_ = other.wallObjects_;
42  A_ = other.A_;
43  AB_ = other.AB_;
44  C_ = other.C_;
45  logger(DEBUG, "IntersectionOfWalls(IntersectionOfWalls&) constructed.");
46 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
std::vector< InfiniteWall > wallObjects_
The wall "segments"/directions that together make up the finite wall.
std::vector< Vec3D > C_
A vector that stores the intersection point of three different InfiniteWall.
std::vector< Vec3D > A_
A vector that stores a point for each intersecting line between two different InfiniteWall.
BaseWall()
Default constructor. It makes an empty BaseWall.
Definition: BaseWall.cc:31
std::vector< Vec3D > AB_
A vector that stores the direction of the intersecting lines between two different InfiniteWall...
IntersectionOfWalls::~IntersectionOfWalls ( )
virtual

Destructor.

Definition at line 48 of file IntersectionOfWalls.cc.

References DEBUG, and logger.

49 {
50  logger(DEBUG, "~IntersectionOfWalls() has been called.");
51 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")

Member Function Documentation

void IntersectionOfWalls::addObject ( Vec3D  normal,
Vec3D  point 
)

Adds a wall to the set of infinite walls, given an outward normal vector s.t. normal*x=normal*point.

Parameters
[in]normalThe normal to this wallObject.
[in]pointOne of the points of the wallObject.

Adds a wall to the set of finite walls, given an outward unit normal vector s.t. normal*x=normal*point for all x of the wallObject. First make the InfiniteWall, then compute all intersections, which are then stored in A_, AB_ and C_.

Definition at line 90 of file IntersectionOfWalls.cc.

References A_, AB_, C_, Vec3D::cross(), Vec3D::dot(), BaseInteractable::getPosition(), logger, Vec3D::normalize(), VERBOSE, wallObjects_, X, Y, and Z.

Referenced by ChuteWithHopper::addHopper(), addObject(), createOpenPrism(), createPrism(), read(), and WallHandler::readOldObject().

91 {
92  normal.normalize();
93 
94  //n is the index of the new wall
95  std::size_t n = wallObjects_.size();
96  wallObjects_.resize(n + 1);
97  wallObjects_[n].set(normal, point);
98 
99  // AB[n*(n-1)/2+m] is the direction of the intersecting line between walls m and n, m<n
100  // A[n*(n-1)/2+m] is a point on the intersecting line between walls m and n, m<n
101  // See http://www.netcomuk.co.uk/~jenolive/vect18d.html for finding the line where two planes meet
102  AB_.resize(n * (n + 1) / 2); //0 + 1 + 2 + ... + indexNew, total number of walls you need
103  A_.resize(n * (n + 1) / 2);
104  for (std::size_t m = 0; m < n; m++)
105  {
106  std::size_t id = (n - 1) * n / 2 + m;
107  //first we cross the wall normals and normalize to obtain AB
108  AB_[id] = Vec3D::cross(wallObjects_[m].getNormal(), wallObjects_[n].getNormal());
109  AB_[id] /= sqrt(AB_[id].getLengthSquared());
110  //then we find a point A (using AB*x=0 as a third plane)
111  Mdouble invdet = 1.0 / (+wallObjects_[n].getNormal().X * (wallObjects_[m].getNormal().Y * AB_[id].Z - AB_[id].Y * wallObjects_[m].getNormal().Z)
112  - wallObjects_[n].getNormal().Y * (wallObjects_[m].getNormal().X * AB_[id].Z - wallObjects_[m].getNormal().Z * AB_[id].X)
113  + wallObjects_[n].getNormal().Z * (wallObjects_[m].getNormal().X * AB_[id].Y - wallObjects_[m].getNormal().Y * AB_[id].X));
114 
115  A_[id] = Vec3D( +(wallObjects_[m].getNormal().Y * AB_[id].Z - AB_[id].Y * wallObjects_[m].getNormal().Z) * Vec3D::dot(wallObjects_[n].getPosition(),wallObjects_[n].getNormal())
116  -(wallObjects_[n].getNormal().Y * AB_[id].Z - wallObjects_[n].getNormal().Z * AB_[id].Y) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
117  +(wallObjects_[n].getNormal().Y * wallObjects_[m].getNormal().Z - wallObjects_[n].getNormal().Z * wallObjects_[m].getNormal().Y) * 0.0,
118  -(wallObjects_[m].getNormal().X * AB_[id].Z - wallObjects_[m].getNormal().Z * AB_[id].X) * Vec3D::dot(wallObjects_[n].getPosition(),wallObjects_[n].getNormal())
119  +(wallObjects_[n].getNormal().X * AB_[id].Z - wallObjects_[n].getNormal().Z * AB_[id].X) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
120  -(wallObjects_[n].getNormal().X * wallObjects_[m].getNormal().Z - wallObjects_[m].getNormal().X * wallObjects_[n].getNormal().Z) * 0.0,
121  +(wallObjects_[m].getNormal().X * AB_[id].Y - AB_[id].X * wallObjects_[m].getNormal().Y) * Vec3D::dot(wallObjects_[n].getPosition(),wallObjects_[n].getNormal())
122  -(wallObjects_[n].getNormal().X * AB_[id].Y - AB_[id].X * wallObjects_[n].getNormal().Y) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
123  +(wallObjects_[n].getNormal().X * wallObjects_[m].getNormal().Y - wallObjects_[m].getNormal().X * wallObjects_[n].getNormal().Y) * 0.0) * invdet;
124  }
125 
126  // C[(n-2)*(n-1)*n/6+(m-1)*m/2+l] is a point intersecting walls l, m and n, l<m<n
127  C_.resize((n - 1) * n * (n + 1) / 6);
128  for (std::size_t m = 0; m < n; m++)
129  {
130  for (std::size_t l = 0; l < m; l++)
131  {
132  std::size_t id = (n - 2) * (n - 1) * n / 6 + (m - 1) * m / 2 + l;
133  Mdouble invdet = 1.0 / (+wallObjects_[n].getNormal().X * (wallObjects_[m].getNormal().Y * wallObjects_[l].getNormal().Z - wallObjects_[l].getNormal().Y * wallObjects_[m].getNormal().Z)
134  - wallObjects_[n].getNormal().Y * (wallObjects_[m].getNormal().X * wallObjects_[l].getNormal().Z - wallObjects_[m].getNormal().Z * wallObjects_[l].getNormal().X)
135  + wallObjects_[n].getNormal().Z * (wallObjects_[m].getNormal().X * wallObjects_[l].getNormal().Y - wallObjects_[m].getNormal().Y * wallObjects_[l].getNormal().X));
136  C_[id] = Vec3D(+(wallObjects_[m].getNormal().Y * wallObjects_[l].getNormal().Z - wallObjects_[l].getNormal().Y * wallObjects_[m].getNormal().Z) * Vec3D::dot(wallObjects_[n].getPosition(),wallObjects_[n].getNormal())
137  - (wallObjects_[n].getNormal().Y * wallObjects_[l].getNormal().Z - wallObjects_[n].getNormal().Z * wallObjects_[l].getNormal().Y) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
138  + (wallObjects_[n].getNormal().Y * wallObjects_[m].getNormal().Z - wallObjects_[n].getNormal().Z * wallObjects_[m].getNormal().Y) * Vec3D::dot(wallObjects_[l].getPosition(),wallObjects_[l].getNormal()),
139  -(wallObjects_[m].getNormal().X * wallObjects_[l].getNormal().Z - wallObjects_[m].getNormal().Z * wallObjects_[l].getNormal().X) * Vec3D::dot(wallObjects_[n].getPosition(),wallObjects_[n].getNormal())
140  + (wallObjects_[n].getNormal().X * wallObjects_[l].getNormal().Z - wallObjects_[n].getNormal().Z * wallObjects_[l].getNormal().X) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
141  - (wallObjects_[n].getNormal().X * wallObjects_[m].getNormal().Z - wallObjects_[m].getNormal().X * wallObjects_[n].getNormal().Z) * Vec3D::dot(wallObjects_[l].getPosition(),wallObjects_[l].getNormal()),
142  +(wallObjects_[m].getNormal().X * wallObjects_[l].getNormal().Y - wallObjects_[l].getNormal().X * wallObjects_[m].getNormal().Y) * Vec3D::dot(wallObjects_[n].getPosition(),wallObjects_[n].getNormal())
143  - (wallObjects_[n].getNormal().X * wallObjects_[l].getNormal().Y - wallObjects_[l].getNormal().X * wallObjects_[n].getNormal().Y) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
144  + (wallObjects_[n].getNormal().X * wallObjects_[m].getNormal().Y - wallObjects_[m].getNormal().X * wallObjects_[n].getNormal().Y) * Vec3D::dot(wallObjects_[l].getPosition(),wallObjects_[l].getNormal())) * invdet;
145  }
146  }
147 
148  logger(VERBOSE, "%", *this);
149  for (InfiniteWall w : wallObjects_)
150  logger(VERBOSE, "wallObject %, %", w.getNormal(), w.getPosition());
151  for (Vec3D v : A_)
152  logger(VERBOSE, "A %", v);
153  for (Vec3D v : AB_)
154  logger(VERBOSE, "AB %", v);
155  for (Vec3D v : C_)
156  logger(VERBOSE, "C %", v);
157 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
void normalize()
Makes this Vec3D unit length.
Definition: Vector.cc:234
double Mdouble
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:187
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
std::vector< InfiniteWall > wallObjects_
The wall "segments"/directions that together make up the finite wall.
std::vector< Vec3D > C_
A vector that stores the intersection point of three different InfiniteWall.
std::vector< Vec3D > A_
A vector that stores a point for each intersecting line between two different InfiniteWall.
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:268
std::vector< Vec3D > AB_
A vector that stores the direction of the intersecting lines between two different InfiniteWall...
This is a class defining walls.
Definition: InfiniteWall.h:43
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:45
void IntersectionOfWalls::addObject ( Vec3D  normal,
Mdouble  position 
)

Adds a wall to the set of finite walls, given an outward normal vector s. t. normal*x=position.

Deprecated:
Don't use this function, instead use the function addObject(Vec3D, Vec3D).
Parameters
[in]normalThe normal to the wallObject.
[in]positionThe position of the wallObject in the direction of the normal vector.

Definition at line 163 of file IntersectionOfWalls.cc.

References addObject(), logger, and WARN.

164 {
165  logger(WARN, "This function is deprecated, use IntersectionOfWalls::addObject(Vec3D, Vec3D) instead.");
166  addObject(normal, position*normal);
167 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
void addObject(Vec3D normal, Vec3D point)
Adds a wall to the set of infinite walls, given an outward normal vector s.t. normal*x=normal*point.
void IntersectionOfWalls::clear ( )
virtual

Removes all parts of the walls.

Deprecated:
Please don't use any clear() anymore, it will be gone soon.

Reimplemented from BaseWall.

Definition at line 74 of file IntersectionOfWalls.cc.

References A_, AB_, C_, and wallObjects_.

Referenced by createOpenPrism().

75 {
76  wallObjects_.clear();
77  A_.clear();
78  AB_.clear();
79  C_.clear();
80 }
std::vector< InfiniteWall > wallObjects_
The wall "segments"/directions that together make up the finite wall.
std::vector< Vec3D > C_
A vector that stores the intersection point of three different InfiniteWall.
std::vector< Vec3D > A_
A vector that stores a point for each intersecting line between two different InfiniteWall.
std::vector< Vec3D > AB_
A vector that stores the direction of the intersecting lines between two different InfiniteWall...
IntersectionOfWalls * IntersectionOfWalls::copy ( ) const
overridevirtual

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

Returns
pointer to a IntersectionOfWalls object allocated using new.

Implements BaseWall.

Definition at line 69 of file IntersectionOfWalls.cc.

References IntersectionOfWalls().

Referenced by operator=().

70 {
71  return new IntersectionOfWalls(*this);
72 }
IntersectionOfWalls()
Default constructor.
void IntersectionOfWalls::createOpenPrism ( std::vector< Vec3D points,
Vec3D  prismAxis 
)

Creates an open prism which is a polygon between the points, except the first and last point, and extends infinitely in the PrismAxis direction.

Parameters
[in]pointsA vector of 3D-vectors which contains the points between which the polygon is drawn.
[in]prismAxisA 3D-vector which represents the direction in which the prism is extended infinitely.

Create an open prism which is a polygon with no connection between the first and last point, and extending infinitely in the other direction, which is defined as PrismAxis. Do this by adding the walls between the consecutive points one by one.

Definition at line 177 of file IntersectionOfWalls.cc.

References addObject(), clear(), and Vec3D::cross().

Referenced by createOpenPrism(), and createPrism().

178 {
179  clear();
180  //note: use i+1 < points.size() instead of i < points.size()-1, otherwise it creates havoc if point has zero entries.
181  for (unsigned int i = 0; i+1 < points.size(); i++)
182  addObject(Vec3D::cross(points[i] - points[i + 1], prismAxis), points[i]);
183 }
void addObject(Vec3D normal, Vec3D point)
Adds a wall to the set of infinite walls, given an outward normal vector s.t. normal*x=normal*point.
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:268
MERCURY_DEPRECATED void clear()
Removes all parts of the walls.
void IntersectionOfWalls::createOpenPrism ( std::vector< Vec3D points)

Creates an open prism which is a polygon between the points, except the first and last point, and extends infinitely in the direction perpendicular to the first and second wall.

Parameters
[in]pointsA vector of 3D-vectors which contains the points between which the polygon is drawn.

Create an open prism which is a polygon with no connection between the first and last point, and extending infinitely in the direction perpendicular to the first and second wall. Do this by first computing in which direction the wall must be extended infinitely, then call createOpenPrism(points, prismAxis).

Definition at line 205 of file IntersectionOfWalls.cc.

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

206 {
207  Vec3D prismAxis = Vec3D::cross(
208  Vec3D::getUnitVector(points[2] - points[0]),
209  Vec3D::getUnitVector(points[1] - points[0]));
210  createOpenPrism(points, prismAxis);
211 }
static Vec3D getUnitVector(const Vec3D &a)
Returns a unit Vec3D based on a.
Definition: Vector.cc:441
void createOpenPrism(std::vector< Vec3D > points, Vec3D prismAxis)
Creates an open prism which is a polygon between the points, except the first and last point...
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:268
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:45
void IntersectionOfWalls::createPrism ( std::vector< Vec3D points,
Vec3D  prismAxis 
)

Creates an open prism which is a polygon between the points and extends infinitely in the PrismAxis direction.

Parameters
[in]pointsA vector of 3D-vectors which contains the points between which the polygon is drawn.
[in]prismAxisA 3D-vector which represents the direction in which the prism is extended infinitely.

Create an open prism which is a polygon and extending infinitely in the other direction, which is defined as PrismAxis. Do this by first creating an open prism and then connect the last and the first point.

Definition at line 192 of file IntersectionOfWalls.cc.

References addObject(), createOpenPrism(), and Vec3D::cross().

Referenced by createPrism().

193 {
194  createOpenPrism(points, prismAxis);
195  addObject(Vec3D::cross(points.back() - points.front(), prismAxis), points.front());
196 }
void addObject(Vec3D normal, Vec3D point)
Adds a wall to the set of infinite walls, given an outward normal vector s.t. normal*x=normal*point.
void createOpenPrism(std::vector< Vec3D > points, Vec3D prismAxis)
Creates an open prism which is a polygon between the points, except the first and last point...
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:268
void IntersectionOfWalls::createPrism ( std::vector< Vec3D points)

Creates an open prism which is a polygon between the points and extends infinitely in the direction perpendicular to the first and second wall.

Parameters
[in]pointsA vector of 3D-vectors which contains the points between which the polygon is drawn.

Create an open prism which is a polygon, and extending infinitely in the direction perpendicular to the first and second wall. Do this by first computing in which direction the wall must be extended infinitely, then call createOpenPrism(points, prismAxis).

Definition at line 220 of file IntersectionOfWalls.cc.

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

221 {
222  Vec3D prismAxis = Vec3D::cross(
223  Vec3D::getUnitVector(points[2] - points[0]),
224  Vec3D::getUnitVector(points[1] - points[0]));
225  createPrism(points, prismAxis);
226 }
static Vec3D getUnitVector(const Vec3D &a)
Returns a unit Vec3D based on a.
Definition: Vector.cc:441
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:268
void createPrism(std::vector< Vec3D > points, Vec3D prismAxis)
Creates an open prism which is a polygon between the points and extends infinitely in the PrismAxis d...
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:45
bool IntersectionOfWalls::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 IntersectionOfWalls. If there is a collision, this function also computes the distance between the BaseParticle and IntersectionOfWalls and the normal of the IntersectionOfWalls at the intersection point. It does this by calling IntersectionOfWalls::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.

Implements BaseWall.

Definition at line 242 of file IntersectionOfWalls.cc.

References BaseParticle::getInteractionRadius(), and BaseInteractable::getPosition().

Referenced by AxisymmetricIntersectionOfWalls::getDistanceAndNormal(), and getInteractionWith().

243 {
244  return getDistanceAndNormal(p.getPosition(), p.getInteractionRadius(), distance, normal_return);
245 }
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...
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Mdouble getInteractionRadius() const
Returns the particle's interaction radius, which might be different from radius_ (e.g., when dealing with wet particles)
bool IntersectionOfWalls::getDistanceAndNormal ( const Vec3D position,
Mdouble  wallInteractionRadius,
Mdouble distance,
Vec3D normal_return 
) const

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

Parameters
[in]positionThe position of the object there is possible an interaction with.
[in]wallInteractionRadiusThe maximum distance between the IntersectionOfWalls and the input argument position for which there is an interaction.
[out]distanceThe distance of the object at position to this wall.
[out]normal_returnIf there was an interaction, the normal vector to this wall will be placed here.
Returns
A boolean which says whether or not there was an interaction.

This function computes whether or not there is an interaction between an object at the given distance and this IntersectionOfWalls. If there is an interaction, this function also computes the distance between the BaseParticle and IntersectionOfWalls and the normal of the IntersectionOfWalls at the intersection point. First check if the distance between the object at position and the IntersectionOfWalls is smaller or greater than the wallInteractionRadius. If there is no interaction, return false, the output parameters then have no meaning. If there is an interaction, find out which (one or more) of the InfiniteWall there is an interaction with. Then compute the distance between the particle and InfiniteWall and the normal to the interaction point. 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.

Definition at line 269 of file IntersectionOfWalls.cc.

References A_, AB_, C_, Vec3D::dot(), Vec3D::getLengthSquared(), and wallObjects_.

270 {
271  distance = -1e20;
272  Mdouble distance2 = -1e20;
273  Mdouble distance3 = -1e20;
274  Mdouble distanceCurrent;
275  unsigned int id=0;
276  unsigned int id2=0;
277  unsigned int id3=0;
278 
279  //The object has to touch each wall each wall (distanceCurrent) and keep the minimum distance (distance) and wall index (id)
280  for (unsigned int i = 0; i < wallObjects_.size(); i++)
281  {
282  // Calculate distance to each wall (distanceCurrent);
283  distanceCurrent = wallObjects_[i].getDistance(position);
284  // The object has to touch each wall (distanceCurrent >= wallInteractionRadius), otherwise return false (i.e. no contact)
285  // This means that for each InfiniteWall in wallObjects_, the particle is either "inside"
286  // the wall or touching it. If not, there is no interaction.
287  if (distanceCurrent >= wallInteractionRadius)
288  return false;
289  // Find out which of the InfiniteWalls is interacting with the particle.
290  // Keep the minimum distance (distance) and wall index (id)
291  // and store up to two walls (id2, id3) and their distances (distance2, distance3),
292  // if the possible contact point is near the intersection between id and id2 (and id3)
293  if (distanceCurrent > distance)
294  {
295  if (distance > -wallInteractionRadius)
296  {
297  if (distance2 > -wallInteractionRadius)
298  {
299  distance3 = distance;
300  id3 = id;
301  }
302  else
303  {
304  distance2 = distance;
305  id2 = id;
306  }
307  }
308  distance = distanceCurrent;
309  id = i;
310  }
311  else if (distanceCurrent > -wallInteractionRadius)
312  {
313  if (distance2 > -wallInteractionRadius)
314  {
315  distance3 = distanceCurrent;
316  id3 = i;
317  }
318  else
319  {
320  distance2 = distanceCurrent;
321  id2 = i;
322  }
323  }
324  }
325 
326  //If we are here, the closest wall is id;
327  //if distance2>-P.Radius (and distance3>-P.Radius), the possible contact point
328  // is near the intersection between id and id2 (and id3)
329  if (distance2 > -wallInteractionRadius)
330  {
331  //D is the point on wall id closest to P
332  Vec3D D = position + wallObjects_[id].getNormal() * distance;
333  //If the distance of D to id2 is positive, the contact is with the intersection
334  bool intersection_with_id2 = (wallObjects_[id2].getDistance(D) > 0.0);
335 
336  if (distance3 > -wallInteractionRadius && (wallObjects_[id3].getDistance(D) > 0.0))
337  {
338  if (intersection_with_id2)
339  {
340  //possible contact is with intersection of id,id2,id3
341  //we know id2<id3
342  unsigned int index =
343  (id < id2) ? ((id3 - 2) * (id3 - 1) * id3 / 6 + (id2 - 1) * id2 / 2 + id) :
344  (id < id3) ? ((id3 - 2) * (id3 - 1) * id3 / 6 + (id - 1) * id / 2 + id2) :
345  ((id - 2) * (id - 1) * id / 6 + (id3 - 1) * id3 / 2 + id2);
346  normal_return = position - C_[index];
347  distance = sqrt(normal_return.getLengthSquared());
348  if (distance >= wallInteractionRadius)
349  return false; //no contact
350  normal_return /= -distance;
351  return true; //contact with id,id2,id3
352  }
353  else
354  {
355  intersection_with_id2 = true;
356  distance2 = distance3;
357  id2 = id3;
358  }
359  }
360 
361  if (intersection_with_id2)
362  { //possible contact is with intersection of id,id2
363  unsigned int index = (id > id2) ? ((id - 1) * id / 2 + id2) : ((id2 - 1) * id2 / 2 + id);
364  Vec3D AC = position - A_[index];
365  normal_return = AC - AB_[index] * Vec3D::dot(AC, AB_[index]);
366  distance = sqrt(normal_return.getLengthSquared());
367  if (distance >= wallInteractionRadius)
368  return false; //no contact
369  normal_return /= -distance;
370  return true; //contact with two walls
371  }
372  }
373  //contact is with id
374  normal_return = wallObjects_[id].getNormal();
375  return true;
376 }
static Mdouble getLengthSquared(const Vec3D &a)
Calculates the squared length of a Vec3D: .
Definition: Vector.cc:304
double Mdouble
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:187
std::vector< InfiniteWall > wallObjects_
The wall "segments"/directions that together make up the finite wall.
std::vector< Vec3D > C_
A vector that stores the intersection point of three different InfiniteWall.
std::vector< Vec3D > A_
A vector that stores a point for each intersecting line between two different InfiniteWall.
std::vector< Vec3D > AB_
A vector that stores the direction of the intersecting lines between two different InfiniteWall...
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:45
BaseInteraction * IntersectionOfWalls::getInteractionWith ( BaseParticle p,
Mdouble  timeStamp,
InteractionHandler interactionHandler 
)
virtual

Get the interaction between this IntersectionOfWalls 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.
Todo:
{DK: What is the contact point for interactions with walls}

Implements BaseInteractable.

Definition at line 454 of file IntersectionOfWalls.cc.

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

455 {
456  Mdouble distance;
457  Vec3D normal;
458 
459  if (getDistanceAndNormal(*p, distance, normal))
460  {
461  BaseInteraction* c = interactionHandler->getInteraction(p, this, timeStamp);
462  c->setNormal(-normal);
463  c->setDistance(distance);
464  c->setOverlap(p->getRadius() - distance);
466  c->setContactPoint(p->getPosition() - (p->getRadius() - 0.5 * c->getOverlap()) * c->getNormal());
467  return c;
468  }
469  else
470  {
471  return nullptr;
472  }
473 }
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.
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...
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 ...
const Vec3D & getNormal() const
Gets the normal vector between the two interacting objects.
Mdouble getRadius() const
Returns the particle's radius_.
Mdouble getOverlap() const
Returns a Mdouble with the current overlap between the two interacting objects.
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:45
std::string IntersectionOfWalls::getName ( ) const
overridevirtual

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

Returns
The string "IntersectionOfWalls".

Implements BaseObject.

Definition at line 442 of file IntersectionOfWalls.cc.

443 {
444  return "IntersectionOfWalls";
445 }
void IntersectionOfWalls::move ( const Vec3D move)
virtual

Move the IntersectionOfWalls 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 InterSectionOfWalls 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 IntersectionOfWalls as a whole.

Todo:
We should use the position_ and orientation_ of the IntersectionOfWalls; that way, IntersectionOfWalls 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 389 of file IntersectionOfWalls.cc.

References A_, C_, BaseInteractable::move(), and wallObjects_.

390 {
392  for(Vec3D& a : A_)
393  {
394  a += move;
395  }
396  for(Vec3D& c : C_)
397  {
398  c += move;
399  }
400  for(InfiniteWall& o : wallObjects_)
401  {
402  o.move(move);
403  }
404 }
void move(const Vec3D &move)
Move the IntersectionOfWalls to a new position, which is a Vec3D from the old position.
std::vector< InfiniteWall > wallObjects_
The wall "segments"/directions that together make up the finite wall.
std::vector< Vec3D > C_
A vector that stores the intersection point of three different InfiniteWall.
std::vector< Vec3D > A_
A vector that stores a point for each intersecting line between two different InfiniteWall.
This is a class defining walls.
Definition: InfiniteWall.h:43
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:45
virtual void move(const Vec3D &move)
Moves this BaseInteractable by adding an amount to the position.
IntersectionOfWalls & IntersectionOfWalls::operator= ( const IntersectionOfWalls other)

Copy assignment operator.

Parameters
[in]otherThe IntersectionOfWalls that must be copied.

Definition at line 56 of file IntersectionOfWalls.cc.

References copy(), DEBUG, and logger.

57 {
58  logger(DEBUG, "IntersectionOfWalls::operator= called.");
59  if (this == &other)
60  {
61  return *this;
62  }
63  return *(other.copy());
64 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
IntersectionOfWalls * copy() const override
Wall copy method. It calls the copy constructor of this Wall, useful for polymorphism.
void IntersectionOfWalls::read ( std::istream &  is)
virtual

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

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

Reimplemented from BaseWall.

Definition at line 409 of file IntersectionOfWalls.cc.

References addObject(), and BaseWall::read().

Referenced by AxisymmetricIntersectionOfWalls::read().

410 {
411  BaseWall::read(is);
412  std::string dummy;
413  int n;
414  is >> dummy >> n;
415 
416  Vec3D normal;
417  Vec3D position;
418  for (int i = 0; i < n; i++)
419  {
420  is >> dummy >> normal >> dummy >> position;
421  addObject(normal, position);
422  }
423 }
void addObject(Vec3D normal, Vec3D point)
Adds a wall to the set of infinite walls, given an outward normal vector s.t. normal*x=normal*point.
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:45
void read(std::istream &is)
Function that reads a BaseWall from an input stream, usually a restart file.
Definition: BaseWall.cc:60
void IntersectionOfWalls::write ( std::ostream &  os) const
virtual

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

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

Reimplemented from BaseWall.

Definition at line 429 of file IntersectionOfWalls.cc.

References wallObjects_, and BaseWall::write().

Referenced by AxisymmetricIntersectionOfWalls::write().

430 {
431  BaseWall::write(os);
432  os << " numIntersectionOfWalls " << wallObjects_.size();
433  for (std::vector<InfiniteWall>::const_iterator it = wallObjects_.begin(); it != wallObjects_.end(); ++it)
434  {
435  os << " normal " << it->getNormal() << " position " << it->getPosition();
436  }
437 }
std::vector< InfiniteWall > wallObjects_
The wall "segments"/directions that together make up the finite wall.
void write(std::ostream &os) const
Function that writes a BaseWall to an output stream, usually a restart file.
Definition: BaseWall.cc:68

Member Data Documentation

std::vector<Vec3D> IntersectionOfWalls::A_
private

A vector that stores a point for each intersecting line between two different InfiniteWall.

A[n*(n-1)/2+m] is a point on the intersecting line between walls m and n, m<n.

Definition at line 160 of file IntersectionOfWalls.h.

Referenced by addObject(), clear(), getDistanceAndNormal(), IntersectionOfWalls(), and move().

std::vector<Vec3D> IntersectionOfWalls::AB_
private

A vector that stores the direction of the intersecting lines between two different InfiniteWall.

AB[n*(n-1)/2+m] is the direction of the intersecting line between walls m and n, m<n.

Definition at line 167 of file IntersectionOfWalls.h.

Referenced by addObject(), clear(), getDistanceAndNormal(), and IntersectionOfWalls().

std::vector<Vec3D> IntersectionOfWalls::C_
private

A vector that stores the intersection point of three different InfiniteWall.

C[(n-2)*(n-1)*n/6+(m-1)*m/2+l] is a point intersecting walls l, m and n, l<m<n

Definition at line 174 of file IntersectionOfWalls.h.

Referenced by addObject(), clear(), getDistanceAndNormal(), IntersectionOfWalls(), and move().

std::vector<InfiniteWall> IntersectionOfWalls::wallObjects_
private

The wall "segments"/directions that together make up the finite wall.

An intersection of walls exists of a number of infinite walls that are cut of at the intersection points. These InfiniteWall are saved in this vector called iWObjects_.

Definition at line 153 of file IntersectionOfWalls.h.

Referenced by addObject(), clear(), getDistanceAndNormal(), IntersectionOfWalls(), move(), and write().


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