PeriodicBoundary Class Reference

Defines a pair of periodic walls. Inherits from BaseBoundary. More...

#include <PeriodicBoundary.h>

+ Inheritance diagram for PeriodicBoundary:

Public Member Functions

 PeriodicBoundary ()
 default constructor More...
 
 ~PeriodicBoundary () override
 destructor More...
 
PeriodicBoundarycopy () const override
 copy method More...
 
 PeriodicBoundary (const PeriodicBoundary &other)
 copy constructor More...
 
void set (Vec3D normal, Mdouble distanceLeft, Mdouble distanceRight)
 Defines a PeriodicBoundary by its normal and positions. More...
 
void set (Vec3D normal, Vec3D positionLeft, Vec3D positionRight)
 As above, but by specifying two positions that the boundaries go through instead of distanceLeft and distanceRight. More...
 
void set (Vec3D normal, Mdouble distanceLeft, Mdouble distanceRight, Vec3D shiftDirection)
 For general parallelogramic domains, the direction of the shift vector can to be set manually. More...
 
void setPlanewiseShift (Vec3D planewiseShift)
 Set the planewise shift (projected onto the planewise direction, and zero by default). More...
 
Vec3D getNormal () const
 returns the vector normal to the periodic boundary More...
 
Mdouble getDistanceLeft () const
 Returns the distance of the left wall to the origin, in normal direction. More...
 
Mdouble getDistanceRight () const
 Returns the distance of the right wall to the origin, in normal direction. More...
 
Vec3D getShift () const
 Returns the vector going from the left to the right side of the periodic boundary. More...
 
void moveLeft (Mdouble distanceLeft)
 Sets the distance from the origin of the 'left' periodic wall. More...
 
void moveRight (Mdouble distanceRight)
 Sets the distance from the origin of the 'right' periodic wall. More...
 
Mdouble getDistance (const BaseParticle &p) const override
 Returns the distance of the edge to the particle. More...
 
Mdouble getDistance (const Vec3D &position) const override
 Returns the distance of the edge to the position. More...
 
virtual void shiftPosition (BaseParticle *p) const override
 shifts the particle More...
 
void shiftPosition (Vec3D &p) const
 
virtual void shiftPositions (Vec3D &postition1, Vec3D &postion2) const
 shifts two positions More...
 
virtual bool isClosestToLeftBoundary (const BaseParticle &p) const
 Returns TRUE if particle checked is closest to the 'left' edge, and FALSE if it is closest to the 'right' edge. More...
 
bool isClosestToLeftBoundary (const Vec3D &p) const override
 
virtual void createPeriodicParticles (ParticleHandler &pH) override
 Checks distance of particle to closer edge and creates a periodic copy if necessary. More...
 
void createGhostParticle (BaseParticle *pReal)
 Creates and adds a ghost particle from a given real particle. More...
 
void createPeriodicParticle (BaseParticle *p, ParticleHandler &pH) override
 Creates a single periodic particle if required from a given particle. More...
 
virtual void checkBoundaryAfterParticlesMove (ParticleHandler &pH) override
 Loops over particles, checks if each particle has crossed either boundary edge, and applies a shift if that is the case. More...
 
virtual void read (std::istream &is) override
 reads boundary properties from istream More...
 
MERCURYDPM_DEPRECATED void oldRead (std::istream &is)
 deprecated version of CubeInsertionBoundary::read(). More...
 
void write (std::ostream &os) const override
 writes boundary properties to ostream More...
 
std::string getName () const override
 Returns the name of the object. More...
 
- Public Member Functions inherited from BasePeriodicBoundary
 BasePeriodicBoundary ()
 default constructor. More...
 
 BasePeriodicBoundary (const BasePeriodicBoundary &b)
 copy constructor More...
 
 ~BasePeriodicBoundary () override
 destructor More...
 
void setPeriodicHandler (PeriodicBoundaryHandler *periodicHandler)
 Sets the periodicBoundaryHandler, required for parallel periodic boundaries. More...
 
PeriodicBoundaryHandlergetPeriodicHandler () const
 Returns the periodic boundary handler. More...
 
void createPeriodicParticles (ParticleHandler &pH) override
 Creates periodic ocpies of given particle in case of periodic boundaries in serial build. More...
 
virtual void modifyPeriodicComplexity (std::vector< int > &complexity, int &totalPeriodicComplexity, BaseParticle *particle, int i) const
 Modifies periodic complexity of a particle if necessary (i.e. maser boundary) More...
 
virtual void performActionsBeforeAddingParticles ()
 Actions that need to be performed before adding new ghost particles. More...
 
- Public Member Functions inherited from BaseBoundary
 BaseBoundary ()
 default constructor. More...
 
 BaseBoundary (const BaseBoundary &b)
 copy constructor More...
 
 ~BaseBoundary () override
 destructor More...
 
virtual void createPeriodicParticle (BaseParticle *p UNUSED, ParticleHandler &pH UNUSED)
 Creates a periodic particle in case of periodic boundaries in serial build. More...
 
virtual void createPeriodicParticles (ParticleHandler &pH UNUSED)
 Creates periodic copies of given particle in case of periodic boundaries. More...
 
virtual void checkBoundaryBeforeTimeStep (DPMBase *md)
 Virtual function that does things before each time step. More...
 
virtual void actionsBeforeTimeLoop ()
 Virtual function that does something after DPMBase::setupInitialConditions but before the first time step. More...
 
virtual void modifyGhostAfterCreation (BaseParticle *particle, int i)
 
virtual void writeVTK (std::fstream &file)
 
void setHandler (BoundaryHandler *handler)
 Sets the boundary's BoundaryHandler. More...
 
BoundaryHandlergetHandler () const
 Returns the boundary's BoundaryHandler. More...
 
- 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
 

Protected Attributes

Vec3D normal_
 outward unit normal vector for right edge More...
 
Mdouble distanceLeft_
 position of left edge, s.t. normal*x = distanceLeft_ More...
 
Mdouble distanceRight_
 position of right edge, s.t. normal*x = distanceRight_ More...
 
Mdouble scaleFactor_
 This is the normal to rescale the normal vector to a unit vectors. More...
 
Vec3D shift_
 shift from left to right boundary More...
 

Detailed Description

Defines a pair of periodic walls. Inherits from BaseBoundary.

The particles are in {x: position_left<=normal*x <position_right}, with normal being the outward unit normal vector of the right wall. If a particle moves outside these boundaries, it will be shifted.

Constructor & Destructor Documentation

◆ PeriodicBoundary() [1/2]

PeriodicBoundary::PeriodicBoundary ( )

default constructor

constructor

39 {
40  distanceLeft_ = std::numeric_limits<double>::quiet_NaN();
41  distanceRight_ = std::numeric_limits<double>::quiet_NaN();
42  scaleFactor_ = std::numeric_limits<double>::quiet_NaN();
43 
44  logger(DEBUG, "PeriodicBoundary::PeriodicBoundary() finished");
45 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ DEBUG
BasePeriodicBoundary()
default constructor.
Definition: BasePeriodicBoundary.cc:35
Mdouble scaleFactor_
This is the normal to rescale the normal vector to a unit vectors.
Definition: PeriodicBoundary.h:211
Mdouble distanceLeft_
position of left edge, s.t. normal*x = distanceLeft_
Definition: PeriodicBoundary.h:201
Mdouble distanceRight_
position of right edge, s.t. normal*x = distanceRight_
Definition: PeriodicBoundary.h:206

References DEBUG, distanceLeft_, distanceRight_, logger, and scaleFactor_.

Referenced by copy().

◆ ~PeriodicBoundary()

PeriodicBoundary::~PeriodicBoundary ( )
override

destructor

destructor

51 {
52  logger(DEBUG, "PeriodicBoundary::~PeriodicBoundary() finished");
53 }

References DEBUG, and logger.

◆ PeriodicBoundary() [2/2]

PeriodicBoundary::PeriodicBoundary ( const PeriodicBoundary other)

copy constructor

Copy constructor

67  : BasePeriodicBoundary(other)
68 {
69  normal_ = other.normal_;
70  scaleFactor_ = other.scaleFactor_;
73  shift_ = other.shift_;
74 }
Vec3D normal_
outward unit normal vector for right edge
Definition: PeriodicBoundary.h:197
Vec3D shift_
shift from left to right boundary
Definition: PeriodicBoundary.h:216

References distanceLeft_, distanceRight_, normal_, scaleFactor_, and shift_.

Member Function Documentation

◆ checkBoundaryAfterParticlesMove()

void PeriodicBoundary::checkBoundaryAfterParticlesMove ( ParticleHandler pH)
overridevirtual

Loops over particles, checks if each particle has crossed either boundary edge, and applies a shift if that is the case.

Loops through all particles to see if they have become ghosts. If that is the case their position is shifted. Note: This is only for a serial build - periodic particles work different in paralle

Parameters
[in]

Reimplemented from BasePeriodicBoundary.

Reimplemented in SubcriticalMaserBoundaryTEST.

375 {
376 #ifdef MERCURYDPM_USE_MPI
377  if (NUMBER_OF_PROCESSORS == 1)
378  {
379 #endif
380  for (auto p = pH.begin(); p != pH.end(); ++p)
381  {
382  if (getDistance((*p)->getPosition()) < 0)
383  {
384  shiftPosition(*p);
386  }
387  }
388 #ifdef MERCURYDPM_USE_MPI
389  }
390 #endif
391 }
#define NUMBER_OF_PROCESSORS
For the MPI communication routines this quantity is often required. defining this macro makes the cod...
Definition: GeneralDefine.h:62
BoundaryHandler * getHandler() const
Returns the boundary's BoundaryHandler.
Definition: BaseBoundary.cc:143
const std::vector< T * >::const_iterator begin() const
Gets the begin of the const_iterator over all Object in this BaseHandler.
Definition: BaseHandler.h:690
const std::vector< T * >::const_iterator end() const
Gets the end of the const_iterator over all BaseBoundary in this BaseHandler.
Definition: BaseHandler.h:704
DPMBase * getDPMBase()
Gets the problem that is solved using this handler.
Definition: BaseHandler.h:725
virtual void hGridUpdateMove(BaseParticle *, Mdouble)
Definition: DPMBase.cc:1933
Mdouble getDistance(const BaseParticle &p) const override
Returns the distance of the edge to the particle.
Definition: PeriodicBoundary.cc:197
virtual void shiftPosition(BaseParticle *p) const override
shifts the particle
Definition: PeriodicBoundary.cc:219
static Mdouble getLengthSquared(const Vec3D &a)
Calculates the squared length of a Vec3D: .
Definition: Vector.h:332

References BaseHandler< T >::begin(), BaseHandler< T >::end(), getDistance(), BaseHandler< T >::getDPMBase(), BaseBoundary::getHandler(), Vec3D::getLengthSquared(), DPMBase::hGridUpdateMove(), NUMBER_OF_PROCESSORS, shift_, and shiftPosition().

◆ copy()

PeriodicBoundary * PeriodicBoundary::copy ( ) const
overridevirtual

copy method

Copy method; creates a copy on the heap and returns its pointer.

Implements BasePeriodicBoundary.

Reimplemented in SubcriticalMaserBoundaryTEST.

59 {
60  return new PeriodicBoundary(*this);
61 }
PeriodicBoundary()
default constructor
Definition: PeriodicBoundary.cc:37

References PeriodicBoundary().

◆ createGhostParticle()

void PeriodicBoundary::createGhostParticle ( BaseParticle pReal)

Creates and adds a ghost particle from a given real particle.

315 {
317 
318  //Step 1: Copy the particle to new ghost particle.
319  BaseParticle* pGhost = pReal->copy();
320 
321  //Step 2: Copy the interactions of the ghost particle.
322  pGhost->copyInteractionsForPeriodicParticles(*pReal);
323 
324  //Step 3: Shift the ghost to the 'reflected' location.
325  shiftPosition(pGhost);
326 
327  //Step 4: If Particle is double shifted, get correct original particle
328  BaseParticle* from = pReal;
329  while (from->getPeriodicFromParticle() != nullptr)
330  from = from->getPeriodicFromParticle();
331  pGhost->setPeriodicFromParticle(from);
332  pGhost->setPeriodicGhostParticle(true);
333 
334  pH.addObject(pGhost);
335 }
void copyInteractionsForPeriodicParticles(const BaseInteractable &p)
Copies interactions to this BaseInteractable whenever a periodic copy made.
Definition: BaseInteractable.cc:386
Definition: BaseParticle.h:54
void setPeriodicFromParticle(BaseParticle *p)
Assigns the pointer to the 'original' particle this one's a periodic copy of (used in periodic bounda...
Definition: BaseParticle.h:441
virtual BaseParticle * copy() const =0
Particle copy method. It calls to copy constructor of this Particle, useful for polymorphism.
void setPeriodicGhostParticle(bool flag)
Flags the status of the particle to be a ghost in periodic boundary or not.
Definition: BaseParticle.cc:302
BaseParticle * getPeriodicFromParticle() const
Returns the 'original' particle this one's a periodic copy of.
Definition: BaseParticle.h:341
ParticleHandler particleHandler
An object of the class ParticleHandler, contains the pointers to all the particles created.
Definition: DPMBase.h:1437
Container to store all BaseParticle.
Definition: ParticleHandler.h:48
void addObject(BaseParticle *P) override
Adds a BaseParticle to the ParticleHandler.
Definition: ParticleHandler.cc:171

References ParticleHandler::addObject(), BaseParticle::copy(), BaseInteractable::copyInteractionsForPeriodicParticles(), BaseHandler< T >::getDPMBase(), BaseBoundary::getHandler(), BaseParticle::getPeriodicFromParticle(), DPMBase::particleHandler, BaseParticle::setPeriodicFromParticle(), BaseParticle::setPeriodicGhostParticle(), and shiftPosition().

Referenced by createPeriodicParticle(), and SubcriticalMaserBoundaryTEST::createPeriodicParticle().

◆ createPeriodicParticle()

void PeriodicBoundary::createPeriodicParticle ( BaseParticle p,
ParticleHandler pH 
)
override

Creates a single periodic particle if required from a given particle.

Checks the distance of given particle to the closest of both periodic walls, and creates a periodic copy of the particle if needed (i.e. if the particle is closer to the periodic wall than the radius of the largest particle in the system).

Parameters
[in]pParticle to be checked and possibly periodically copied
[in,out]pHSystem's ParticleHandler, (1) from which the interaction radius of its largest particle is retrieved to determine the maximum distance from the wall at which a particle should still have a periodic copy created, and (2) to which a possible periodic copy of the particle will be added
305 {
306  //note that getDistance sets closestToLeftBoundary_ to true or false depending on which side is closest
308  if ((getDistance(*p) < maxDistance)&&(!p->isClump()))
309  {
311  }
312 }
double Mdouble
Definition: GeneralDefine.h:34
bool isClump() const
Checks if particle is a clump (container)
Definition: BaseParticle.h:664
Mdouble getMaxInteractionRadius() const
Returns the particle's interaction radius, which might be different from radius_ (e....
Definition: BaseParticle.h:362
BaseParticle * getLargestParticle() const
Returns the pointer of the largest particle in the particle handler. When mercury is running in paral...
Definition: ParticleHandler.cc:534
void createGhostParticle(BaseParticle *pReal)
Creates and adds a ghost particle from a given real particle.
Definition: PeriodicBoundary.cc:314

References createGhostParticle(), getDistance(), ParticleHandler::getLargestParticle(), BaseParticle::getMaxInteractionRadius(), and BaseParticle::isClump().

Referenced by createPeriodicParticles().

◆ createPeriodicParticles()

void PeriodicBoundary::createPeriodicParticles ( ParticleHandler pH)
overridevirtual

Checks distance of particle to closer edge and creates a periodic copy if necessary.

Checks the distance of given particle to the closest of both periodic walls, and creates a periodic copy of the particle if needed (i.e. if the particle is closer to the periodic wall than the radius of the largest particle in the system).

Parameters
[in]pParticle to be checked and possibly periodically copied
[in,out]pHSystem's ParticleHandler, (1) from which the interaction radius of its largest particle is retrieved to determine the maximum distance from the wall at which a particle should still have a periodic copy created, and (2) to which a possible periodic copy of the particle will be added
350 {
351 #ifdef MERCURYDPM_USE_MPI
352  if (NUMBER_OF_PROCESSORS == 1)
353  {
354 #endif
355  unsigned numberOfParticles = pH.getSize();
356 
357  for (unsigned i = 0; i < numberOfParticles; i++)
358  {
360  }
361 #ifdef MERCURYDPM_USE_MPI
362  }
363 #endif
364 }
unsigned int getSize() const
Gets the size of the particleHandler (including mpi and periodic particles)
Definition: BaseHandler.h:655
T * getObject(const unsigned int id)
Gets a pointer to the Object at the specified index in the BaseHandler.
Definition: BaseHandler.h:613
void createPeriodicParticle(BaseParticle *p, ParticleHandler &pH) override
Creates a single periodic particle if required from a given particle.
Definition: PeriodicBoundary.cc:304
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51

References createPeriodicParticle(), BaseHandler< T >::getObject(), BaseHandler< T >::getSize(), constants::i, and NUMBER_OF_PROCESSORS.

◆ getDistance() [1/2]

Mdouble PeriodicBoundary::getDistance ( const BaseParticle p) const
overridevirtual

Returns the distance of the edge to the particle.

Returns the distance to the closest edge of the boundary to the particle. Since this function should be called before calculating any Particle-Wall interactions, it can also be used to set the shift vector in case of curved walls. Positive means that the particle is insid-> the periodic domain, negative means that it is outsid-> the periodic domain.

Parameters
[in]pA reference to the particle which distance to the periodic boundary is calculated

Implements BasePeriodicBoundary.

198 {
199  return getDistance(p.getPosition());
200 }
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:218

References BaseInteractable::getPosition().

Referenced by ChuteWithPeriodicInflow::Check_and_Duplicate_Periodic_Particle(), checkBoundaryAfterParticlesMove(), DPMBase::checkParticleForInteractionLocalPeriodic(), createPeriodicParticle(), SubcriticalMaserBoundaryTEST::getDistance(), and ChuteWithPeriodicInflow::integrateBeforeForceComputation().

◆ getDistance() [2/2]

Mdouble PeriodicBoundary::getDistance ( const Vec3D position) const
overridevirtual

Returns the distance of the edge to the position.

Returns the distance to the edge closest to the position

Parameters
[in]positionA reference to the position which distance to the periodic boundary is to be calculated

Implements BasePeriodicBoundary.

Reimplemented in SubcriticalMaserBoundaryTEST.

208 {
209  Mdouble distanceFromPlaneThroughOrigin = Vec3D::dot(position, normal_);
210  return std::min(distanceFromPlaneThroughOrigin - distanceLeft_,
211  distanceRight_ - distanceFromPlaneThroughOrigin);
212 }
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:76

References distanceLeft_, distanceRight_, Vec3D::dot(), and normal_.

◆ getDistanceLeft()

Mdouble PeriodicBoundary::getDistanceLeft ( ) const

Returns the distance of the left wall to the origin, in normal direction.

Returns
The distance of the left wall to the origin, in normal direction
144 {
145  return distanceLeft_;
146 }

References distanceLeft_.

Referenced by ConstantMassFlowMaserBoundary::ConstantMassFlowMaserBoundary(), and SubcriticalMaserBoundary::SubcriticalMaserBoundary().

◆ getDistanceRight()

Mdouble PeriodicBoundary::getDistanceRight ( ) const

Returns the distance of the right wall to the origin, in normal direction.

Returns
The distance of the left wall to the origin, in normal direction
152 {
153  return distanceRight_;
154 }

References distanceRight_.

Referenced by ConstantMassFlowMaserBoundary::ConstantMassFlowMaserBoundary(), and SubcriticalMaserBoundary::SubcriticalMaserBoundary().

◆ getName()

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

Returns the name of the object.

Returns the name of the object class

Returns
the object's class' name, i.e. 'CubeInsertionBoundary'

Implements BaseObject.

Reimplemented in SubcriticalMaserBoundaryTEST.

441 {
442  return "PeriodicBoundary";
443 }

◆ getNormal()

Vec3D PeriodicBoundary::getNormal ( ) const

returns the vector normal to the periodic boundary

Returns
The vector perpendicular to the periodic boundary
136 {
137  return normal_;
138 }

References normal_.

Referenced by ConstantMassFlowMaserBoundary::ConstantMassFlowMaserBoundary(), and SubcriticalMaserBoundary::SubcriticalMaserBoundary().

◆ getShift()

Vec3D PeriodicBoundary::getShift ( ) const

Returns the vector going from the left to the right side of the periodic boundary.

Returns
The vector going from the left to the right sid-> of the periodic boundary
160 {
161  return shift_;
162 }

References shift_.

Referenced by ConstantMassFlowMaserBoundary::ConstantMassFlowMaserBoundary(), and SubcriticalMaserBoundary::SubcriticalMaserBoundary().

◆ isClosestToLeftBoundary() [1/2]

bool PeriodicBoundary::isClosestToLeftBoundary ( const BaseParticle p) const
virtual

Returns TRUE if particle checked is closest to the 'left' edge, and FALSE if it is closest to the 'right' edge.

276 {
278 }
virtual bool isClosestToLeftBoundary(const BaseParticle &p) const
Returns TRUE if particle checked is closest to the 'left' edge, and FALSE if it is closest to the 'ri...
Definition: PeriodicBoundary.cc:275

References BaseInteractable::getPosition().

Referenced by ChuteWithPeriodicInflow::integrateBeforeForceComputation(), shiftPosition(), and shiftPositions().

◆ isClosestToLeftBoundary() [2/2]

bool PeriodicBoundary::isClosestToLeftBoundary ( const Vec3D p) const
overridevirtual

Returns TRUE if position checked is closest to the 'left' edge, and FALSE if it is closest to the 'right' edge`.

Implements BasePeriodicBoundary.

287 {
288  const Mdouble distance = Vec3D::dot(p, normal_);
289  return (distanceRight_ - distance > distance - distanceLeft_);
290 }

References distanceLeft_, distanceRight_, Vec3D::dot(), and normal_.

◆ moveLeft()

void PeriodicBoundary::moveLeft ( Mdouble  distanceLeft)

Sets the distance from the origin of the 'left' periodic wall.

Allows the left periodic boundary to be moved to a new position and automatically changes its shift value

Parameters
[in]distanceLeftThe distance (from the origin) to which the left boundary is moved
171 {
172  distanceLeft_ = distanceLeft * scaleFactor_;
174 }

References distanceLeft_, distanceRight_, normal_, scaleFactor_, and shift_.

◆ moveRight()

void PeriodicBoundary::moveRight ( Mdouble  distanceRight)

Sets the distance from the origin of the 'right' periodic wall.

Allows the right periodic wall to be moved to a new position and automatically changes its shift value

Parameters
[in]distanceRightThe distance (from the origin) to which the right boundary is moved
183 {
184  distanceRight_ = distanceRight * scaleFactor_;
186 }

References distanceLeft_, distanceRight_, normal_, scaleFactor_, and shift_.

◆ oldRead()

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

deprecated version of CubeInsertionBoundary::read().

Deprecated version of read().

Deprecated:
Should be gone by Mercury 2.0. Instead, use CubeInsertionBoundary::read().
413 {
414  std::string dummy;
415  is >> dummy >> normal_
416  >> dummy >> scaleFactor_ // JMFT: I'd like to deprecate scaleFactor_
417  >> dummy >> distanceLeft_
418  >> dummy >> distanceRight_
419  >> dummy >> shift_;
420 }

References distanceLeft_, distanceRight_, normal_, scaleFactor_, and shift_.

◆ read()

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

reads boundary properties from istream

Reads the boundary properties from an istream

Parameters
[in]isthe istream

Reimplemented from BasePeriodicBoundary.

Reimplemented in SubcriticalMaserBoundaryTEST.

398 {
400  std::string dummy;
401  is >> dummy >> normal_
402  >> dummy >> scaleFactor_ // JMFT: I'd like to deprecate scaleFactor_
403  >> dummy >> distanceLeft_
404  >> dummy >> distanceRight_
405  >> dummy >> shift_;
406 }
void read(std::istream &is) override
Reads the object's id_ from given istream.
Definition: BasePeriodicBoundary.cc:71

References distanceLeft_, distanceRight_, normal_, BasePeriodicBoundary::read(), scaleFactor_, and shift_.

Referenced by SubcriticalMaserBoundaryTEST::read().

◆ set() [1/3]

void PeriodicBoundary::set ( Vec3D  normal,
Mdouble  distanceLeft,
Mdouble  distanceRight 
)

Defines a PeriodicBoundary by its normal and positions.

Defines the boundary, given a normal vector such that all particles are within {x: position_left<=normal*x<position_right}. The shift vector is set assuming that the domain is rectangular (shift parallel to normal).

Parameters
[in]normalThe normal vector pointing from the left wall into the domain
[in]distanceLeftThe (signed) distance between the left wall and the origin
[in]distanceRightThe (signed) distance between the right wall and the origin
85 {
86  // factor is used to set normal to unit length
87  scaleFactor_ = 1. / std::sqrt(Vec3D::dot(normal, normal));
88  normal_ = normal * scaleFactor_;
89  distanceLeft_ = distanceLeft * scaleFactor_;
90  distanceRight_ = distanceRight * scaleFactor_;
91  logger.assert_always(distanceRight_ > distanceLeft_,
92  "PeriodicBoundary::set: left distance needs to be smaller than right distance");
94 }

References distanceLeft_, distanceRight_, Vec3D::dot(), logger, normal_, scaleFactor_, and shift_.

Referenced by Chutebelt::actionsAfterTimeStep(), AngledPeriodicBoundaryUnitTest::AngledPeriodicBoundaryUnitTest(), ChutePeriodicDemo::ChutePeriodicDemo(), ChuteWithContraction::ChuteWithContraction(), ChuteWithPeriodicInflowAndContraction::ChuteWithPeriodicInflowAndContraction(), ClosedCSCWalls::ClosedCSCWalls(), ContractionWithPeriodicInflow::ContractionWithPeriodicInflow(), Funnel::create_walls(), SegregationPeriodic::createWalls(), CSCWalls::CSCWalls(), ChuteWithPeriodicInflow::ExtendInWidth(), FluxAndPeriodicBoundarySelfTest::FluxAndPeriodicBoundarySelfTest(), InitialConditions< SpeciesType >::InitialConditions(), main(), BoundaryHandler::readOldObject(), DPMBase::readParAndIniFiles(), StressStrainControlBoundary::set(), set(), ChutePeriodic::setupInitialConditions(), MercuryLogo::setupInitialConditions(), my_problem::setupInitialConditions(), RotatingDrumWet::setupInitialConditions(), MaserRepeatedOutInMPI2Test::setupInitialConditions(), MpiPeriodicBoundaryUnitTest::setupInitialConditions(), PeriodicBounaryEnteringMPIDomainTest::setupInitialConditions(), CubicCell::setupInitialConditions(), FreeCooling2DinWalls::setupInitialConditions(), FreeCooling3DDemoProblem::setupInitialConditions(), FreeCoolingDemoProblem::setupInitialConditions(), MinimalExampleDrum::setupInitialConditions(), Cstatic2d::setupInitialConditions(), SilbertPeriodic::setupInitialConditions(), SegregationPeriodic::setupInitialConditions(), Chutebelt::setupInitialConditions(), TriangulatedScrewSelfTest::setupInitialConditions(), ScalingTestInitialConditionsRelax::setupInitialConditions(), Tutorial6::setupInitialConditions(), Packing::setupInitialConditions(), FullRestartTest::setupInitialConditions(), MD_demo::setupInitialConditions(), MpiMaserChuteTest::setupInitialConditions(), PeriodicWalls::setupInitialConditions(), PeriodicWallsWithSlidingFrictionUnitTest::setupInitialConditions(), ChuteBottom::setupInitialConditions(), and Chute::setupSideWalls().

◆ set() [2/3]

void PeriodicBoundary::set ( Vec3D  normal,
Mdouble  distanceLeft,
Mdouble  distanceRight,
Vec3D  shiftDirection 
)

For general parallelogramic domains, the direction of the shift vector can to be set manually.

like PeriodicBoundary::set(normal, distanceLeft, distanceRight), but including the possibility of setting the shift direction vector.

Parameters
[in]normalThe normal vector pointing from the left wall into the domain
[in]distanceLeftThe (signed) distance between the left wall and the origin
[in]distanceRightThe (signed) distance between the right wall and the origin
[in]shiftDirectionThe vector over which particles will be shifted when moving through the PeriodicBoundary
111 {
112  // factor is used to set normal to unit length
113  scaleFactor_ = 1. / std::sqrt(Vec3D::dot(normal, normal));
114  normal_ = normal * scaleFactor_;
115  distanceLeft_ = distanceLeft * scaleFactor_;
116  distanceRight_ = distanceRight * scaleFactor_;
117  // factor is used to set shift vector to correct length
119  shift_ = shiftDirection * scaleFactor_;
120 }

References distanceLeft_, distanceRight_, Vec3D::dot(), normal_, scaleFactor_, and shift_.

◆ set() [3/3]

void PeriodicBoundary::set ( Vec3D  normal,
Vec3D  positionLeft,
Vec3D  positionRight 
)

As above, but by specifying two positions that the boundaries go through instead of distanceLeft and distanceRight.

97 {
98  set(normal, Vec3D::dot(positionLeft, normal), Vec3D::dot(positionRight, normal));
99 }
void set(Vec3D normal, Mdouble distanceLeft, Mdouble distanceRight)
Defines a PeriodicBoundary by its normal and positions.
Definition: PeriodicBoundary.cc:84

References Vec3D::dot(), and set().

◆ setPlanewiseShift()

void PeriodicBoundary::setPlanewiseShift ( Vec3D  planewiseShift)

Set the planewise shift (projected onto the planewise direction, and zero by default).

Sets the shift_ vector through setting the planewise shift. We delete the component of planewiseShift that is parallel to normal_.

127 {
128  planewiseShift -= Vec3D::dot(planewiseShift, normal_) / Vec3D::dot(normal_, normal_) * normal_;
129  shift_ = normal_ * (distanceRight_ - distanceLeft_) + planewiseShift;
130 }

References distanceLeft_, distanceRight_, Vec3D::dot(), normal_, and shift_.

◆ shiftPosition() [1/2]

void PeriodicBoundary::shiftPosition ( BaseParticle p) const
overridevirtual

shifts the particle

Parameters
[in]pA pointer to the particle which will be shifted.

Shifts the particle either to the left or right, using the method isClosestToLeftBoundary to determine which sid-> it should be shifted to.

Parameters
[in]pA pointer to the particle which will be shifted.

Implements BasePeriodicBoundary.

220 {
221  if (isClosestToLeftBoundary(*p))
222  {
223  p->move(shift_);
224  }
225  else
226  {
227  p->move(-shift_);
228  }
229 }
virtual void move(const Vec3D &move)
Moves this BaseInteractable by adding an amount to the position.
Definition: BaseInteractable.cc:215

References isClosestToLeftBoundary(), BaseInteractable::move(), and shift_.

Referenced by ChuteWithPeriodicInflow::Check_and_Duplicate_Periodic_Particle(), SubcriticalMaserBoundaryTEST::checkBoundaryAfterParticleMoved(), checkBoundaryAfterParticlesMove(), Mercury3Dclump::checkClumpForInteractionPeriodic(), DPMBase::checkParticleForInteractionLocalPeriodic(), createGhostParticle(), and ChuteWithPeriodicInflow::integrateBeforeForceComputation().

◆ shiftPosition() [2/2]

void PeriodicBoundary::shiftPosition ( Vec3D p) const

Shifts the particle either to the left or right, using the method isClosestToLeftBoundary to determine which sid-> it should be shifted to.

Parameters
[in]pA pointer to the particle which will be shifted.
237 {
239  {
240  p += shift_;
241  }
242  else
243  {
244  p -= shift_;
245  }
246 }

References isClosestToLeftBoundary(), and shift_.

◆ shiftPositions()

void PeriodicBoundary::shiftPositions ( Vec3D position1,
Vec3D position2 
) const
virtual

shifts two positions

Shifts two given positions by the shift_ vector.

Parameters
[in]position1The first position to be shifted
[in]position2The second position to be shifted
Todo:
(AT) see toDo of PeriodicBoundary::shiftPosition().
256 {
257  if (isClosestToLeftBoundary(position1))
258  {
259  position1 += shift_;
260  position2 += shift_;
261  }
262  else
263  {
264  position1 -= shift_;
265  position2 -= shift_;
266  }
267 }

References isClosestToLeftBoundary(), and shift_.

◆ write()

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

writes boundary properties to ostream

Writes boundary's properties to an ostream

Parameters
[in]osthe ostream

Reimplemented from BasePeriodicBoundary.

Reimplemented in SubcriticalMaserBoundaryTEST.

427 {
429  os << " normal " << normal_
430  << " scaleFactor " << scaleFactor_ // JMFT: I'd like to deprecate scaleFactor_
431  << " distanceLeft " << distanceLeft_
432  << " distanceRight " << distanceRight_
433  << " shift " << shift_;
434 }
void write(std::ostream &os) const override
Adds object's id_ to given ostream.
Definition: BasePeriodicBoundary.cc:80

References distanceLeft_, distanceRight_, normal_, scaleFactor_, shift_, and BasePeriodicBoundary::write().

Referenced by SubcriticalMaserBoundaryTEST::write().

Member Data Documentation

◆ distanceLeft_

Mdouble PeriodicBoundary::distanceLeft_
protected

◆ distanceRight_

◆ normal_

◆ scaleFactor_

Mdouble PeriodicBoundary::scaleFactor_
protected

This is the normal to rescale the normal vector to a unit vectors.

Referenced by moveLeft(), moveRight(), oldRead(), PeriodicBoundary(), read(), set(), and write().

◆ shift_


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