MercuryDPM  Alpha
 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:

Classes

struct  normalAndPosition
 

Public Member Functions

 IntersectionOfWalls ()
 Default constructor. More...
 
 IntersectionOfWalls (const IntersectionOfWalls &other)
 Copy constructor. More...
 
 IntersectionOfWalls (std::vector< normalAndPosition > walls, const ParticleSpecies *species)
 Constructor setting values. 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 () override
 Removes all parts of the walls. More...
 
void setSpecies (const ParticleSpecies *species)
 sets species of subwalls as well More...
 
void setHandler (WallHandler *wallHandler) override
 A function which sets the WallHandler for this BaseWall. 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) override
 Move the IntersectionOfWalls to a new position, which is a Vec3D from the old position. More...
 
void read (std::istream &is) override
 Reads an IntersectionOfWalls from an input stream, for example a restart file. More...
 
void write (std::ostream &os) const override
 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...
 
void writeVTK (VTKContainer &vtk) const override
 
- 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...
 
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
 
virtual std::vector
< BaseInteraction * > 
getInteractionWith (BaseParticle *p, Mdouble timeStamp, InteractionHandler *interactionHandler)
 
- 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...
 

Protected Attributes

std::vector< InfiniteWallwallObjects_
 The wall "segments"/directions that together make up the finite wall. More...
 

Private Attributes

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

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 IntersectionOfWalls is convex polygon defined as an intersection of InfiniteWall's.

It 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 each wall, specify a normal and position vector
w.addObject(Vec3D(-1.0,0.0,0.0),Vec3D(1.0,0.0,0.0));
w.addObject(Vec3D(1.0,0.0,0.0),Vec3D(0.0,0.0,0.0));
w.addObject(Vec3D(0.0,-1.0,0.0),Vec3D(0.0,1.0,0.0));
w.addObject(Vec3D(0.0,1.0,0.0),Vec3D(0.0,0.0,0.0));
wallHandler.copyAndAddObject(w);

A particle of radius r and position x touches an InfiniteWall with normal n and position p if $p-n\cdot x\leq r$ (note 'touching particles' also includes particles that are completely enclosed inside the wall). A particle touches an IntersectionOfWalls if it touches all InfiniteWall objects (shown in the image below).

For a demonstration on how to use this class, see T8: Motion of a particle in a box with an obstacle and Flow through a 3D hourglass/silo (shown in the image below).

Definition at line 55 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, BaseWall::getHandler(), logger, and wallObjects_.

39  : BaseWall(other)
40 {
41  wallObjects_ = other.wallObjects_;
42  for (auto& wall : wallObjects_)
43  {
44  if (getHandler() != nullptr)
45  wall.setHandler(getHandler());
46  }
47  A_ = other.A_;
48  AB_ = other.AB_;
49  C_ = other.C_;
50  logger(DEBUG, "IntersectionOfWalls(IntersectionOfWalls&) constructed.");
51 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
std::vector< InfiniteWall > wallObjects_
The wall "segments"/directions that together make up the finite wall.
WallHandler * getHandler() const
A function which returns the WallHandler that handles this BaseWall.
Definition: BaseWall.cc:85
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 ( std::vector< normalAndPosition walls,
const ParticleSpecies species 
)

Constructor setting values.

Definition at line 53 of file IntersectionOfWalls.cc.

References addObject(), and setSpecies().

54 {
55  setSpecies(species);
56  for (auto wall : walls)
57  {
58  addObject(wall.normal,wall.position);
59  }
60 }
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 setSpecies(const ParticleSpecies *species)
sets species of subwalls as well
IntersectionOfWalls::~IntersectionOfWalls ( )
virtual

Destructor.

Definition at line 63 of file IntersectionOfWalls.cc.

References DEBUG, and logger.

64 {
65  logger(DEBUG, "~IntersectionOfWalls() has been called.");
66 }
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 123 of file IntersectionOfWalls.cc.

References A_, AB_, C_, Vec3D::cross(), Vec3D::dot(), InfiniteWall::getNormal(), BaseInteractable::getPosition(), BaseInteractable::getSpecies(), logger, Vec3D::normalize(), InfiniteWall::set(), BaseWall::setSpecies(), VERBOSE, wallObjects_, X, Y, and Z.

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

124 {
125  normal.normalize();
126 
127  //n is the index of the new wall
128  std::size_t n = wallObjects_.size();
129  InfiniteWall w;
130  if (getSpecies()) w.setSpecies(getSpecies());
131  w.set(normal, point);
132  wallObjects_.push_back(w);
133 
134  // AB[n*(n-1)/2+m] is the direction of the intersecting line between walls m and n, m<n
135  // A[n*(n-1)/2+m] is a point on the intersecting line between walls m and n, m<n
136  // See http://www.netcomuk.co.uk/~jenolive/vect18d.html for finding the line where two planes meet
137  AB_.resize(n * (n + 1) / 2); //0 + 1 + 2 + ... + indexNew, total number of walls you need
138  A_.resize(n * (n + 1) / 2);
139  for (std::size_t m = 0; m < n; m++)
140  {
141  std::size_t id = (n - 1) * n / 2 + m;
142  //first we cross the wall normals and normalize to obtain AB
143  AB_[id] = Vec3D::cross(wallObjects_[m].getNormal(), wallObjects_[n].getNormal());
144  AB_[id] /= sqrt(AB_[id].getLengthSquared());
145  //then we find a point A (using AB*x=0 as a third plane)
146  Mdouble invdet = 1.0 / (+wallObjects_[n].getNormal().X * (wallObjects_[m].getNormal().Y * AB_[id].Z - AB_[id].Y * wallObjects_[m].getNormal().Z)
147  - wallObjects_[n].getNormal().Y * (wallObjects_[m].getNormal().X * AB_[id].Z - wallObjects_[m].getNormal().Z * AB_[id].X)
148  + wallObjects_[n].getNormal().Z * (wallObjects_[m].getNormal().X * AB_[id].Y - wallObjects_[m].getNormal().Y * AB_[id].X));
149 
150  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())
151  -(wallObjects_[n].getNormal().Y * AB_[id].Z - wallObjects_[n].getNormal().Z * AB_[id].Y) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
152  +(wallObjects_[n].getNormal().Y * wallObjects_[m].getNormal().Z - wallObjects_[n].getNormal().Z * wallObjects_[m].getNormal().Y) * 0.0,
153  -(wallObjects_[m].getNormal().X * AB_[id].Z - wallObjects_[m].getNormal().Z * AB_[id].X) * Vec3D::dot(wallObjects_[n].getPosition(),wallObjects_[n].getNormal())
154  +(wallObjects_[n].getNormal().X * AB_[id].Z - wallObjects_[n].getNormal().Z * AB_[id].X) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
155  -(wallObjects_[n].getNormal().X * wallObjects_[m].getNormal().Z - wallObjects_[m].getNormal().X * wallObjects_[n].getNormal().Z) * 0.0,
156  +(wallObjects_[m].getNormal().X * AB_[id].Y - AB_[id].X * wallObjects_[m].getNormal().Y) * Vec3D::dot(wallObjects_[n].getPosition(),wallObjects_[n].getNormal())
157  -(wallObjects_[n].getNormal().X * AB_[id].Y - AB_[id].X * wallObjects_[n].getNormal().Y) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
158  +(wallObjects_[n].getNormal().X * wallObjects_[m].getNormal().Y - wallObjects_[m].getNormal().X * wallObjects_[n].getNormal().Y) * 0.0) * invdet;
159  }
160 
161  // C[(n-2)*(n-1)*n/6+(m-1)*m/2+l] is a point intersecting walls l, m and n, l<m<n
162  C_.resize((n - 1) * n * (n + 1) / 6);
163  for (std::size_t m = 0; m < n; m++)
164  {
165  for (std::size_t l = 0; l < m; l++)
166  {
167  std::size_t id = (n - 2) * (n - 1) * n / 6 + (m - 1) * m / 2 + l;
168  Mdouble invdet = 1.0 / (+wallObjects_[n].getNormal().X * (wallObjects_[m].getNormal().Y * wallObjects_[l].getNormal().Z - wallObjects_[l].getNormal().Y * wallObjects_[m].getNormal().Z)
169  - wallObjects_[n].getNormal().Y * (wallObjects_[m].getNormal().X * wallObjects_[l].getNormal().Z - wallObjects_[m].getNormal().Z * wallObjects_[l].getNormal().X)
170  + wallObjects_[n].getNormal().Z * (wallObjects_[m].getNormal().X * wallObjects_[l].getNormal().Y - wallObjects_[m].getNormal().Y * wallObjects_[l].getNormal().X));
171  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())
172  - (wallObjects_[n].getNormal().Y * wallObjects_[l].getNormal().Z - wallObjects_[n].getNormal().Z * wallObjects_[l].getNormal().Y) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
173  + (wallObjects_[n].getNormal().Y * wallObjects_[m].getNormal().Z - wallObjects_[n].getNormal().Z * wallObjects_[m].getNormal().Y) * Vec3D::dot(wallObjects_[l].getPosition(),wallObjects_[l].getNormal()),
174  -(wallObjects_[m].getNormal().X * wallObjects_[l].getNormal().Z - wallObjects_[m].getNormal().Z * wallObjects_[l].getNormal().X) * Vec3D::dot(wallObjects_[n].getPosition(),wallObjects_[n].getNormal())
175  + (wallObjects_[n].getNormal().X * wallObjects_[l].getNormal().Z - wallObjects_[n].getNormal().Z * wallObjects_[l].getNormal().X) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
176  - (wallObjects_[n].getNormal().X * wallObjects_[m].getNormal().Z - wallObjects_[m].getNormal().X * wallObjects_[n].getNormal().Z) * Vec3D::dot(wallObjects_[l].getPosition(),wallObjects_[l].getNormal()),
177  +(wallObjects_[m].getNormal().X * wallObjects_[l].getNormal().Y - wallObjects_[l].getNormal().X * wallObjects_[m].getNormal().Y) * Vec3D::dot(wallObjects_[n].getPosition(),wallObjects_[n].getNormal())
178  - (wallObjects_[n].getNormal().X * wallObjects_[l].getNormal().Y - wallObjects_[l].getNormal().X * wallObjects_[n].getNormal().Y) * Vec3D::dot(wallObjects_[m].getPosition(),wallObjects_[m].getNormal())
179  + (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;
180  }
181  }
182 
183  logger(VERBOSE, "%", *this);
184  for (InfiniteWall w : wallObjects_)
185  logger(VERBOSE, "wallObject %, %", w.getNormal(), w.getPosition());
186  for (Vec3D v : A_)
187  logger(VERBOSE, "A %", v);
188  for (Vec3D v : AB_)
189  logger(VERBOSE, "AB %", v);
190  for (Vec3D v : C_)
191  logger(VERBOSE, "C %", v);
192 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
void normalize()
Makes this Vec3D unit length.
Definition: Vector.cc:214
const ParticleSpecies * getSpecies() const
Returns a pointer to the species of this BaseInteractable.
double Mdouble
Vec3D getNormal() const
Access function for normal.
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:167
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:255
std::vector< Vec3D > AB_
A vector that stores the direction of the intersecting lines between two different InfiniteWall...
void set(Vec3D normal, Vec3D point)
Defines a standard wall, given an outward normal vector s.t. normal*x=normal*point for all x of the w...
Definition: InfiniteWall.cc:79
This is a class defining walls.
Definition: InfiniteWall.h:47
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:45
void setSpecies(const ParticleSpecies *species)
Define the species of this wall.
Definition: BaseWall.cc:113
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 198 of file IntersectionOfWalls.cc.

References addObject(), logger, and WARN.

199 {
200  logger(WARN, "This function is deprecated, use IntersectionOfWalls::addObject(Vec3D, Vec3D) instead.");
201  addObject(normal, position*normal);
202 }
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 ( )
overridevirtual

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 98 of file IntersectionOfWalls.cc.

References A_, AB_, C_, and wallObjects_.

Referenced by createOpenPrism().

99 {
100  wallObjects_.clear();
101  A_.clear();
102  AB_.clear();
103  C_.clear();
104 }
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 93 of file IntersectionOfWalls.cc.

References IntersectionOfWalls().

Referenced by operator=().

94 {
95  return new IntersectionOfWalls(*this);
96 }
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 212 of file IntersectionOfWalls.cc.

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

Referenced by createOpenPrism(), and createPrism().

213 {
214  clear();
215  //note: use i+1 < points.size() instead of i < points.size()-1, otherwise it creates havoc if point has zero entries.
216  for (unsigned int i = 0; i+1 < points.size(); i++)
217  addObject(Vec3D::cross(points[i] - points[i + 1], prismAxis), points[i]);
218 }
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.
MERCURY_DEPRECATED void clear() override
Removes all parts of the walls.
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:255
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 240 of file IntersectionOfWalls.cc.

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

241 {
242  Vec3D prismAxis = Vec3D::cross(
243  Vec3D::getUnitVector(points[1] - points[0]),
244  Vec3D::getUnitVector(points[2] - points[0]));
245  createOpenPrism(points, prismAxis);
246 }
static Vec3D getUnitVector(const Vec3D &a)
Returns a unit Vec3D based on a.
Definition: Vector.cc:428
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:255
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 227 of file IntersectionOfWalls.cc.

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

Referenced by createPrism().

228 {
229  createOpenPrism(points, prismAxis);
230  addObject(Vec3D::cross(points.back() - points.front(), prismAxis), points.front());
231 }
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:255
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 255 of file IntersectionOfWalls.cc.

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

256 {
257  Vec3D prismAxis = Vec3D::cross(
258  Vec3D::getUnitVector(points[1] - points[0]),
259  Vec3D::getUnitVector(points[2] - points[0]));
260  createPrism(points, prismAxis);
261 }
static Vec3D getUnitVector(const Vec3D &a)
Returns a unit Vec3D based on a.
Definition: Vector.cc:428
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:255
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 277 of file IntersectionOfWalls.cc.

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

Referenced by AxisymmetricIntersectionOfWalls::getDistanceAndNormal().

278 {
279  return getDistanceAndNormal(p.getPosition(), p.getInteractionRadius(), distance, normal_return);
280 }
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 304 of file IntersectionOfWalls.cc.

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

305 {
306  if (wallObjects_.size()==0)
307  {
308  logger(DEBUG,"Empty IntersectionOfWalls");
309  return false;
310  }
311 
312  distance = -1e20;
313  Mdouble distance2 = -1e20;
314  Mdouble distance3 = -1e20;
315  Mdouble distanceCurrent;
316  unsigned int id=0;
317  unsigned int id2=0;
318  unsigned int id3=0;
319 
320  //The object has to touch each wall each wall (distanceCurrent) and keep the minimum distance (distance) and wall index (id)
321  for (unsigned int i = 0; i < wallObjects_.size(); i++)
322  {
323  // Calculate distance to each wall (distanceCurrent);
324  distanceCurrent = wallObjects_[i].getDistance(position);
325  // The object has to touch each wall (distanceCurrent >= wallInteractionRadius), otherwise return false (i.e. no contact)
326  // This means that for each InfiniteWall in wallObjects_, the particle is either "inside"
327  // the wall or touching it. If not, there is no interaction.
328  if (distanceCurrent >= wallInteractionRadius)
329  return false;
330  // Find out which of the InfiniteWalls is interacting with the particle.
331  // Keep the minimum distance (distance) and wall index (id)
332  // and store up to two walls (id2, id3) and their distances (distance2, distance3),
333  // if the possible contact point is near the intersection between id and id2 (and id3)
334  if (distanceCurrent > distance)
335  {
336  if (distance > -wallInteractionRadius)
337  {
338  if (distance2 > -wallInteractionRadius)
339  {
340  distance3 = distance;
341  id3 = id;
342  }
343  else
344  {
345  distance2 = distance;
346  id2 = id;
347  }
348  }
349  distance = distanceCurrent;
350  id = i;
351  }
352  else if (distanceCurrent > -wallInteractionRadius)
353  {
354  if (distance2 > -wallInteractionRadius)
355  {
356  distance3 = distanceCurrent;
357  id3 = i;
358  }
359  else
360  {
361  distance2 = distanceCurrent;
362  id2 = i;
363  }
364  }
365  }
366 
367  //If we are here, the closest wall is id;
368  //if distance2>-P.Radius (and distance3>-P.Radius), the possible contact point
369  // is near the intersection between id and id2 (and id3)
370  if (distance2 > -wallInteractionRadius)
371  {
372  //D is the point on wall id closest to P
373  Vec3D D = position + wallObjects_[id].getNormal() * distance;
374  //If the distance of D to id2 is positive, the contact is with the intersection
375  bool intersection_with_id2 = (wallObjects_[id2].getDistance(D) > 0.0);
376 
377  if (distance3 > -wallInteractionRadius && (wallObjects_[id3].getDistance(D) > 0.0))
378  {
379  if (intersection_with_id2)
380  {
381  //possible contact is with intersection of id,id2,id3
382  //we know id2<id3
383  unsigned int index =
384  (id < id2) ? ((id3 - 2) * (id3 - 1) * id3 / 6 + (id2 - 1) * id2 / 2 + id) :
385  (id < id3) ? ((id3 - 2) * (id3 - 1) * id3 / 6 + (id - 1) * id / 2 + id2) :
386  ((id - 2) * (id - 1) * id / 6 + (id3 - 1) * id3 / 2 + id2);
387  normal_return = position - C_[index];
388  distance = sqrt(normal_return.getLengthSquared());
389  if (distance >= wallInteractionRadius)
390  return false; //no contact
391  normal_return /= -distance;
392  return true; //contact with id,id2,id3
393  }
394  else
395  {
396  intersection_with_id2 = true;
397  distance2 = distance3;
398  id2 = id3;
399  }
400  }
401 
402  if (intersection_with_id2)
403  { //possible contact is with intersection of id,id2
404  unsigned int index = (id > id2) ? ((id - 1) * id / 2 + id2) : ((id2 - 1) * id2 / 2 + id);
405  Vec3D AC = position - A_[index];
406  normal_return = AC - AB_[index] * Vec3D::dot(AC, AB_[index]);
407  distance = sqrt(normal_return.getLengthSquared());
408  if (distance >= wallInteractionRadius)
409  return false; //no contact
410  normal_return /= -distance;
411  return true; //contact with two walls
412  }
413  }
414  //contact is with id
415  normal_return = wallObjects_[id].getNormal();
416  return true;
417 }
static Mdouble getLengthSquared(const Vec3D &a)
Calculates the squared length of a Vec3D: .
Definition: Vector.cc:291
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
double Mdouble
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:167
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
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 483 of file IntersectionOfWalls.cc.

484 {
485  return "IntersectionOfWalls";
486 }
void IntersectionOfWalls::move ( const Vec3D move)
overridevirtual

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 430 of file IntersectionOfWalls.cc.

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

431 {
433  for(Vec3D& a : A_)
434  {
435  a += move;
436  }
437  for(Vec3D& c : C_)
438  {
439  c += move;
440  }
441  for(InfiniteWall& o : wallObjects_)
442  {
443  o.move(move);
444  }
445 }
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.
void move(const Vec3D &move) override
Move the IntersectionOfWalls to a new position, which is a Vec3D from the old position.
This is a class defining walls.
Definition: InfiniteWall.h:47
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 80 of file IntersectionOfWalls.cc.

References copy(), DEBUG, and logger.

81 {
82  logger(DEBUG, "IntersectionOfWalls::operator= called.");
83  if (this == &other)
84  {
85  return *this;
86  }
87  return *(other.copy());
88 }
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)
overridevirtual

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 450 of file IntersectionOfWalls.cc.

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

Referenced by AxisymmetricIntersectionOfWalls::read().

451 {
452  BaseWall::read(is);
453  std::string dummy;
454  int n;
455  is >> dummy >> n;
456 
457  Vec3D normal;
458  Vec3D position;
459  for (int i = 0; i < n; i++)
460  {
461  is >> dummy >> normal >> dummy >> position;
462  addObject(normal, position);
463  }
464 }
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::setHandler ( WallHandler handler)
overridevirtual

A function which sets the WallHandler for this BaseWall.

Parameters
[in]handlerA pointer to the BaseHandler that handles this wall.

Reimplemented from BaseWall.

Definition at line 106 of file IntersectionOfWalls.cc.

References BaseWall::setHandler(), and wallObjects_.

107 {
108  BaseWall::setHandler(wallHandler);
109  for (InfiniteWall w : wallObjects_)
110  {
111  BaseWall::setHandler(wallHandler);
112  }
113 }
virtual void setHandler(WallHandler *handler)
A function which sets the WallHandler for this BaseWall.
Definition: BaseWall.cc:76
std::vector< InfiniteWall > wallObjects_
The wall "segments"/directions that together make up the finite wall.
This is a class defining walls.
Definition: InfiniteWall.h:47
void IntersectionOfWalls::setSpecies ( const ParticleSpecies species)

sets species of subwalls as well

Definition at line 68 of file IntersectionOfWalls.cc.

References BaseWall::setSpecies(), and wallObjects_.

Referenced by ChuteWithHopper::addHopper(), IntersectionOfWalls(), WallHandler::readObject(), and WallHandler::readOldObject().

69 {
70  BaseWall::setSpecies(species);
71  for (auto wall : wallObjects_)
72  {
73  wall.setSpecies(species);
74  }
75 }
std::vector< InfiniteWall > wallObjects_
The wall "segments"/directions that together make up the finite wall.
void setSpecies(const ParticleSpecies *species)
Define the species of this wall.
Definition: BaseWall.cc:113
void IntersectionOfWalls::write ( std::ostream &  os) const
overridevirtual

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 470 of file IntersectionOfWalls.cc.

References wallObjects_, and BaseWall::write().

Referenced by AxisymmetricIntersectionOfWalls::write().

471 {
472  BaseWall::write(os);
473  os << " numIntersectionOfWalls " << wallObjects_.size();
474  for (std::vector<InfiniteWall>::const_iterator it = wallObjects_.begin(); it != wallObjects_.end(); ++it)
475  {
476  os << " normal " << it->getNormal() << " position " << it->getPosition();
477  }
478 }
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
void IntersectionOfWalls::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 488 of file IntersectionOfWalls.cc.

References BaseWall::intersectVTK(), and wallObjects_.

489 {
490  for (auto wall=wallObjects_.begin(); wall!=wallObjects_.end(); wall++)
491  {
492  std::vector<Vec3D> points;
493  wall->createVTK (points);
494  for (auto other=wallObjects_.begin(); other!=wallObjects_.end(); other++)
495  {
496  if (other!=wall)
497  {
498  intersectVTK(points, -other->getNormal(), other->getPosition());
499  }
500  }
501  wall->addToVTK (points, vtk);
502  }
503 }
void intersectVTK(std::vector< Vec3D > &points, const Vec3D normal, const Vec3D position) const
Definition: BaseWall.cc:163
std::vector< InfiniteWall > wallObjects_
The wall "segments"/directions that together make up the finite wall.

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 184 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 191 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 198 of file IntersectionOfWalls.h.

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

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

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 176 of file IntersectionOfWalls.h.

Referenced by addObject(), clear(), getDistanceAndNormal(), IntersectionOfWalls(), move(), setHandler(), setSpecies(), write(), AxisymmetricIntersectionOfWalls::writeVTK(), and writeVTK().


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