revision: v0.14
BasicIntersectionOfWalls Class Reference

Restriction of a wall to the intersection with another wall. More...

#include <BasicIntersectionOfWalls.h>

+ Inheritance diagram for BasicIntersectionOfWalls:

Public Member Functions

 BasicIntersectionOfWalls ()
 Default constructor, the normal is infinitely long. More...
 
 BasicIntersectionOfWalls (const BasicIntersectionOfWalls &w)
 Copy constructor, copy the given wall. More...
 
 ~BasicIntersectionOfWalls () override
 Default destructor. More...
 
BasicIntersectionOfWallscopy () const override
 Wall copy method. It calls the copy constructor of this Wall, useful for polymorphism. More...
 
unsigned long getNumberOfObjects ()
 Returns the number of objects. More...
 
void add (BaseWall &wall)
 Defines a standard wall, given an outward normal vector s.t. normal*x=normal*point for all x of the 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...
 
void read (std::istream &is) override
 Reads BasicIntersectionOfWalls from a restart file. More...
 
void oldRead (std::istream &is)
 Reads BasicIntersectionOfWalls from an old-style restart file. More...
 
void write (std::ostream &os) const override
 Writes the BasicIntersectionOfWalls to an output stream, usually a restart file. More...
 
std::string getName () const override
 Returns the name of the object, in this case the string "BasicIntersectionOfWalls". More...
 
void getVTK (std::vector< Vec3D > &points, std::vector< std::vector< double >> &triangleStrips)
 
BaseWallgetObject (unsigned i)
 
- Public Member Functions inherited from BaseWall
 BaseWall ()
 Default constructor. More...
 
 BaseWall (const BaseWall &w)
 Copy constructor. More...
 
 ~BaseWall () override
 Default destructor. More...
 
virtual bool getDistanceNormalOverlap (const BaseParticle &P, Mdouble &distance, Vec3D &normal_return, Mdouble &overlap) const
 
virtual bool getDistanceNormalOverlapSuperquadric (const SuperQuadricParticle &p, Mdouble &distance, Vec3D &normal_return, Mdouble &overlap) const
 
virtual Vec3D getFurthestPointSuperQuadric (const Vec3D &normalBodyFixed, const Vec3D &axes, Mdouble eps1, Mdouble eps2) const
 
virtual void setHandler (WallHandler *handler)
 A function which sets the WallHandler for this BaseWall. More...
 
WallHandlergetHandler () const
 A function which returns the WallHandler that handles this BaseWall. More...
 
void setIndSpecies (unsigned int indSpecies) override
 Define the species of this wall using the index of the species in the SpeciesHandler in this DPMBase. More...
 
void setSpecies (const ParticleSpecies *species)
 Defines the species of the current wall. More...
 
bool isFixed () const override
 
void setForceControl (Vec3D forceGoal, Vec3D gainFactor, Vec3D baseVelocity={0, 0, 0})
 Slowly adjusts the force on a wall towards a specified goal, by adjusting (prescribing) the velocity of the wall. More...
 
virtual bool isLocal (Vec3D &min, Vec3D &max) const
 
bool getLinePlaneIntersect (Vec3D &intersect, const Vec3D &p0, const Vec3D &p1, const Vec3D &n, const Vec3D &p)
 
bool isInsideWallVTK (const Vec3D &point, const Vec3D &normal, const Vec3D &position) const
 
void projectOntoWallVTK (Vec3D &point0, const Vec3D &point1, const Vec3D &normal, const Vec3D &position) const
 
void intersectVTK (std::vector< Vec3D > &points, Vec3D normal, Vec3D position) const
 
virtual BaseInteractiongetInteractionWithSuperQuad (SuperQuadricParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler)
 
virtual void writeVTK (VTKContainer &vtk) const
 
void getVTK (std::vector< Vec3D > &points, std::vector< std::vector< double >> &triangleStrips)
 
const Vec3D getAxis () const
 
BaseInteractiongetInteractionWith (BaseParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler) override
 Returns the interaction between this wall and a given particle, nullptr if there is no interaction. More...
 
virtual void actionsOnRestart ()
 No implementation but can be overidden in its derived classes. More...
 
virtual void actionsAfterParticleGhostUpdate ()
 No implementation but can be overidden in its derived classes. More...
 
virtual void handleParticleAddition (unsigned int id, BaseParticle *p)
 Handles the addition of particles to the particleHandler. More...
 
virtual void handleParticleRemoval (unsigned int id)
 Handles the addition of particles to the particleHandler. More...
 
virtual void checkInteractions (InteractionHandler *interactionHandler, unsigned int timeStamp)
 Check if all interactions are valid. More...
 
bool getVTKVisibility () const
 
void setVTKVisibility (bool vtkVisibility)
 
void addRenderedWall (BaseWall *w)
 
BaseWallgetRenderedWall (size_t i) const
 
void removeRenderedWalls ()
 
void renderWall (VTKContainer &vtk)
 
void addParticlesAtWall (unsigned numElements=50)
 
void setVelocityControl (Vec3D forceGoal, Vec3D gainFactor, Vec3D baseVelocity)
 
- Public Member Functions inherited from BaseInteractable
 BaseInteractable ()
 Default BaseInteractable constructor. More...
 
 BaseInteractable (const BaseInteractable &p)
 Copy constructor. More...
 
 ~BaseInteractable () override
 Destructor, it simply destructs the BaseInteractable and all the objects it contains. More...
 
unsigned int getIndSpecies () const
 Returns the index of the species associated with the interactable object. More...
 
const ParticleSpeciesgetSpecies () const
 Returns a pointer to the species of this BaseInteractable. More...
 
void setSpecies (const ParticleSpecies *species)
 Sets the species of this BaseInteractable. More...
 
const Vec3DgetForce () const
 Returns the force on this BaseInteractable. More...
 
const Vec3DgetTorque () const
 Returns the torque on this BaseInteractable. More...
 
void setForce (const Vec3D &force)
 Sets the force on this BaseInteractable. More...
 
void setTorque (const Vec3D &torque)
 Sets the torque on this BaseInteractable. More...
 
void addForce (const Vec3D &addForce)
 Adds an amount to the force on this BaseInteractable. More...
 
void addTorque (const Vec3D &addTorque)
 Adds an amount to the torque on this BaseInteractable. More...
 
void resetForceTorque (int numberOfOMPthreads)
 
void sumForceTorqueOMP ()
 
const Vec3DgetPosition () const
 Returns the position of this BaseInteractable. More...
 
const QuaterniongetOrientation () const
 Returns the orientation of this BaseInteractable. More...
 
void setPosition (const Vec3D &position)
 Sets the position of this BaseInteractable. More...
 
void setOrientationViaNormal (Vec3D normal)
 Sets the orientation of this BaseInteractable by defining the vector that results from the rotation of the (1,0,0) vector. More...
 
void setOrientationViaEuler (Vec3D eulerAngle)
 Sets the orientation of this BaseInteractable by defining the euler angles. More...
 
void setOrientation (const Quaternion &orientation)
 Sets the orientation of this BaseInteractable. More...
 
virtual void move (const Vec3D &move)
 Moves this BaseInteractable by adding an amount to the position. More...
 
virtual void rotate (const Vec3D &angularVelocityDt)
 Rotates this BaseInteractable. More...
 
const std::vector< BaseInteraction * > & getInteractions () const
 Returns a list of interactions which belong to this interactable. More...
 
void addInteraction (BaseInteraction *I)
 Adds an interaction to this BaseInteractable. More...
 
bool removeInteraction (BaseInteraction *I)
 Removes an interaction from this BaseInteractable. More...
 
void copyInteractionsForPeriodicParticles (const BaseInteractable &p)
 Copies interactions to this BaseInteractable whenever a periodic copy made. More...
 
void setVelocity (const Vec3D &velocity)
 set the velocity of the BaseInteractable. More...
 
void setAngularVelocity (const Vec3D &angularVelocity)
 set the angular velocity of the BaseInteractble. More...
 
void addVelocity (const Vec3D &velocity)
 adds an increment to the velocity. More...
 
void addAngularVelocity (const Vec3D &angularVelocity)
 add an increment to the angular velocity. More...
 
virtual const Vec3DgetVelocity () const
 Returns the velocity of this interactable. More...
 
virtual const Vec3DgetAngularVelocity () const
 Returns the angular velocity of this interactable. More...
 
void setPrescribedPosition (const std::function< Vec3D(double)> &prescribedPosition)
 Allows the position of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedPosition (double time)
 Computes the position from the user defined prescribed position function. More...
 
void setPrescribedVelocity (const std::function< Vec3D(double)> &prescribedVelocity)
 Allows the velocity of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedVelocity (double time)
 Computes the velocity from the user defined prescribed velocity function. More...
 
void setPrescribedOrientation (const std::function< Quaternion(double)> &prescribedOrientation)
 Allows the orientation of the infinite mass interactbale to be prescribed. More...
 
void applyPrescribedOrientation (double time)
 Computes the orientation from the user defined prescribed orientation function. More...
 
void setPrescribedAngularVelocity (const std::function< Vec3D(double)> &prescribedAngularVelocity)
 Allows the angular velocity of the infinite mass interactable to be prescribed. More...
 
void applyPrescribedAngularVelocity (double time)
 Computes the angular velocity from the user defined prescribed angular velocity. More...
 
virtual const Vec3D getVelocityAtContact (const Vec3D &contact) const
 Returns the velocity at the contact point, use by many force laws. More...
 
void integrateBeforeForceComputation (double time, double timeStep)
 This is part of integrate routine for objects with infinite mass. More...
 
void integrateAfterForceComputation (double time, double timeStep)
 This is part of the integration routine for objects with infinite mass. More...
 
virtual Mdouble getInvMass () const
 
virtual Mdouble getCurvature (const Vec3D &labFixedCoordinates) const
 
- Public Member Functions inherited from BaseObject
 BaseObject ()=default
 Default constructor. More...
 
 BaseObject (const BaseObject &p)=default
 Copy constructor, copies all the objects BaseObject contains. More...
 
virtual ~BaseObject ()=default
 virtual destructor More...
 
virtual void moveInHandler (unsigned int index)
 Except that it is virtual, it does the same thing as setIndex() does. More...
 
void setIndex (unsigned int index)
 Allows one to assign an index to an object in the handler/container. More...
 
void setId (unsigned long id)
 Assigns a unique identifier to each object in the handler (container) which remains constant even after the object is deleted from the container/handler. More...
 
unsigned int getIndex () const
 Returns the index of the object in the handler. More...
 
unsigned int getId () const
 Returns the unique identifier of any particular object. More...
 
void setGroupId (unsigned groupId)
 
unsigned getGroupId () const
 

Private Attributes

std::vector< BaseWall * > walls_
 

Additional Inherited Members

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

Detailed Description

Restriction of a wall to the intersection with another wall.

This is a class defining walls. It defines the interaction of regular walls and periodic walls with particles as defined in Particle Modifications:

Constructor & Destructor Documentation

◆ BasicIntersectionOfWalls() [1/2]

BasicIntersectionOfWalls::BasicIntersectionOfWalls ( )

Default constructor, the normal is infinitely long.

36 {
37  logger(DEBUG, "BasicIntersectionOfWalls::BasicIntersectionOfWalls ) finished");
38 }

References FATAL, and logger.

Referenced by copy().

◆ BasicIntersectionOfWalls() [2/2]

BasicIntersectionOfWalls::BasicIntersectionOfWalls ( const BasicIntersectionOfWalls b)

Copy constructor, copy the given wall.

Parameters
[in]wBasicIntersectionOfWalls that has to be copied.

First copy the attributes of the BaseWall, then copy the ones that are specific for the BasicIntersectionOfWalls.

46  : BaseWall(b)
47 {
48  for (auto& w : b.walls_)
49  {
50  walls_.push_back(w->copy());
51  }
52  logger(DEBUG, "BasicIntersectionOfWalls::BasicIntersectionOfWalls(const BasicIntersectionOfWalls &p) finished");
53 }

References FATAL, logger, and walls_.

◆ ~BasicIntersectionOfWalls()

BasicIntersectionOfWalls::~BasicIntersectionOfWalls ( )
override

Default destructor.

56 {
57  for (auto& w : walls_)
58  {
59  delete w;
60  }
61  logger(DEBUG, "BasicIntersectionOfWalls::~BasicIntersectionOfWalls finished");
62 }

References FATAL, logger, and walls_.

Member Function Documentation

◆ add()

void BasicIntersectionOfWalls::add ( BaseWall wall)

Defines a standard wall, given an outward normal vector s.t. normal*x=normal*point for all x of the wall.

Todo:
TW maybe the Restricted wall should be templated with the wall type such that we don't need to use new and delete.
91 {
92  walls_.push_back(wall.copy());
93  walls_.back()->setId(walls_.size());
94 }

References BaseWall::copy(), and walls_.

Referenced by UnionOfWalls::setupInitialConditions().

◆ copy()

BasicIntersectionOfWalls * BasicIntersectionOfWalls::copy ( ) const
overridevirtual

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

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

Implements BaseWall.

68 {
69  return new BasicIntersectionOfWalls(*this);
70 }

References BasicIntersectionOfWalls().

◆ getDistanceAndNormal()

bool BasicIntersectionOfWalls::getDistanceAndNormal ( const BaseParticle p,
Mdouble distance,
Vec3D normal 
) 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 for which the distance to the wall must be computed.
[out]distanceDistance between the particle and the wall.
[out]normalThe normal of this wall, will only be set if there is a collision.
Returns
A boolean value for whether or not there is a collision.

Return the smallest distance and normal

Implements BaseWall.

104 {
105  if (walls_.empty())
106  {
107  logger(DEBUG, "Empty BasicIntersectionOfWalls");
108  return false;
109  }
110 
111  distance = -constants::inf; //distance of the closest wall
112  Mdouble distance2 = -constants::inf; //distance of the second-closest wall
113  Mdouble distance3 = -constants::inf; //distance of the third-closest wall
114  Mdouble distanceCurrent; //distance of teh current wall
115  Vec3D normal2; //normal of the second-closest wall
116  Vec3D normal3; //normal of the third-closest wall
117  Vec3D normalCurrent; //normal of the current wall
118  unsigned int id = static_cast<unsigned int>(-1); //id of the closest wall
119  unsigned int id2 = static_cast<unsigned int>(-1); //id of the second-closest wall
120  unsigned int id3 = static_cast<unsigned int>(-1); //id of the third-closest wall
121  Mdouble wallInteractionRadius = p.getWallInteractionRadius(this);
122  Vec3D position = p.getPosition() - getPosition();
123  getOrientation().rotateBack(position);
124  SphericalParticle shifted;
125  shifted.setSpecies(p.getSpecies());
126  shifted.setPosition(position);
127  shifted.setRadius(p.getRadius());
128 
129  //The object has to touch each wall each wall (distanceCurrent) and keep the minimum distance (distance) and wall index (id)
130  for (unsigned int i = 0; i < walls_.size(); i++)
131  {
132  //immediately stop if one of teh walls is not in interaction distance
133  if (!walls_[i]->getDistanceAndNormal(shifted, distanceCurrent, normalCurrent))
134  return false;
135  // Find out which of the walls is interacting with the particle.
136  // Keep the minimum distance (distance) and wall index (id)
137  // and store up to two walls (id2, id3) and their distances (distance2, distance3),
138  // if the possible contact point is near the intersection between id and id2 (and id3)
139  if (distanceCurrent > distance)
140  {
141  if (distance > -wallInteractionRadius)
142  {
143  if (distance2 > -wallInteractionRadius)
144  {
145  distance3 = distance;
146  normal3 = normal;
147  id3 = id;
148  }
149  else
150  {
151  distance2 = distance;
152  normal2 = normal;
153  id2 = id;
154  }
155  }
156  distance = distanceCurrent;
157  normal = normalCurrent;
158  id = i;
159  }
160  else if (distanceCurrent > -wallInteractionRadius)
161  {
162  if (distance2 > -wallInteractionRadius)
163  {
164  distance3 = distanceCurrent;
165  normal3 = normalCurrent;
166  id3 = i;
167  }
168  else
169  {
170  distance2 = distanceCurrent;
171  normal2 = normalCurrent;
172  id2 = i;
173  }
174  }
175  }
176  //logger(INFO, "particle %: contact with % % %, n % % %", p.getId(), distance, distance2, distance3, normal, normal2, normal3);
177 
178 
179  //If we are here, the closest wall is id;
180  //if distance2>-P.Radius, it's a possible face or vertex contact
181  //if distance3>-P.Radius, it's a possible vertex contact
182  if (distance2 > -wallInteractionRadius)
183  {
184  //contact point on wall id
185  SphericalParticle contact;
186  contact.setSpecies(p.getSpecies());
187  contact.setPosition(position + distance * normal);
188  contact.setRadius(0);
189 
190  //If the distance of D to id2 is positive, the contact is with the intersection
191  bool contactPointOutsideID2 = !walls_[id2]->getDistanceAndNormal(contact, distanceCurrent, normalCurrent);
192 
193  if (distance3 > -wallInteractionRadius &&
194  !walls_[id3]->getDistanceAndNormal(contact, distanceCurrent, normalCurrent))
195  {
196  if (contactPointOutsideID2)
197  {
198  //possible contact is with intersection of id,id2,id3
199  //we know id2<id3
200  unsigned int index =
201  (id < id2) ? ((id3 - 2) * (id3 - 1) * id3 / 6 + (id2 - 1) * id2 / 2 + id) :
202  (id < id3) ? ((id3 - 2) * (id3 - 1) * id3 / 6 + (id - 1) * id / 2 + id2) :
203  ((id - 2) * (id - 1) * id / 6 + (id3 - 1) * id3 / 2 + id2);
204  //find vertex C
205  Matrix3D N(normal.X, normal.Y, normal.Z, normal2.X, normal2.Y, normal2.Z, normal3.X, normal3.Y,
206  normal3.Z);
207  Vec3D P(Vec3D::dot(position, normal) + distance, Vec3D::dot(position, normal2) + distance2,
208  Vec3D::dot(position, normal3) + distance3);
209  Vec3D vertex = N.ldivide(P);
210  normal = position - vertex;
211  distance = sqrt(normal.getLengthSquared());
212  if (distance >= wallInteractionRadius)
213  return false; //no contact
214  normal /= -distance;
215  //logger(INFO, "vertex contact with % % %: pos %, n %, d %", id, id2, id3, p.getPosition(), normal, distance);
216  getOrientation().rotate(normal);
217  return true;
218  }
219  else
220  {
221  contactPointOutsideID2 = true;
222  distance2 = distance3;
223  id2 = id3;
224  }
225  }
226 
227  if (contactPointOutsideID2)
228  {
229  //find contact point C
230  normal3 = Vec3D::cross(normal, normal2);
231  Matrix3D N(normal.X, normal.Y, normal.Z, normal2.X, normal2.Y, normal2.Z, normal3.X, normal3.Y,
232  normal3.Z);
233  Vec3D P(Vec3D::dot(position, normal) + distance, Vec3D::dot(position, normal2) + distance2,
234  Vec3D::dot(position, normal3));
235  Vec3D C = N.ldivide(P);
236  normal = position - C;
237  distance = sqrt(normal.getLengthSquared());
238  if (distance >= wallInteractionRadius)
239  return false; //no contact
240  normal /= -distance;
241  //logger(INFO, "edge contact with % %: normal %", id, id2, normal);
242  getOrientation().rotate(normal);
243  return true;
244  }
245  }
246  //logger(INFO, "face contact with %", id);
247  getOrientation().rotate(normal);
248  return true;
249 }

References Vec3D::cross(), Vec3D::dot(), FATAL, BaseInteractable::getOrientation(), BaseInteractable::getPosition(), BaseParticle::getRadius(), BaseInteractable::getSpecies(), BaseParticle::getWallInteractionRadius(), constants::i, constants::inf, Matrix3D::ldivide(), logger, Quaternion::rotate(), Quaternion::rotateBack(), BaseInteractable::setPosition(), BaseParticle::setRadius(), BaseParticle::setSpecies(), walls_, Vec3D::X, Vec3D::Y, and Vec3D::Z.

Referenced by getVTK().

◆ getName()

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

Returns the name of the object, in this case the string "BasicIntersectionOfWalls".

Returns
The string "BasicIntersectionOfWalls", which is the name of this class.

Implements BaseObject.

292 {
293  return "BasicIntersectionOfWalls";
294 }

◆ getNumberOfObjects()

unsigned long BasicIntersectionOfWalls::getNumberOfObjects ( )

Returns the number of objects.

\detail Suppose your simulation adds to an BasicIntersectionOfWalls after a certain time or condition is met. Checking the number of objects is useful for checking if this has happened yet, when restarting.

78 {
79  return walls_.size();
80 }

References walls_.

◆ getObject()

BaseWall * BasicIntersectionOfWalls::getObject ( unsigned  i)
351 {
352  logger.assert_always(walls_.size() > i, "Index % exceeds number of walls %", i, walls_.size());
353  return walls_[i];
354 }

References constants::i, logger, and walls_.

◆ getVTK()

void BasicIntersectionOfWalls::getVTK ( std::vector< Vec3D > &  points,
std::vector< std::vector< double >> &  triangleStrips 
)
Todo:
change getVTK to writeVTK
Parameters
points
triangleStrips
Todo:
this function could be improved; might not plot full wall
297 {
298 // for (auto w : walls_) {
299 // w->writeVTK (points, triangleStrips);
300 // }
301  //writes points and strips for all walls; points are added to the global point vector, but the strips are held back
302  std::vector<std::vector<double>> myTriangleStrips;
303  unsigned long n = points.size();
304  for (auto w : walls_)
305  {
306  w->getVTK(points, myTriangleStrips);
307  }
308  //add position of the BasicIntersectionOfWalls to the point
309  for (std::vector<Vec3D>::iterator p = points.begin() + n; p != points.end(); p++)
310  {
311  getOrientation().rotate(*p);
312  }
313  //create a vector which points are in the wall (actually, only the new points are necessary)
314  std::vector<bool> pointInWall;
315  pointInWall.reserve(points.size());
316  SphericalParticle particle;
317  particle.setSpecies(getSpecies());
318  particle.setRadius(1e-10); //points within that distance are declared part of the wall
319  Mdouble distance;
320  Vec3D normal;
321  for (auto p : points)
322  {
323  particle.setPosition(p);
324  pointInWall.push_back(getDistanceAndNormal(particle, distance, normal));
325  }
326  //now loop through myTriangleStrips to find the strip parts that are fully inside the wall
327  std::vector<double> strip;
328  for (auto t : myTriangleStrips)
329  {
330  for (unsigned i : t)
331  {
332  if (pointInWall[i] == true)
333  {
334  strip.push_back(i);
335  }
336  else
337  {
338  if (strip.size() > 2)
339  triangleStrips.push_back(strip);
340  strip.clear();
341  }
342  }
343  if (strip.size() > 2)
344  triangleStrips.push_back(strip);
345  strip.clear();
346  }
348 }

References getDistanceAndNormal(), BaseInteractable::getOrientation(), BaseInteractable::getSpecies(), constants::i, n, Quaternion::rotate(), BaseInteractable::setPosition(), BaseParticle::setRadius(), BaseParticle::setSpecies(), and walls_.

◆ oldRead()

void BasicIntersectionOfWalls::oldRead ( std::istream &  is)

Reads BasicIntersectionOfWalls from an old-style restart file.

◆ read()

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

Reads BasicIntersectionOfWalls from a restart file.

Parameters
[in]isThe input stream from which the BasicIntersectionOfWalls is read.

Reimplemented from BaseWall.

255 {
256  BaseWall::read(is);
257  std::string dummy;
258  unsigned size;
259  // read e.g. "numIntersectionOfWalls 2"
260  is >> dummy >> size;
261  for (unsigned i = 0; i < size; i++)
262  {
263  std::string type;
264  // read e.g. "IntersectionOfWalls"
265  is >> type;
266  BaseWall* wall = getHandler()->createObject(type);
267  wall->setHandler(getHandler());
268  wall->read(is);
269  walls_.push_back(wall);
270  walls_.back()->setId(walls_.size());
271  }
272 }

References WallHandler::createObject(), BaseWall::getHandler(), constants::i, BaseWall::read(), BaseWall::setHandler(), and walls_.

◆ write()

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

Writes the BasicIntersectionOfWalls to an output stream, usually a restart file.

Parameters
[in]osThe output stream the BasicIntersectionOfWalls is written to.

Reimplemented from BaseWall.

278 {
279  BaseWall::write(os);
280  os << " numIntersectionOfWalls " << walls_.size();
281  for (auto w : walls_)
282  {
283  os << " ";
284  w->write(os);
285  }
286 }

References walls_, and BaseWall::write().

Member Data Documentation

◆ walls_

std::vector<BaseWall*> BasicIntersectionOfWalls::walls_
private

Outward normal vector. This does not have to be a unit vector.

Referenced by add(), BasicIntersectionOfWalls(), getDistanceAndNormal(), getNumberOfObjects(), getObject(), getVTK(), read(), write(), and ~BasicIntersectionOfWalls().


The documentation for this class was generated from the following files:
BaseWall::BaseWall
BaseWall()
Default constructor.
Definition: BaseWall.cc:36
BaseInteractable::getSpecies
const ParticleSpecies * getSpecies() const
Returns a pointer to the species of this BaseInteractable.
Definition: BaseInteractable.h:108
BaseInteractable::setPosition
void setPosition(const Vec3D &position)
Sets the position of this BaseInteractable.
Definition: BaseInteractable.h:239
BasicIntersectionOfWalls::walls_
std::vector< BaseWall * > walls_
Definition: BasicIntersectionOfWalls.h:114
BaseWall
Basic class for walls.
Definition: BaseWall.h:48
logger
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
Vec3D::X
Mdouble X
the vector components
Definition: Vector.h:65
Vec3D::dot
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:76
BaseWall::write
void write(std::ostream &os) const override
Function that writes a BaseWall to an output stream, usually a restart file.
Definition: BaseWall.cc:102
BaseInteractable::getOrientation
const Quaternion & getOrientation() const
Returns the orientation of this BaseInteractable.
Definition: BaseInteractable.h:230
BaseParticle::setRadius
virtual void setRadius(Mdouble radius)
Sets the particle's radius_ (and adjusts the mass_ accordingly, based on the particle's species)
Definition: BaseParticle.cc:542
Vec3D
Definition: Vector.h:50
Mdouble
double Mdouble
Definition: GeneralDefine.h:34
BaseWall::copy
virtual BaseWall * copy() const =0
Pure virtual function that can be overwritten in inherited classes in order to copy a BaseWall.
Matrix3D
Implementation of a 3D matrix.
Definition: Matrix.h:38
BaseParticle::getRadius
Mdouble getRadius() const
Returns the particle's radius.
Definition: BaseParticle.h:348
BaseParticle::getWallInteractionRadius
Mdouble getWallInteractionRadius(const BaseWall *wall) const
returns the radius plus the interactionDistance
Definition: BaseParticle.h:386
BaseParticle::setSpecies
void setSpecies(const ParticleSpecies *species)
Definition: BaseParticle.cc:804
Log::FATAL
@ FATAL
SphericalParticle
A spherical particle is the most simple particle used in MercuryDPM.
Definition: SphericalParticle.h:37
BaseWall::setHandler
virtual void setHandler(WallHandler *handler)
A function which sets the WallHandler for this BaseWall.
Definition: BaseWall.cc:127
BasicIntersectionOfWalls::BasicIntersectionOfWalls
BasicIntersectionOfWalls()
Default constructor, the normal is infinitely long.
Definition: BasicIntersectionOfWalls.cc:35
Vec3D::Y
Mdouble Y
Definition: Vector.h:65
constants::i
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51
BasicIntersectionOfWalls::getDistanceAndNormal
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....
Definition: BasicIntersectionOfWalls.cc:103
Quaternion::rotateBack
void rotateBack(Vec3D &position) const
Definition: Quaternion.cc:592
n
const unsigned n
Definition: CG3DPackingUnitTest.cpp:32
Vec3D::cross
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:163
BaseWall::getHandler
WallHandler * getHandler() const
A function which returns the WallHandler that handles this BaseWall.
Definition: BaseWall.cc:136
Quaternion::rotate
void rotate(Vec3D &position) const
Definition: Quaternion.cc:563
WallHandler::createObject
static BaseWall * createObject(const std::string &type)
Create a new wall, with the type given as a string (required for restarting).
Definition: WallHandler.cc:123
BaseWall::read
void read(std::istream &is) override
Function that reads a BaseWall from an input stream, usually a restart file.
Definition: BaseWall.cc:80
Vec3D::Z
Mdouble Z
Definition: Vector.h:65
constants::inf
const Mdouble inf
Definition: GeneralDefine.h:44
BaseInteractable::getPosition
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:218