MercuryDPM  Alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
BaseWall Class Referenceabstract

Basic class for walls. More...

#include <BaseWall.h>

+ Inheritance diagram for BaseWall:

Public Member Functions

 BaseWall ()
 Default constructor. It makes an empty BaseWall. More...
 
 BaseWall (const BaseWall &w)
 Copy constructor. More...
 
virtual ~BaseWall ()
 Default destructor. More...
 
virtual BaseWallcopy () const =0
 Pure virtual function that copies a BaseWall. More...
 
void read (std::istream &is)
 Function that reads a BaseWall from an input stream, usually a restart file. More...
 
void write (std::ostream &os) const
 Function that writes a BaseWall to an output stream, usually a restart file. More...
 
virtual MERCURY_DEPRECATED void clear ()
 A function that removes all data from this BaseWall, so sets handler_ to nullptr. More...
 
virtual bool getDistanceAndNormal (const BaseParticle &P, Mdouble &distance, Vec3D &normal_return) const =0
 Pure virtual function that computes the distance of a BaseParticle to this wall and returns the normal of this wall if there is a collision. More...
 
virtual void setHandler (WallHandler *handler)
 A function which sets the WallHandler for this BaseWall. More...
 
WallHandlergetHandler () const
 A function which returns the WallHandler that handles this BaseWall. More...
 
void setIndSpecies (unsigned int indSpecies)
 Define the species of this wall using the index of the species in the SpeciesHandler in this DPMBase. More...
 
void setSpecies (const ParticleSpecies *species)
 Define the species of this wall. More...
 
bool getLinePlaneIntersect (Vec3D &intersect, const Vec3D &p0, const Vec3D &p1, const Vec3D &n, const Vec3D &p)
 
bool isInsideWallVTK (const Vec3D &point, const Vec3D &normal, const Vec3D &position) const
 
void projectOntoWallVTK (Vec3D &point0, const Vec3D &point1, const Vec3D &normal, const Vec3D &position) const
 
void intersectVTK (std::vector< Vec3D > &points, const Vec3D normal, const Vec3D position) const
 
virtual void writeVTK (VTKContainer &vtk) 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...
 
virtual void move (const Vec3D &move)
 Moves this BaseInteractable by adding an amount to the position. 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 std::string getName () const =0
 A purely virtual function. 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...
 

Static Public Member Functions

static void addToVTK (const std::vector< Vec3D > &points, VTKContainer &vtk)
 

Public Attributes

WallHandlerhandler_
 

Detailed Description

Basic class for walls.

Class from which all walls inherit. Please note the getVelocity can for some walls be dependent on which point on the wall is meant.

Definition at line 44 of file BaseWall.h.

Constructor & Destructor Documentation

BaseWall::BaseWall ( )

Default constructor. It makes an empty BaseWall.

Definition at line 31 of file BaseWall.cc.

References DEBUG, handler_, and logger.

32 {
33  handler_ = nullptr;
34  logger(DEBUG, "BaseWall::BaseWall() finished");
35 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
WallHandler * handler_
Definition: BaseWall.h:154
BaseWall::BaseWall ( const BaseWall w)

Copy constructor.

Parameters
[in]wWall that must be copied.

Definition at line 40 of file BaseWall.cc.

References DEBUG, handler_, and logger.

42 {
43  handler_ = w.handler_;
44  logger(DEBUG, "BaseWall::BaseWall(const BaseWall &p) finished");
45 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
WallHandler * handler_
Definition: BaseWall.h:154
BaseInteractable()
Default BaseInteractable constructor, it simply creates an empty BaseInteractable.
BaseWall::~BaseWall ( )
virtual

Default destructor.

Definition at line 47 of file BaseWall.cc.

References DEBUG, and logger.

48 {
49  logger(DEBUG, "BaseWall::~BaseWall() finished");
50 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")

Member Function Documentation

void BaseWall::addToVTK ( const std::vector< Vec3D > &  points,
VTKContainer vtk 
)
static

Takes the points provided and adds a triangle strip connecting these points to the vtk container

Definition at line 289 of file BaseWall.cc.

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

Referenced by InfiniteWall::writeVTK().

290 {
291  if (points.size()!=0) {
292  //all all values in myPoints to points
293  vtk.points.insert(vtk.points.end(), points.begin(), points.end());
294 
295  // create one cell object containing all indices of the added points (created a triangle strip connecting these points)
296  std::vector<double> cell;
297  cell.reserve(vtk.points.size()+1);
298  cell.push_back(vtk.points.size()-1);
299  for (unsigned i=vtk.points.size()-points.size(); i<vtk.points.size(); i++) {
300  cell.push_back(i);
301  }
302 
303  //add this triangle strip to the vtk file
304  vtk.triangleStrips.push_back(cell);
305  }
306 }
std::vector< std::vector< double > > triangleStrips
Definition: BaseWall.h:36
std::vector< Vec3D > points
Definition: BaseWall.h:35
void BaseWall::clear ( )
virtual

A function that removes all data from this BaseWall, so sets handler_ to nullptr.

Deprecated:
Please don't use any clear() anymore, it will be gone soon.
Todo:
TW Why is this function deprecated? How else do I reset the wall properties e.g. of an intersection of walls

Reimplemented in IntersectionOfWalls, InfiniteWallWithHole, and CylindricalWall.

Definition at line 52 of file BaseWall.cc.

References DEBUG, and logger.

53 {
54  logger(DEBUG, "BaseWall::clear(), this function shouldn't be called");
55 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
virtual BaseWall* BaseWall::copy ( ) const
pure virtual

Pure virtual function that copies a BaseWall.

Returns
A pointer to the new BaseWall.

Implemented in TriangulatedWall, AxisymmetricIntersectionOfWalls, IntersectionOfWalls, InfiniteWall, SphericalWall, Coil, RestrictedWall, Screw, InfiniteWallWithHole, and CylindricalWall.

Referenced by RestrictedWall::set().

virtual bool BaseWall::getDistanceAndNormal ( const BaseParticle P,
Mdouble distance,
Vec3D normal_return 
) const
pure virtual

Pure virtual function that computes the distance of a BaseParticle to this wall and returns the normal of this wall if there is a collision.

Beware, the distance and normal are output parameters, not return values!

Parameters
[in]PReference to the BaseParticle we want to compute the distance to the BaseWall of.
[out]distanceDistance of the BaseParticle to the BaseWall.
[out]normal_returnThe normal of the wall. Is only given if there is a collision.
Returns
A boolean which indicates if there is a collision between the BaseParticle and the wall.

Implemented in IntersectionOfWalls, InfiniteWall, TriangulatedWall, InfiniteWallWithHole, AxisymmetricIntersectionOfWalls, SphericalWall, RestrictedWall, Coil, CylindricalWall, and Screw.

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

WallHandler * BaseWall::getHandler ( ) const

A function which returns the WallHandler that handles this BaseWall.

Returns
A pointer to the WallHandler that manages this BaseWall.

Definition at line 85 of file BaseWall.cc.

References handler_.

Referenced by InfiniteWall::createVTK(), IntersectionOfWalls::IntersectionOfWalls(), setHandler(), setIndSpecies(), setSpecies(), and AxisymmetricIntersectionOfWalls::writeVTK().

86 {
87  return handler_;
88 }
WallHandler * handler_
Definition: BaseWall.h:154
std::vector< BaseInteraction * > BaseWall::getInteractionWith ( BaseParticle p,
Mdouble  timeStamp,
InteractionHandler interactionHandler 
)
virtual
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 BaseWall and the BaseParticle at the timeStamp.
Todo:
{DK: What is the contact point for interactions with walls}

Implements BaseInteractable.

Reimplemented in TriangulatedWall, and RestrictedWall.

Definition at line 266 of file BaseWall.cc.

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

267 {
268  Mdouble distance;
269  Vec3D normal;
270  std::vector<BaseInteraction*> interactions;
271  if (getDistanceAndNormal(*p,distance,normal))
272  {
273  BaseInteraction* c = interactionHandler->getInteraction(p, this, timeStamp);
274  c->setNormal(-normal);
275  c->setDistance(distance);
276  c->setOverlap(p->getRadius() - distance);
278  c->setContactPoint(p->getPosition()-(p->getRadius() - 0.5 * c->getOverlap()) * c->getNormal());
279  interactions.push_back(c);
280  }
281  return interactions;
282 }
void setNormal(Vec3D normal)
Sets the normal vector between the two interacting objects.
void setOverlap(Mdouble overlap)
Set the overlap between the two interacting object.
BaseInteraction * getInteraction(BaseInteractable *P, BaseInteractable *I, Mdouble timeStamp)
Returns the Interaction between the BaseInteractable's P and I.
double Mdouble
void setContactPoint(Vec3D contactPoint)
Set the location of the contact point between the two interacting objects.
void setDistance(Mdouble distance)
Sets the interaction distance between the two interacting objects.
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Stores information about interactions between two interactable objects; often particles but could be ...
virtual bool getDistanceAndNormal(const BaseParticle &P, Mdouble &distance, Vec3D &normal_return) const =0
Pure virtual function that computes the distance of a BaseParticle to this wall and returns the norma...
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
bool BaseWall::getLinePlaneIntersect ( Vec3D intersect,
const Vec3D p0,
const Vec3D p1,
const Vec3D n,
const Vec3D p 
)

Definition at line 128 of file BaseWall.cc.

References Vec3D::dot().

129 {
130  // t = (p-p0).n / (p1-p0).n
131  //first compute the denominator
132  Mdouble denominator = Vec3D::dot(p1-p0,n);
133  if (fabs(denominator)>=1e-10)
134  {
135  Mdouble t = Vec3D::dot(p-p0,n)/denominator;
136  if (t<1+1e-12&&t>-1e-12)
137  intersect = p0 + t * (p1-p0);
138  return true;
139  }
140  return false;
141 }
double Mdouble
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:167
void BaseWall::intersectVTK ( std::vector< Vec3D > &  points,
const Vec3D  normal,
const Vec3D  position 
) const

Intersects a of set VTK points representing a wall with a half-space defined by position and normal; This is used in createVTK to restrict the VTK points representing a wall to the inside of the domain.

Checks if a set of VTK points is inside a half-space defined by position and normal; all points outside the half-space are projected onto the half-space boundary. Thus, if the old set of points represented a wall object, the new set of points is represents the intersection of the wall with the half-space.

Definition at line 163 of file BaseWall.cc.

References DEBUG, isInsideWallVTK(), logger, and projectOntoWallVTK().

Referenced by InfiniteWall::createVTK(), AxisymmetricIntersectionOfWalls::writeVTK(), and IntersectionOfWalls::writeVTK().

164 {
165  // find first point in Wall
166  std::vector<Vec3D>::iterator firstIn=points.begin();
167  for (;firstIn!=points.end(); firstIn++) {
168  //stop if points[first] is in domain
169  if (isInsideWallVTK(*firstIn,normal,position)) {
170  break;
171  }
172  }
173 
174  //if all points are out of the wall
175  if (firstIn==points.end()) {
176  logger(DEBUG, "BaseWall::intersectVTK: all points out of wall");
177  return;
178  }
179 
180  // find first point out of the wall after firstIn
181  std::vector<Vec3D>::iterator firstOut=firstIn+1;
182  for (;firstOut!=points.end(); firstOut++) {
183  if (!isInsideWallVTK(*firstOut,normal,position)) {
184  break;
185  }
186  }
187 
188  //if all points are in the wall
189  if (firstOut==points.end()&&firstIn==points.begin()) {
190  logger(DEBUG, "BaseWall::intersectVTK: points completely in wall; removing points");
191  points.clear();
192  return;
193  }
194 
195  //if the sequence starts with a point out of the wall
196  //Several cases have to be distinguished, multiple points in the wall the last point in or out of the wall: ooiiioo, ooioo, ooiii, or ooi
197  //In addition, we add the case iiioo and ioo
198  if (firstIn!=points.begin() || !isInsideWallVTK(points.back(),normal,position)) {
199  // remove unnessesary points in the wall
200  // ooiiioo -> ooiioo, ooiii -> ooii, iiioo -> iioo
201  if (firstOut-firstIn>2) {
202  logger(DEBUG, "BaseWall::intersectVTK: remove unnessesary points in the wall");
203  //necessary reset of firstOut, as erase invalidates iterators after erase point
204  points.erase(firstIn+1,firstOut-1); //note: erase does not delete the last point
205  firstOut = firstIn+2;
206  }
207 
208  // if there is only one point in the wall, make it two
209  // ooioo -> ooiioo, ooi -> ooii
210  if (firstOut==firstIn+1) {
211  logger(DEBUG, "BaseWall::intersectVTK: there is only one point in the wall, make it two");
212  //necessary reset of firstIn, firstOut, as insert invalidates iterators
213  unsigned in = firstIn-points.begin();
214  points.insert(firstIn+1,*firstIn);
215  firstIn = points.begin()+in;//necessary, unless capacity is set right
216  firstOut = firstIn+2;
217  }
218 
219  // three cases remain: ooiioo, ooii, iioo
220 
221  //move both points onto the surface of the wall
222  if (firstIn!=points.begin()) {
223  logger(DEBUG, "BaseWall::intersectVTK: move first point onto the surface of the wall");
224  projectOntoWallVTK(*firstIn, *(firstIn-1), normal, position);
225  } else {
226  logger(DEBUG, "BaseWall::intersectVTK: move first point (at the beginning of the list) onto the surface of the wall");
227  projectOntoWallVTK(points.front(), points.back(), normal, position);
228  }
229 
230  if (firstOut!=points.end()) {
231  logger(DEBUG, "BaseWall::intersectVTK: move second point onto the surface of the wall");
232  projectOntoWallVTK(*(firstOut-1), *firstOut, normal, position);
233  } else {
234  logger(DEBUG, "BaseWall::intersectVTK: move second point (at the end of the list) onto the surface of the wall");
235  projectOntoWallVTK(points.back(), points.front(), normal, position);
236  }
237  //if sequence starts and ends with a point in the wall: iiiooiii
238  } else {
239  logger(DEBUG, "BaseWall::intersectVTK: sequence starts and ends with a point in the wall");
240 
241  // find first point in wall after firstOut
242  for (firstIn=firstOut+1;firstIn!=points.end(); firstIn++) {
243  if (isInsideWallVTK(*firstIn,normal,position)) {
244  break;
245  }
246  }
247 
248  // remove unnessesary points in the wall
249  // iiiooiii -> iooi
250  points.erase(firstIn+1,points.end());
251  points.erase(points.begin(),firstOut-1); //note: erase does not delete the last point //note iterators are invalid now
252 
253  //move both points onto the surface of the wall: iooi
254  projectOntoWallVTK(points.front(), *(points.begin()+1), normal, position);
255  projectOntoWallVTK(points.back(), *(points.end()-2), normal, position);
256  }
257 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
bool isInsideWallVTK(const Vec3D &point, const Vec3D &normal, const Vec3D &position) const
Definition: BaseWall.cc:144
void projectOntoWallVTK(Vec3D &point0, const Vec3D &point1, const Vec3D &normal, const Vec3D &position) const
Definition: BaseWall.cc:154
bool BaseWall::isInsideWallVTK ( const Vec3D point,
const Vec3D normal,
const Vec3D position 
) const

Checks if point is in wall (if close to the wall, the point is assumed out of the wall)

Definition at line 144 of file BaseWall.cc.

References Vec3D::dot().

Referenced by intersectVTK().

144  {
145  return Vec3D::dot(position-point,normal)<-1e-12;
146 }
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:167
void BaseWall::projectOntoWallVTK ( Vec3D point0,
const Vec3D point1,
const Vec3D normal,
const Vec3D position 
) const

Moves point0 onto the surface of the wall such that the direction of the edge from point0 to point1 is preserved.

intersectVTK = point[i] + t*(point[i-1]-point[i]) and (intersectVTK - position) normal_=0. => (position-point[i]-t(point[i-1]-point[i]))*normal_=0 t=(position-point[i]))*normal_ / (point[i-1]-point[i]))*normal_

Definition at line 154 of file BaseWall.cc.

References Vec3D::dot().

Referenced by intersectVTK().

154  {
155  Vec3D dPoint = point1-point0;
156  point0 += Vec3D::dot(position-point0,normal)/Vec3D::dot(dPoint,normal)*dPoint;
157 }
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:167
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:45
void BaseWall::read ( std::istream &  is)
virtual

Function that reads a BaseWall from an input stream, usually a restart file.

Parameters
[in]isInput stream from which the BaseWall is read.

Implements BaseInteractable.

Reimplemented in IntersectionOfWalls, TriangulatedWall, InfiniteWall, InfiniteWallWithHole, SphericalWall, Coil, RestrictedWall, Screw, and CylindricalWall.

Definition at line 60 of file BaseWall.cc.

References BaseInteractable::read().

Referenced by CylindricalWall::read(), Screw::read(), RestrictedWall::read(), Coil::read(), SphericalWall::read(), InfiniteWallWithHole::read(), InfiniteWall::read(), and IntersectionOfWalls::read().

61 {
63 }
virtual void read(std::istream &is)=0
Reads a BaseInteractable from an input stream.
void BaseWall::setHandler ( WallHandler handler)
virtual

A function which sets the WallHandler for this BaseWall.

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

Reimplemented in IntersectionOfWalls.

Definition at line 76 of file BaseWall.cc.

References getHandler(), BaseInteractable::getIndSpecies(), BaseHandler< T >::getObject(), handler_, and setSpecies().

Referenced by WallHandler::addObject(), IntersectionOfWalls::setHandler(), and setSpecies().

77 {
78  handler_ = handler;
79  setSpecies(getHandler()->getDPMBase()->speciesHandler.getObject(getIndSpecies()));
80 }
WallHandler * handler_
Definition: BaseWall.h:154
WallHandler * getHandler() const
A function which returns the WallHandler that handles this BaseWall.
Definition: BaseWall.cc:85
T * getObject(const unsigned int id)
Gets a pointer to the Object at the specified index in the BaseHandler.
Definition: BaseHandler.h:451
unsigned int getIndSpecies() const
Returns the index of the Species of this BaseInteractable.
void setSpecies(const ParticleSpecies *species)
Define the species of this wall.
Definition: BaseWall.cc:113
void BaseWall::setIndSpecies ( unsigned int  indSpecies)
virtual

Define the species of this wall using the index of the species in the SpeciesHandler in this DPMBase.

Deprecated:
TW: this function should be taken out and replaced by setSpecies
Parameters
[in]indSpeciesThe index of the species of this BaseWall in the SpeciesHandler.

Reimplemented from BaseInteractable.

Definition at line 93 of file BaseWall.cc.

References ERROR, getHandler(), BaseObject::getId(), BaseInteractable::getIndSpecies(), BaseHandler< T >::getObject(), handler_, logger, BaseInteractable::setIndSpecies(), and setSpecies().

94 {
95  if (handler_ != nullptr)
96  {
97  setSpecies(getHandler()->getDPMBase()->speciesHandler.getObject(getIndSpecies()));
98  }
99  else
100  {
102  logger(ERROR, "setIndSpecies called on a particle with no particle handler.\n"
103  "Therefore I can't request the given species from the species handler.\n"
104  " PartID = %", getId());
105  }
106 }
unsigned int getId() const
Returns the unique identifier of any particular object.
Definition: BaseObject.cc:116
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
WallHandler * handler_
Definition: BaseWall.h:154
WallHandler * getHandler() const
A function which returns the WallHandler that handles this BaseWall.
Definition: BaseWall.cc:85
T * getObject(const unsigned int id)
Gets a pointer to the Object at the specified index in the BaseHandler.
Definition: BaseHandler.h:451
virtual void setIndSpecies(unsigned int indSpecies)
Sets the index of the Species of this BaseInteractable.
unsigned int getIndSpecies() const
Returns the index of the Species of this BaseInteractable.
void setSpecies(const ParticleSpecies *species)
Define the species of this wall.
Definition: BaseWall.cc:113
void BaseWall::setSpecies ( const ParticleSpecies species)

Define the species of this wall.

Todo:
TW: this function should also check if the particle is the correct particle for the species type.
Parameters
[in]speciesThe species this BaseWall is made of.

In addition to the functionality of BaseInteractable::setSpecies, this function sets the pointer to the wallHandler, which is needed to retrieve species information.

Definition at line 113 of file BaseWall.cc.

References BaseHandler< T >::getDPMBase(), BaseSpecies::getHandler(), getHandler(), setHandler(), BaseInteractable::setSpecies(), and DPMBase::wallHandler.

Referenced by IntersectionOfWalls::addObject(), Chute::createBottom(), InfiniteWall::InfiniteWall(), helpers::loadingTest(), helpers::normalAndTangentialLoadingTest(), FileReader::read(), WallHandler::readObject(), WallHandler::readOldObject(), DPMBase::readParAndIniFiles(), RestrictedWall::set(), setHandler(), setIndSpecies(), IntersectionOfWalls::setSpecies(), ChuteBottom::setupInitialConditions(), Chute::setupSideWalls(), and TriangulatedWall::TriangulatedWall().

114 {
116 
117  //set pointer to the handler, which is needed to retrieve species information
118  if (getHandler() == nullptr)
119  {
120  SpeciesHandler* sH = species->getHandler();
121  DPMBase* dB = sH->getDPMBase();
122  if (dB != nullptr)
123  setHandler(&dB->wallHandler);
124  }
125 }
Container to store all ParticleSpecies.
The DPMBase header includes quite a few header files, defining all the handlers, which are essential...
Definition: DPMBase.h:65
virtual void setHandler(WallHandler *handler)
A function which sets the WallHandler for this BaseWall.
Definition: BaseWall.cc:76
SpeciesHandler * getHandler() const
Returns the pointer to the handler to which this species belongs.
Definition: BaseSpecies.cc:74
WallHandler * getHandler() const
A function which returns the WallHandler that handles this BaseWall.
Definition: BaseWall.cc:85
void setSpecies(const ParticleSpecies *species)
Sets the species of this BaseInteractable.
WallHandler wallHandler
An object of the class WallHandler. Contains pointers to all the walls created.
Definition: DPMBase.h:1006
DPMBase * getDPMBase()
Gets the problem that is solved using this handler.
Definition: BaseHandler.h:543
void BaseWall::write ( std::ostream &  os) const
virtual

Function that writes a BaseWall to an output stream, usually a restart file.

Parameters
[in]osOutput stream the BaseWall has to be written to.

Implements BaseInteractable.

Reimplemented in IntersectionOfWalls, InfiniteWall, TriangulatedWall, InfiniteWallWithHole, SphericalWall, Coil, RestrictedWall, Screw, and CylindricalWall.

Definition at line 68 of file BaseWall.cc.

References BaseInteractable::write().

Referenced by CylindricalWall::write(), Screw::write(), RestrictedWall::write(), Coil::write(), SphericalWall::write(), InfiniteWallWithHole::write(), InfiniteWall::write(), and IntersectionOfWalls::write().

69 {
71 }
virtual void write(std::ostream &os) const =0
Write a BaseInteractable to an output stream.
void BaseWall::writeVTK ( VTKContainer vtk) const
virtual

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 in IntersectionOfWalls, InfiniteWall, AxisymmetricIntersectionOfWalls, RestrictedWall, and TriangulatedWall.

Definition at line 284 of file BaseWall.cc.

References BaseObject::getIndex(), BaseObject::getName(), logger, and WARN.

Referenced by RestrictedWall::writeVTK().

285 {
286  logger(WARN,"Wall % (%) cannot has no vtk writer defined",getIndex(),getName());
287 }
unsigned int getIndex() const
Returns the index of the object in the handler.
Definition: BaseObject.cc:108
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
virtual std::string getName() const =0
A purely virtual function.

Member Data Documentation

WallHandler* BaseWall::handler_

A pointer to the WallHandler that handles this BaseWall.

Definition at line 154 of file BaseWall.h.

Referenced by BaseWall(), getHandler(), setHandler(), and setIndSpecies().


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