ChuteWithPeriodicInflow Class Reference

Particles of a single Species. More...

#include <ChuteWithPeriodicInflow.h>

+ Inheritance diagram for ChuteWithPeriodicInflow:

Public Member Functions

 ChuteWithPeriodicInflow (std::string restart_file)
 
 ChuteWithPeriodicInflow (std::string restart_file, int numRepetitions)
 
 ChuteWithPeriodicInflow (std::string restart_file, int numRepetitions, int numRepetitionsInWidth)
 
double getInfo (const BaseParticle &P) const override
 A virtual function that returns some user-specified information about a particle. More...
 
void loadPeriodicBox (std::string const restart_file)
 loads periodic chute data from restart file More...
 
void AddContinuingBottom (int numRepetitions)
 
void ExtendInWidth (int numRepetitionsInWidth)
 
void actionsBeforeTimeStep () override
 Do not add, only remove particles. More...
 
void cleanChute ()
 Remove particles if they fall below a certain height (allows them to become supercritical) More...
 
void setupInitialConditions () override
 Do not add bottom. More...
 
void integrateBeforeForceComputation ()
 Update particles' and walls' positions and velocities before force computation. More...
 
void printTime () const override
 add some particular output More...
 
void computeInternalForces (BaseParticle *P1, BaseParticle *P2)
 
int Check_and_Duplicate_Periodic_Particle (BaseParticle *i, int nWallPeriodic)
 
void writeXBallsScript () const
 This writes a script which can be used to load the xballs problem to display the data just generated. More...
 
void outputXBallsDataParticlee (const unsigned int i, const unsigned int format, std::ostream &os) const
 
double get_PeriodicBoxLength ()
 
void set_PeriodicBoxLength (double new_)
 
int get_PeriodicBoxNSpecies ()
 
void set_PeriodicBoxNSpecies (int new_)
 
bool InPeriodicBox (BaseParticle *P)
 
- Public Member Functions inherited from Chute
 Chute ()
 This is the default constructor. All it does is set sensible defaults. More...
 
 Chute (const DPMBase &other)
 Copy constructor, converts an existing DPMBase problem into a Chute problem. More...
 
 Chute (const MercuryBase &other)
 Copy constructor, converts an existing MercuryBase problem into a Chute problem. More...
 
 Chute (const Mercury3D &other)
 Copy constructor, converts an existing Mercury3D problem into a Chute problem. More...
 
 Chute (const Chute &other)
 Default copy constructor. More...
 
void constructor ()
 This is the actual constructor METHOD; it is called by all constructors above (except the default copy constructor). More...
 
bool readNextArgument (int &i, int argc, char *argv[]) override
 This method can be used for reading object properties from a string. More...
 
void setupSideWalls ()
 Creates chute side walls (either solid or periodic) More...
 
void makeChutePeriodic ()
 This makes the chute periodic in Y. More...
 
bool getIsPeriodic () const
 Returns whether the chute is periodic in Y. More...
 
void setupInitialConditions () override
 Creates bottom, side walls and a particle insertion boundary. More...
 
void read (std::istream &is, ReadOptions opt=ReadOptions::ReadAll) override
 Reads all chute properties from an istream. More...
 
void write (std::ostream &os, bool writeAllParticles=true) const override
 This function writes the Chute properties to an ostream, and adds the properties of ALL chute particles as well. More...
 
void setFixedParticleRadius (Mdouble fixedParticleRadius)
 Sets the particle radius of the fixed particles which constitute the (rough) chute bottom. More...
 
Mdouble getFixedParticleRadius () const
 Returns the particle radius of the fixed particles which constitute the (rough) chute bottom. More...
 
void setFixedParticleSpacing (Mdouble fixedParticleSpacing)
 Sets the spacing of the fixed particles which constitute the (rough) chute bottom; used in triangular packing only. More...
 
Mdouble getFixedParticleSpacing () const
 Returns the particle radius of the fixed particles which constitute the (rough) chute bottom; used in triangular packing only. More...
 
void setRoughBottomType (RoughBottomType roughBottomType)
 Sets the type of rough bottom of the chute. More...
 
void setRoughBottomType (std::string roughBottomTypeString)
 Sets the type of rough bottom of the chute, using a string with the EXACT enum type as input. More...
 
RoughBottomType getRoughBottomType () const
 Returns the type of (rough) bottom of the chute. More...
 
void setChuteAngle (Mdouble chuteAngle)
 Sets gravity vector according to chute angle (in degrees) More...
 
void setChuteAngleAndMagnitudeOfGravity (Mdouble chuteAngle, Mdouble gravity)
 Sets gravity vector according to chute angle (in degrees) More...
 
Mdouble getChuteAngle () const
 Returns the chute angle (in radians) More...
 
Mdouble getChuteAngleDegrees () const
 Returns the chute angle (in degrees) More...
 
void setMaxFailed (unsigned int maxFailed)
 Sets the number of times a particle will be tried to be added to the insertion boundary. More...
 
unsigned int getMaxFailed () const
 Returns the number of times a particle will be tried to be added to the insertion boundary. More...
 
void setInflowParticleRadius (Mdouble inflowParticleRadius)
 Sets the radius of the inflow particles to a single one (i.e. ensures a monodisperse inflow). More...
 
void setInflowParticleRadius (Mdouble minInflowParticleRadius, Mdouble maxInflowParticleRadius)
 Sets the minimum and maximum radius of the inflow particles. More...
 
void setMinInflowParticleRadius (Mdouble minInflowParticleRadius)
 sets the minimum radius of inflow particles More...
 
void setMaxInflowParticleRadius (Mdouble maxInflowParticleRadius)
 Sets the maximum radius of inflow particles. More...
 
Mdouble getInflowParticleRadius () const
 Returns the average radius of inflow particles. More...
 
Mdouble getMinInflowParticleRadius () const
 returns the minimum radius of inflow particles More...
 
Mdouble getMaxInflowParticleRadius () const
 Returns the maximum radius of inflow particles. More...
 
void setInflowHeight (Mdouble inflowHeight)
 Sets maximum inflow height (Z-direction) More...
 
Mdouble getInflowHeight () const
 Returns the maximum inflow height (Z-direction) More...
 
void setInflowVelocity (Mdouble inflowVelocity)
 Sets the average inflow velocity. More...
 
Mdouble getInflowVelocity () const
 Returns the average inflow velocity. More...
 
void setInflowVelocityVariance (Mdouble inflowVelocityVariance)
 Sets the inflow velocity variance. More...
 
Mdouble getInflowVelocityVariance () const
 Returns the inflow velocity variance. More...
 
void setChuteWidth (Mdouble chuteWidth)
 Sets the chute width (Y-direction) More...
 
Mdouble getChuteWidth () const
 Returns the chute width (Y-direction) More...
 
virtual void setChuteLength (Mdouble chuteLength)
 Sets the chute length (X-direction) More...
 
Mdouble getChuteLength () const
 Returns the chute length (X-direction) More...
 
void setInsertionBoundary (InsertionBoundary *insertionBoundary)
 Sets the chute insertion boundary. More...
 
- Public Member Functions inherited from Mercury3D
 Mercury3D ()
 This is the default constructor. All it does is set sensible defaults. More...
 
 Mercury3D (const DPMBase &other)
 Copy-constructor for creates an Mercury3D problem from an existing MD problem. More...
 
 Mercury3D (const Mercury3D &other)
 Copy-constructor. More...
 
void constructor ()
 Function that sets the SystemDimension and ParticleDimension to 3. More...
 
std::vector< BaseParticle * > hGridFindParticleContacts (const BaseParticle *obj) override
 Returns all particles that have a contact with a given particle. More...
 
- Public Member Functions inherited from MercuryBase
 MercuryBase ()
 This is the default constructor. It sets sensible defaults. More...
 
 ~MercuryBase () override
 This is the default destructor. More...
 
 MercuryBase (const MercuryBase &mercuryBase)
 Copy-constructor. More...
 
void constructor ()
 This is the actual constructor, it is called do both constructors above. More...
 
void hGridActionsBeforeTimeLoop () override
 This sets up the broad phase information, has to be done at this stage because it requires the particle size. More...
 
void hGridActionsBeforeTimeStep () override
 Performs all necessary actions before a time-step, like updating the particles and resetting all the bucket information, etc. More...
 
void read (std::istream &is, ReadOptions opt=ReadOptions::ReadAll) override
 Reads the MercuryBase from an input stream, for example a restart file. More...
 
void write (std::ostream &os, bool writeAllParticles=true) const override
 Writes all data into a restart file. More...
 
Mdouble getHGridCurrentMaxRelativeDisplacement () const
 Returns hGridCurrentMaxRelativeDisplacement_. More...
 
Mdouble getHGridTotalCurrentMaxRelativeDisplacement () const
 Returns hGridTotalCurrentMaxRelativeDisplacement_. More...
 
void setHGridUpdateEachTimeStep (bool updateEachTimeStep)
 Sets whether or not the HGrid must be updated every time step. More...
 
bool getHGridUpdateEachTimeStep () const final
 Gets whether or not the HGrid is updated every time step. More...
 
void setHGridMaxLevels (unsigned int HGridMaxLevels)
 Sets the maximum number of levels of the HGrid in this MercuryBase. More...
 
unsigned int getHGridMaxLevels () const
 Gets the maximum number of levels of the HGrid in this MercuryBase. More...
 
HGridMethod getHGridMethod () const
 Gets whether the HGrid in this MercuryBase is BOTTOMUP or TOPDOWN. More...
 
void setHGridMethod (HGridMethod hGridMethod)
 Sets the HGridMethod to either BOTTOMUP or TOPDOWN. More...
 
HGridDistribution getHGridDistribution () const
 Gets how the sizes of the cells of different levels are distributed. More...
 
void setHGridDistribution (HGridDistribution hGridDistribution)
 Sets how the sizes of the cells of different levels are distributed. More...
 
Mdouble getHGridCellOverSizeRatio () const
 Gets the ratio of the smallest cell over the smallest particle. More...
 
void setHGridCellOverSizeRatio (Mdouble cellOverSizeRatio)
 Sets the ratio of the smallest cell over the smallest particle. More...
 
bool hGridNeedsRebuilding ()
 Gets if the HGrid needs rebuilding before anything else happens. More...
 
virtual unsigned int getHGridTargetNumberOfBuckets () const
 Gets the desired number of buckets, which is the maximum of the number of particles and 10. More...
 
virtual Mdouble getHGridTargetMinInteractionRadius () const
 Gets the desired size of the smallest cells of the HGrid. More...
 
virtual Mdouble getHGridTargetMaxInteractionRadius () const
 Gets the desired size of the largest cells of the HGrid. More...
 
bool checkParticleForInteraction (const BaseParticle &P) final
 Checks if given BaseParticle has an interaction with a BaseWall or other BaseParticle. More...
 
bool checkParticleForInteractionLocal (const BaseParticle &P) final
 Checks if the given BaseParticle has an interaction with a BaseWall or other BaseParticles in a local domain. More...
 
virtual Mdouble userHGridCellSize (unsigned int level)
 Virtual function that enables inheriting classes to implement a function to let the user set the cell size of the HGrid. More...
 
void hGridInfo (std::ostream &os=std::cout) const
 Writes the info of the HGrid to the screen in a nice format. More...
 
- Public Member Functions inherited from DPMBase
void constructor ()
 A function which initialises the member variables to default values, so that the problem can be solved off the shelf; sets up a basic two dimensional problem which can be solved off the shelf. It is called in the constructor DPMBase(). More...
 
 DPMBase ()
 Constructor that calls the "void constructor()". More...
 
 DPMBase (const DPMBase &other)
 Copy constructor type-2. More...
 
virtual ~DPMBase ()
 virtual destructor More...
 
void autoNumber ()
 The autoNumber() function calls three functions: setRunNumber(), readRunNumberFromFile() and incrementRunNumberInFile(). More...
 
std::vector< intget1DParametersFromRunNumber (int size_x) const
 This turns a counter into 1 index, which is a useful feature for performing 1D parameter study. The index run from 1:size_x, while the study number starts at 0 (initially the counter=1 in COUNTER_DONOTDEL) More...
 
std::vector< intget2DParametersFromRunNumber (int size_x, int size_y) const
 This turns a counter into 2 indices which is a very useful feature for performing a 2D study. The indices run from 1:size_x and 1:size_y, while the study number starts at 0 ( initially the counter=1 in COUNTER_DONOTDEL) More...
 
std::vector< intget3DParametersFromRunNumber (int size_x, int size_y, int size_z) const
 This turns a counter into 3 indices, which is a useful feature for performing a 3D parameter study. The indices run from 1:size_x, 1:size_y and 1:size_z, while the study number starts at 0 ( initially the counter=1 in COUNTER_DONOTDEL) More...
 
int launchNewRun (const char *name, bool quick=false)
 This launches a code from within this code. Please pass the name of the code to run. More...
 
void setRunNumber (int runNumber)
 This sets the counter/Run number, overriding the defaults. More...
 
int getRunNumber () const
 This returns the current value of the counter (runNumber_) More...
 
virtual void decompose ()
 Sends particles from processorId to the root processor. More...
 
void solve ()
 The work horse of the code. More...
 
void initialiseSolve ()
 Beginning of the solve routine, before time stepping. More...
 
void finaliseSolve ()
 End of the solve routine, after time stepping. More...
 
virtual void computeOneTimeStep ()
 Performs everything needed for one time step, used in the time-loop of solve(). More...
 
void checkSettings ()
 Checks if the essentials are set properly to go ahead with solving the problem. More...
 
void forceWriteOutputFiles ()
 Writes output files immediately, even if the current time step was not meant to be written. Also resets the last saved time step. More...
 
virtual void writeOutputFiles ()
 Writes simulation data to all the main Mercury files: .data, .ene, .fstat, .xballs and .restart (see the Mercury website for more details regarding these files). More...
 
void solve (int argc, char *argv[])
 The work horse of the code. Can handle flags from the command line. More...
 
ParticleVtkWritergetVtkWriter () const
 
virtual void writeRestartFile ()
 Stores all the particle data for current save time step to a "restart" file, which is a file simply intended to store all the information necessary to "restart" a simulation from a given time step (see also MercuryDPM.org for more information on restart files). More...
 
void writeDataFile ()
 
void writeEneFile ()
 
void writeFStatFile ()
 
void fillDomainWithParticles (unsigned N=50)
 
bool readRestartFile (ReadOptions opt=ReadOptions::ReadAll)
 Reads all the particle data corresponding to a given, existing . restart file (for more details regarding restart files, refer to the training materials on the MercuryDPM website).Returns true if it is successful, false otherwise. More...
 
int readRestartFile (std::string fileName, ReadOptions opt=ReadOptions::ReadAll)
 The same as readRestartFile(bool), but also reads all the particle data corresponding to the current saved time step. More...
 
virtual BaseWallreadUserDefinedWall (const std::string &type) const
 Allows you to read in a wall defined in a Driver directory; see USER/Luca/ScrewFiller. More...
 
virtual void readOld (std::istream &is)
 Reads all data from a restart file, e.g. domain data and particle data; old version. More...
 
bool readDataFile (std::string fileName="", unsigned int format=0)
 This allows particle data to be reloaded from data files. More...
 
bool readParAndIniFiles (std::string fileName)
 Allows the user to read par.ini files (useful to read files produced by the MDCLR simulation code - external to MercuryDPM) More...
 
bool readNextDataFile (unsigned int format=0)
 Reads the next data file with default format=0. However, one can modify the format based on whether the particle data corresponds to 3D or 2D data- see Visualising data in xballs. More...
 
void readNextFStatFile ()
 Reads the next fstat file. More...
 
bool findNextExistingDataFile (Mdouble tMin, bool verbose=true)
 Finds and opens the next data file, if such a file exists. More...
 
bool readArguments (int argc, char *argv[])
 Can interpret main function input arguments that are passed by the driver codes. More...
 
bool checkParticleForInteractionLocalPeriodic (const BaseParticle &P)
 
void readSpeciesFromDataFile (bool read=true)
 
void importParticlesAs (ParticleHandler &particleHandler, InteractionHandler &interactionHandler, const ParticleSpecies *species)
 Copies particles, interactions assigning species from a local simulation to a global one. Useful for the creation of a cluster. More...
 
MERCURYDPM_DEPRECATED FilegetDataFile ()
 The non const version. Allows one to edit the File::dataFile. More...
 
MERCURYDPM_DEPRECATED FilegetEneFile ()
 The non const version. Allows to edit the File::eneFile. More...
 
MERCURYDPM_DEPRECATED FilegetFStatFile ()
 The non const version. Allows to edit the File::fStatFile. More...
 
MERCURYDPM_DEPRECATED FilegetRestartFile ()
 The non const version. Allows to edit the File::restartFile. More...
 
MERCURYDPM_DEPRECATED FilegetStatFile ()
 The non const version. Allows to edit the File::statFile. More...
 
FilegetInteractionFile ()
 Return a reference to the file InteractionsFile. More...
 
MERCURYDPM_DEPRECATED const FilegetDataFile () const
 The const version. Does not allow for any editing of the File::dataFile. More...
 
MERCURYDPM_DEPRECATED const FilegetEneFile () const
 The const version. Does not allow for any editing of the File::eneFile. More...
 
MERCURYDPM_DEPRECATED const FilegetFStatFile () const
 The const version. Does not allow for any editing of the File::fStatFile. More...
 
MERCURYDPM_DEPRECATED const FilegetRestartFile () const
 The const version. Does not allow for any editing of the File::restartFile. More...
 
MERCURYDPM_DEPRECATED const FilegetStatFile () const
 The const version. Does not allow for any editing of the File::statFile. More...
 
const FilegetInteractionFile () const
 Returns a constant reference to an Interactions file. More...
 
const std::string & getName () const
 Returns the name of the file. Does not allow to change it though. More...
 
void setName (const std::string &name)
 Allows to set the name of all the files (ene, data, fstat, restart, stat) More...
 
void setName (const char *name)
 Calls setName(std::string) More...
 
void setSaveCount (unsigned int saveCount)
 Sets File::saveCount_ for all files (ene, data, fstat, restart, stat) More...
 
void setFileType (FileType fileType)
 Sets File::fileType_ for all files (ene, data, fstat, restart, stat) More...
 
void setOpenMode (std::fstream::openmode openMode)
 Sets File::openMode_ for all files (ene, data, fstat, restart, stat) More...
 
void resetFileCounter ()
 Resets the file counter for each file i.e. for ene, data, fstat, restart, stat) More...
 
void closeFiles ()
 Closes all files (ene, data, fstat, restart, stat) that were opened to read or write. More...
 
void setVTKOutputDirectory (const std::string &dir)
 Sets the output directory of the VTK files. More...
 
void setLastSavedTimeStep (unsigned int nextSavedTimeStep)
 Sets the next time step for all the files (ene, data, fstat, restart, stat) at which the data is to be written or saved. More...
 
Mdouble getTime () const
 Returns the current simulation time. More...
 
Mdouble getNextTime () const
 Returns the current simulation time. More...
 
unsigned int getNumberOfTimeSteps () const
 Returns the current counter of time-steps, i.e. the number of time-steps that the simulation has undergone so far. More...
 
void setTime (Mdouble time)
 Sets a new value for the current simulation time. More...
 
void setTimeMax (Mdouble newTMax)
 Sets a new value for the maximum simulation duration. More...
 
Mdouble getTimeMax () const
 Returns the maximum simulation duration. More...
 
void setLogarithmicSaveCount (Mdouble logarithmicSaveCountBase)
 Sets File::logarithmicSaveCount_ for all files (ene, data, fstat, restart, stat) More...
 
void setNToWrite (int nToWrite)
 set the number of elements to write to the screen More...
 
int getNToWrite () const
 get the number of elements to write to the More...
 
void setRotation (bool rotation)
 Sets whether particle rotation is enabled or disabled. More...
 
bool getRotation () const
 Indicates whether particle rotation is enabled or disabled. More...
 
MERCURYDPM_DEPRECATED void setWallsWriteVTK (FileType writeWallsVTK)
 Sets whether walls are written into a VTK file. More...
 
MERCURYDPM_DEPRECATED void setWallsWriteVTK (bool)
 Sets whether walls are written into a VTK file. More...
 
MERCURYDPM_DEPRECATED void setInteractionsWriteVTK (bool)
 Sets whether interactions are written into a VTK file. More...
 
void setParticlesWriteVTK (bool writeParticlesVTK)
 Sets whether particles are written in a VTK file. More...
 
void setSuperquadricParticlesWriteVTK (bool writeSuperquadricParticlesVTK)
 
MERCURYDPM_DEPRECATED FileType getWallsWriteVTK () const
 Returns whether walls are written in a VTK file. More...
 
bool getParticlesWriteVTK () const
 Returns whether particles are written in a VTK file. More...
 
bool getSuperquadricParticlesWriteVTK () const
 
Mdouble getXMin () const
 If the length of the problem domain in x-direction is XMax - XMin, then getXMin() returns XMin. More...
 
Mdouble getXMax () const
 If the length of the problem domain in x-direction is XMax - XMin, then getXMax() returns XMax. More...
 
Mdouble getYMin () const
 If the length of the problem domain in y-direction is YMax - YMin, then getYMin() returns YMin. More...
 
Mdouble getYMax () const
 If the length of the problem domain in y-direction is YMax - YMin, then getYMax() returns XMax. More...
 
Mdouble getZMin () const
 If the length of the problem domain in z-direction is ZMax - ZMin, then getZMin() returns ZMin. More...
 
Mdouble getZMax () const
 If the length of the problem domain in z-direction is ZMax - ZMin, then getZMax() returns ZMax. More...
 
Mdouble getXCenter () const
 
Mdouble getYCenter () const
 
Mdouble getZCenter () const
 
Vec3D getCenter () const
 get center of domain More...
 
Vec3D getMin () const
 Returns the minimum coordinates of the problem domain. More...
 
Vec3D getMax () const
 Returns the maximum coordinates of the problem domain. More...
 
void setXMin (Mdouble newXMin)
 Sets the value of XMin, the lower bound of the problem domain in the x-direction. More...
 
void setYMin (Mdouble newYMin)
 Sets the value of YMin, the lower bound of the problem domain in the y-direction. More...
 
void setZMin (Mdouble newZMin)
 Sets the value of ZMin, the lower bound of the problem domain in the z-direction. More...
 
void setXMax (Mdouble newXMax)
 Sets the value of XMax, the upper bound of the problem domain in the x-direction. More...
 
void setYMax (Mdouble newYMax)
 Sets the value of YMax, the upper bound of the problem domain in the y-direction. More...
 
void setZMax (Mdouble newZMax)
 Sets the value of ZMax, the upper bound of the problem domain in the z-direction. More...
 
void setMax (const Vec3D &max)
 Sets the maximum coordinates of the problem domain. More...
 
void setMax (Mdouble, Mdouble, Mdouble)
 Sets the maximum coordinates of the problem domain. More...
 
void setDomain (const Vec3D &min, const Vec3D &max)
 Sets the minimum coordinates of the problem domain. More...
 
void setMin (const Vec3D &min)
 Sets the minimum coordinates of the problem domain. More...
 
void setMin (Mdouble, Mdouble, Mdouble)
 Sets the minimum coordinates of the problem domain. More...
 
void setTimeStep (Mdouble newDt)
 Sets a new value for the simulation time step. More...
 
Mdouble getTimeStep () const
 Returns the simulation time step. More...
 
void setNumberOfOMPThreads (int numberOfOMPThreads)
 Sets the number of omp threads. More...
 
int getNumberOfOMPThreads () const
 Returns the number of omp threads. More...
 
void setXBallsColourMode (int newCMode)
 Set the xballs output mode. More...
 
int getXBallsColourMode () const
 Get the xballs colour mode (CMode). More...
 
void setXBallsVectorScale (double newVScale)
 Set the scale of vectors in xballs. More...
 
double getXBallsVectorScale () const
 Returns the scale of vectors used in xballs. More...
 
void setXBallsAdditionalArguments (std::string newXBArgs)
 Set the additional arguments for xballs. More...
 
std::string getXBallsAdditionalArguments () const
 Returns the additional arguments for xballs. More...
 
void setXBallsScale (Mdouble newScale)
 Sets the scale of the view (either normal, zoom in or zoom out) to display in xballs. The default is fit to screen. More...
 
double getXBallsScale () const
 Returns the scale of the view in xballs. More...
 
void setGravity (Vec3D newGravity)
 Sets a new value for the gravitational acceleration. More...
 
Vec3D getGravity () const
 Returns the gravitational acceleration. More...
 
void setBackgroundDrag (Mdouble backgroundDrag)
 Simple access function to turn on a background drag. The force of particleVelocity*drag is applied (note, it allowed to be negative i.e. create energy) More...
 
const Mdouble getBackgroundDrag () const
 Return the background drag. More...
 
void setDimension (unsigned int newDim)
 Sets both the system dimensions and the particle dimensionality. More...
 
void setSystemDimensions (unsigned int newDim)
 Sets the system dimensionality. More...
 
unsigned int getSystemDimensions () const
 Returns the system dimensionality. More...
 
void setParticleDimensions (unsigned int particleDimensions)
 Sets the particle dimensionality. More...
 
unsigned int getParticleDimensions () const
 Returns the particle dimensionality. More...
 
std::string getRestartVersion () const
 This is to take into account for different Mercury versions. Returns the version of the restart file. More...
 
void setRestartVersion (std::string newRV)
 Sets restart_version. More...
 
bool getRestarted () const
 Returns the flag denoting if the simulation was restarted or not. More...
 
void setRestarted (bool newRestartedFlag)
 Allows to set the flag stating if the simulation is to be restarted or not. More...
 
bool getAppend () const
 Returns whether the "append" option is on or off. More...
 
void setAppend (bool newAppendFlag)
 Sets whether the "append" option is on or off. More...
 
Mdouble getElasticEnergy () const
 Returns the global elastic energy within the system. More...
 
Mdouble getKineticEnergy () const
 Returns the global kinetic energy stored in the system. More...
 
Mdouble getGravitationalEnergy (Vec3D origin={0, 0, 0}) const
 Returns the global gravitational potential energy stored relative to a user-defined origin in the system. Note, if no origin is specified the real origin i.e. (0,0,0) is taken. More...
 
Mdouble getRotationalEnergy () const
 Returns the global rotational energy stored in the system. More...
 
Mdouble getTotalEnergy () const
 
Mdouble getTotalMass () const
 JMFT: Return the total mass of the system, excluding fixed particles. More...
 
Vec3D getCentreOfMass () const
 JMFT: Return the centre of mass of the system, excluding fixed particles. More...
 
Vec3D getTotalMomentum () const
 JMFT: Return the total momentum of the system, excluding fixed particles. More...
 
double getCPUTime ()
 
double getWallTime ()
 
virtual void hGridInsertParticle (BaseParticle *obj UNUSED)
 
virtual void hGridUpdateParticle (BaseParticle *obj UNUSED)
 
virtual void hGridRemoveParticle (BaseParticle *obj UNUSED)
 
bool mpiIsInCommunicationZone (BaseParticle *particle)
 Checks if the position of the particle is in an mpi communication zone or not. More...
 
bool mpiInsertParticleCheck (BaseParticle *P)
 Function that checks if the mpi particle should really be inserted by the current domain. More...
 
void insertGhostParticle (BaseParticle *P)
 This function inserts a particle in the mpi communication boundaries. More...
 
void updateGhostGrid (BaseParticle *P)
 Checks if the Domain/periodic interaction distance needs to be updated and updates it accordingly. More...
 
virtual void gatherContactStatistics (unsigned int index1, int index2, Vec3D Contact, Mdouble delta, Mdouble ctheta, Mdouble fdotn, Mdouble fdott, Vec3D P1_P2_normal_, Vec3D P1_P2_tangential)
 //Not unsigned index because of possible wall collisions. More...
 
void setNumberOfDomains (std::vector< unsigned > direction)
 Sets the number of domains in x-,y- and z-direction. Required for parallel computations. More...
 
void splitDomain (DomainSplit domainSplit)
 
std::vector< unsignedgetNumberOfDomains ()
 returns the number of domains More...
 
DomaingetCurrentDomain ()
 Function that returns a pointer to the domain corresponding to the processor. More...
 
void removeOldFiles () const
 
void setMeanVelocity (Vec3D V_mean_goal)
 This function will help you set a fixed kinetic energy and mean velocity in your system. More...
 
void setMeanVelocityAndKineticEnergy (Vec3D V_mean_goal, Mdouble Ek_goal)
 This function will help you set a fixed kinetic energy and mean velocity in your system. More...
 
Mdouble getTotalVolume () const
 Get the total volume of the cuboid system. More...
 
Matrix3D getKineticStress () const
 Calculate the kinetic stress tensor in the system averaged over the whole volume. More...
 
Matrix3D getStaticStress () const
 Calculate the static stress tensor in the system averaged over the whole volume. More...
 
Matrix3D getTotalStress () const
 Calculate the total stress tensor in the system averaged over the whole volume. More...
 
virtual void handleParticleRemoval (unsigned int id)
 Handles the removal of particles from the particleHandler. More...
 
virtual void handleParticleAddition (unsigned int id, BaseParticle *p)
 Handles the addition of particles to the particleHandler. More...
 
void writePythonFileForVTKVisualisation () const
 writes .py file for ParaView More...
 
void setWritePythonFileForVTKVisualisation (bool forceWritePythonFileForVTKVisualisation)
 
bool getWritePythonFileForVTKVisualisation () const
 
WallVTKWritergetWallVTKWriter ()
 

Public Attributes

SphericalParticle inflowParticle_
 
- Public Attributes inherited from DPMBase
SpeciesHandler speciesHandler
 A handler to that stores the species type i.e. LinearViscoelasticSpecies, etc. More...
 
RNG random
 This is a random generator, often used for setting up the initial conditions etc... More...
 
ParticleHandler particleHandler
 An object of the class ParticleHandler, contains the pointers to all the particles created. More...
 
ParticleHandler paoloParticleHandler
 Fake particleHandler created by Paolo needed temporary by just Paolo. More...
 
WallHandler wallHandler
 An object of the class WallHandler. Contains pointers to all the walls created. More...
 
BoundaryHandler boundaryHandler
 An object of the class BoundaryHandler which concerns insertion and deletion of particles into or from regions. More...
 
PeriodicBoundaryHandler periodicBoundaryHandler
 Internal handler that deals with periodic boundaries, especially in a parallel build. More...
 
DomainHandler domainHandler
 An object of the class DomainHandler which deals with parallel code. More...
 
InteractionHandler interactionHandler
 An object of the class InteractionHandler. More...
 
CGHandler cgHandler
 Object of the class cgHandler. More...
 
File dataFile
 An instance of class File to handle in- and output into a .data file. More...
 
File fStatFile
 An instance of class File to handle in- and output into a .fstat file. More...
 
File eneFile
 An instance of class File to handle in- and output into a .ene file. More...
 
File restartFile
 An instance of class File to handle in- and output into a .restart file. More...
 
File statFile
 An instance of class File to handle in- and output into a .stat file. More...
 
File interactionFile
 File class to handle in- and output into .interactions file. This file hold information about interactions. More...
 
Time clock_
 record when the simulation started More...
 

Private Attributes

double PeriodicBoxLength
 stores the length of the periodic box More...
 
int PeriodicBoxNSpecies
 stores the number of species in the periodic box More...
 
bool FillChute
 

Additional Inherited Members

- Public Types inherited from DPMBase
enum class  ReadOptions : int { ReadAll , ReadNoInteractions , ReadNoParticlesAndInteractions }
 
enum class  DomainSplit {
  X , Y , Z , XY ,
  XZ , YZ , XYZ
}
 
- Static Public Member Functions inherited from DPMBase
static void incrementRunNumberInFile ()
 Increment the run Number (counter value) stored in the file_counter (COUNTER_DONOTDEL) by 1 and store the new value in the counter file. More...
 
static int readRunNumberFromFile ()
 Read the run number or the counter from the counter file (COUNTER_DONOTDEL) More...
 
static bool areInContact (const BaseParticle *pI, const BaseParticle *pJ)
 Checks if two particle are in contact or is there any positive overlap. More...
 
- Protected Member Functions inherited from Chute
void actionsBeforeTimeStep () override
 Calls Chute::cleanChute(). More...
 
void cleanChute ()
 Deletes all outflow particles once every 100 time steps. More...
 
virtual void createBottom ()
 Creates the chute bottom, which can be either flat or one of three flavours of rough. More...
 
virtual void addFlowParticlesCompactly ()
 Add initial flow particles in a dense packing. More...
 
virtual SphericalParticle createFlowParticle ()
 
void printTime () const override
 prints time, max time and number of particles More...
 
- Protected Member Functions inherited from Mercury3D
void hGridFindContactsWithinTargetCell (int x, int y, int z, unsigned int l)
 Finds contacts between particles in the target cell. More...
 
void hGridFindContactsWithTargetCell (int x, int y, int z, unsigned int l, BaseParticle *obj)
 Finds contacts between the BaseParticle and the target cell. More...
 
void computeWallForces (BaseWall *w) override
 Compute contacts with a wall. More...
 
void hGridFindParticlesWithTargetCell (int x, int y, int z, unsigned int l, BaseParticle *obj, std::vector< BaseParticle * > &list)
 Finds particles within target cell and stores them in a list. More...
 
void hGridGetInteractingParticleList (BaseParticle *obj, std::vector< BaseParticle * > &list) override
 Obtains all neighbour particles of a given object, obtained from the hgrid. More...
 
void computeInternalForces (BaseParticle *obj) override
 Finds contacts with the BaseParticle; avoids multiple checks. More...
 
bool hGridHasContactsInTargetCell (int x, int y, int z, unsigned int l, const BaseParticle *obj) const
 Tests if the BaseParticle has contacts with other Particles in the target cell. More...
 
bool hGridHasParticleContacts (const BaseParticle *obj) override
 Tests if a BaseParticle has any contacts in the HGrid. More...
 
void hGridRemoveParticle (BaseParticle *obj) override
 Removes a BaseParticle from the HGrid. More...
 
void hGridUpdateParticle (BaseParticle *obj) override
 Updates the cell (not the level) of a BaseParticle. More...
 
- Protected Member Functions inherited from MercuryBase
void hGridRebuild ()
 This sets up the parameters required for the contact model. More...
 
void hGridInsertParticle (BaseParticle *obj) final
 Inserts a single Particle to current grid. More...
 
void hGridUpdateMove (BaseParticle *iP, Mdouble move) final
 Computes the relative displacement of the given BaseParticle and updates the currentMaxRelativeDisplacement_ accordingly. More...
 
void hGridActionsBeforeIntegration () override
 Resets the currentMaxRelativeDisplacement_ to 0. More...
 
void hGridActionsAfterIntegration () override
 This function has to be called before integrateBeforeForceComputation. More...
 
HGridgetHGrid ()
 Gets the HGrid used by this problem. More...
 
const HGridgetHGrid () const
 Gets the HGrid used by this problem, const version. More...
 
bool readNextArgument (int &i, int argc, char *argv[]) override
 Reads the next command line argument. More...
 
- Protected Member Functions inherited from DPMBase
virtual void computeAllForces ()
 Computes all the forces acting on the particles using the BaseInteractable::setForce() and BaseInteractable::setTorque() More...
 
virtual void computeInternalForce (BaseParticle *, BaseParticle *)
 Computes the forces between two particles (internal in the sense that the sum over all these forces is zero i.e. fully modelled forces) More...
 
virtual void computeExternalForces (BaseParticle *)
 Computes the external forces, such as gravity, acting on particles. More...
 
virtual void computeForcesDueToWalls (BaseParticle *, BaseWall *)
 Computes the forces on the particles due to the walls (normals are outward normals) More...
 
virtual void actionsOnRestart ()
 A virtual function where the users can add extra code which is executed only when the code is restarted. More...
 
virtual void actionsBeforeTimeLoop ()
 A virtual function. Allows one to carry out any operations before the start of the time loop. More...
 
virtual void computeAdditionalForces ()
 A virtual function which allows to define operations to be executed prior to the OMP force collect. More...
 
virtual void actionsAfterSolve ()
 A virtual function which allows to define operations to be executed after the solve(). More...
 
virtual void actionsAfterTimeStep ()
 A virtual function which allows to define operations to be executed after time step. More...
 
void writeVTKFiles () const
 
virtual void outputXBallsData (std::ostream &os) const
 This function writes the location of the walls and particles in a format the XBalls program can read. For more information on the XBalls program, see Visualising data in xballs. More...
 
virtual void outputXBallsDataParticle (unsigned int i, unsigned int format, std::ostream &os) const
 This function writes out the particle locations into an output stream in a format the XBalls program can read. For more information on the XBalls program, see Visualising data in xballs. More...
 
virtual void writeEneHeader (std::ostream &os) const
 Writes a header with a certain format for ENE file. More...
 
virtual void writeFstatHeader (std::ostream &os) const
 Writes a header with a certain format for FStat file. More...
 
virtual void writeEneTimeStep (std::ostream &os) const
 Write the global kinetic, potential energy, etc. in the system. More...
 
virtual void initialiseStatistics ()
 
virtual void outputStatistics ()
 
void gatherContactStatistics ()
 
virtual void processStatistics (bool)
 
virtual void finishStatistics ()
 
virtual void integrateAfterForceComputation ()
 Update particles' and walls' positions and velocities after force computation. More...
 
virtual void checkInteractionWithBoundaries ()
 There are a range of boundaries one could implement depending on ones' problem. This methods checks for interactions between particles and such range of boundaries. See BaseBoundary.h and all the boundaries in the Boundaries folder. More...
 
void setFixedParticles (unsigned int n)
 Sets a number, n, of particles in the particleHandler as "fixed particles". More...
 
virtual bool continueSolve () const
 A virtual function for deciding whether to continue the simulation, based on a user-specified criterion. More...
 
void outputInteractionDetails () const
 Displays the interaction details corresponding to the pointer objects in the interaction handler. More...
 
bool isTimeEqualTo (Mdouble time) const
 Checks whether the input variable "time" is the current time in the simulation. More...
 
void removeDuplicatePeriodicParticles ()
 Removes periodic duplicate Particles. More...
 
void checkAndDuplicatePeriodicParticles ()
 For simulations using periodic boundaries, checks and adds particles when necessary into the particle handler. See DPMBase.cc and PeriodicBoundary.cc for more details. More...
 
void performGhostParticleUpdate ()
 When the Verlet scheme updates the positions and velocities of particles, ghost particles will need an update as wel. Their status will also be updated accordingly. More...
 
void deleteGhostParticles (std::set< BaseParticle * > &particlesToBeDeleted)
 
void synchroniseParticle (BaseParticle *, unsigned fromProcessor=0)
 
void performGhostVelocityUpdate ()
 updates the final time-step velocity of the ghost particles More...
 
void disableSoftStop ()
 This prevents the initialisation of the soft stop feature. More...
 
void discontinueSolve ()
 This discontinues the solve loop after the current time step is finished. More...
 

Detailed Description

Particles of a single Species.

Constructor & Destructor Documentation

◆ ChuteWithPeriodicInflow() [1/3]

ChuteWithPeriodicInflow::ChuteWithPeriodicInflow ( std::string  restart_file)
inline
20  {
21  loadPeriodicBox(restart_file);
22  }
void loadPeriodicBox(std::string const restart_file)
loads periodic chute data from restart file
Definition: ChuteWithPeriodicInflow.h:43

References loadPeriodicBox().

◆ ChuteWithPeriodicInflow() [2/3]

ChuteWithPeriodicInflow::ChuteWithPeriodicInflow ( std::string  restart_file,
int  numRepetitions 
)
inline
25  {
26  loadPeriodicBox(restart_file);
27  AddContinuingBottom(numRepetitions);
28  }
void AddContinuingBottom(int numRepetitions)
Definition: ChuteWithPeriodicInflow.h:70

References AddContinuingBottom(), and loadPeriodicBox().

◆ ChuteWithPeriodicInflow() [3/3]

ChuteWithPeriodicInflow::ChuteWithPeriodicInflow ( std::string  restart_file,
int  numRepetitions,
int  numRepetitionsInWidth 
)
inline
31  {
32  loadPeriodicBox(restart_file);
33  AddContinuingBottom(numRepetitions);
34  ExtendInWidth(numRepetitionsInWidth);
35  }
void ExtendInWidth(int numRepetitionsInWidth)
Definition: ChuteWithPeriodicInflow.h:96

References AddContinuingBottom(), ExtendInWidth(), and loadPeriodicBox().

Member Function Documentation

◆ actionsBeforeTimeStep()

void ChuteWithPeriodicInflow::actionsBeforeTimeStep ( )
inlineoverridevirtual

Do not add, only remove particles.

Reimplemented from DPMBase.

121  {
122  cleanChute();
123  }
void cleanChute()
Remove particles if they fall below a certain height (allows them to become supercritical)
Definition: ChuteWithPeriodicInflow.h:126

References cleanChute().

◆ AddContinuingBottom()

void ChuteWithPeriodicInflow::AddContinuingBottom ( int  numRepetitions)
inline
71  {
72  // creates bottom outside periodic chute of species 1
73  double lengthPeriodicChute = getXMax();
75  particleHandler.setStorageCapacity(N * (2. + numRepetitions));
76  //allows for a gap between the infow and flow regimes
77  double Gap = 0.;
78 
79  for (int j = 1; j < numRepetitions; j++)
80  {
81  for (int i = 0; i < N; i++)
82  {
84  {
87  particleHandler.getLastObject()->move(Vec3D(Gap + lengthPeriodicChute * j, 0., 0.));
88  }
89  }
90  }
91 
92  setXMax(numRepetitions * lengthPeriodicChute);
93  //setHGridNumberOfBucketsToPower(particleHandler.getStorageCapacity());
94  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
void setStorageCapacity(const unsigned int N)
Sets the storage capacity of this BaseHandler.
Definition: BaseHandler.h:677
T * getObject(const unsigned int id)
Gets a pointer to the Object at the specified index in the BaseHandler.
Definition: BaseHandler.h:621
T * getLastObject()
Gets a pointer to the last Object in this BaseHandler.
Definition: BaseHandler.h:642
virtual void move(const Vec3D &move)
Moves this BaseInteractable by adding an amount to the position.
Definition: BaseInteractable.cc:193
bool isFixed() const override
Is fixed Particle function. It returns whether a Particle is fixed or not, by checking its inverse Ma...
Definition: BaseParticle.h:72
virtual BaseParticle * copy() const =0
Particle copy method. It calls to copy constructor of this Particle, useful for polymorphism.
virtual void setSpecies(const ParticleSpecies *species)
Definition: BaseParticle.cc:798
bool FillChute
Definition: ChuteWithPeriodicInflow.h:880
Mdouble getXMax() const
If the length of the problem domain in x-direction is XMax - XMin, then getXMax() returns XMax.
Definition: DPMBase.h:610
SpeciesHandler speciesHandler
A handler to that stores the species type i.e. LinearViscoelasticSpecies, etc.
Definition: DPMBase.h:1433
void setXMax(Mdouble newXMax)
Sets the value of XMax, the upper bound of the problem domain in the x-direction.
Definition: DPMBase.cc:1156
ParticleHandler particleHandler
An object of the class ParticleHandler, contains the pointers to all the particles created.
Definition: DPMBase.h:1443
void addObject(BaseParticle *P) override
Adds a BaseParticle to the ParticleHandler.
Definition: ParticleHandler.cc:150
unsigned int getNumberOfObjects() const override
Returns the number of objects in the container. In parallel code this practice is forbidden to avoid ...
Definition: ParticleHandler.cc:1323
Definition: Kernel/Math/Vector.h:30
@ N
Definition: constructor.cpp:22
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References ParticleHandler::addObject(), BaseParticle::copy(), FillChute, BaseHandler< T >::getLastObject(), ParticleHandler::getNumberOfObjects(), BaseHandler< T >::getObject(), DPMBase::getXMax(), i, BaseParticle::isFixed(), j, BaseInteractable::move(), N, DPMBase::particleHandler, BaseParticle::setSpecies(), BaseHandler< T >::setStorageCapacity(), DPMBase::setXMax(), and DPMBase::speciesHandler.

Referenced by ChuteWithPeriodicInflow().

◆ Check_and_Duplicate_Periodic_Particle()

int ChuteWithPeriodicInflow::Check_and_Duplicate_Periodic_Particle ( BaseParticle i,
int  nWallPeriodic 
)
inline

assumes first periodic wall is in x direction

728  {
729  int C = 0; //Number of particles created
730  for (int k = nWallPeriodic; k < getIsPeriodic(); k++)
731  { //Loop over all still posible walls
734  if ((k ? true : InPeriodicBox(i)) && (perw->getDistance(*i) < i->getRadius() + getMaxInflowParticleRadius()))
735  {
736  SphericalParticle F0 = *i;
737  perw->shiftPosition(i);
738 
739  //If Particle is Mdouble shifted, get correct original particle
740  BaseParticle* From = i;
741  while (From->getPeriodicFromParticle() != NULL)
742  From = From->getPeriodicFromParticle();
743  F0.setPeriodicFromParticle(From);
744 
746  C++;
747 
748  //Check for Mdouble shifted particles
750  }
751  }
752  return (C);
753  }
std::enable_if<!std::is_pointer< U >::value, U * >::type copyAndAddObject(const U &object)
Creates a copy of a Object and adds it to the BaseHandler.
Definition: BaseHandler.h:360
Definition: BaseParticle.h:33
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:416
BaseParticle * getPeriodicFromParticle() const
Returns the 'original' particle this one's a periodic copy of.
Definition: BaseParticle.h:324
int Check_and_Duplicate_Periodic_Particle(BaseParticle *i, int nWallPeriodic)
Definition: ChuteWithPeriodicInflow.h:727
bool InPeriodicBox(BaseParticle *P)
Definition: ChuteWithPeriodicInflow.h:867
Mdouble getMaxInflowParticleRadius() const
Returns the maximum radius of inflow particles.
Definition: Chute.cc:926
bool getIsPeriodic() const
Returns whether the chute is periodic in Y.
Definition: Chute.cc:621
BoundaryHandler boundaryHandler
An object of the class BoundaryHandler which concerns insertion and deletion of particles into or fro...
Definition: DPMBase.h:1458
Defines a pair of periodic walls. Inherits from BaseBoundary.
Definition: PeriodicBoundary.h:20
Mdouble getDistance(const BaseParticle &p) const override
Returns the distance of the edge to the particle.
Definition: PeriodicBoundary.cc:176
virtual void shiftPosition(BaseParticle *p) const override
shifts the particle
Definition: PeriodicBoundary.cc:198
A spherical particle is the most simple particle used in MercuryDPM.
Definition: SphericalParticle.h:16
Definition: matrices.h:74
char char char int int * k
Definition: level2_impl.h:374

References DPMBase::boundaryHandler, BaseHandler< T >::copyAndAddObject(), PeriodicBoundary::getDistance(), Chute::getIsPeriodic(), BaseHandler< T >::getLastObject(), Chute::getMaxInflowParticleRadius(), BaseHandler< T >::getObject(), BaseParticle::getPeriodicFromParticle(), i, InPeriodicBox(), k, DPMBase::particleHandler, BaseParticle::setPeriodicFromParticle(), and PeriodicBoundary::shiftPosition().

◆ cleanChute()

void ChuteWithPeriodicInflow::cleanChute ( )
inline

Remove particles if they fall below a certain height (allows them to become supercritical)

127  {
128  //clean outflow every 100 time steps
129  static int count = 0, maxcount = 100;
130  if (count > maxcount)
131  {
132  count = 0;
133  // delete all outflowing particles
134  for (unsigned int i = 0; i < particleHandler.getNumberOfObjects();)
135  {
136  //~ if (particleHandler.getObject(i)->getPosition().Z<getZMin()-10*getInflowParticleRadius()){
138  {
139  //~ cout << "Remove particle" << endl;
141  }
142  else
143  i++;
144  }
145  }
146  else
147  count++;
148  }
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:197
void removeObject(unsigned int index) override
Removes a BaseParticle from the ParticleHandler.
Definition: ParticleHandler.cc:388
Mdouble X
the vector components
Definition: Kernel/Math/Vector.h:45

References ParticleHandler::getNumberOfObjects(), BaseHandler< T >::getObject(), BaseInteractable::getPosition(), DPMBase::getXMax(), i, DPMBase::particleHandler, ParticleHandler::removeObject(), and Vec3D::X.

Referenced by actionsBeforeTimeStep().

◆ computeInternalForces()

void ChuteWithPeriodicInflow::computeInternalForces ( BaseParticle P1,
BaseParticle P2 
)
inline

Tangential spring information is always store in the real particle with highest index When a Periodic contact is encountered it is always encoutered twice (for normal periodic walls), but the force is only applied to the real particle The tangential spring information for periodic particles is stored in the normal particle (and thus twice for normal periodic interactions) When a Particle is removed the tangential spring information has to be moved

Both options are up to first order the same (the first one is nicer becaus it always keeps the spring tangential, whereas the second one is in a nicer intergration form)

Todo:
{The spring should be cut back such that fdott=mu*fdotn. This is simple for getSlidingDissipation()=0; we have to think about what happens in the sliding case with tang. dissipation; same for walls; Update Dec-2-2011: fixed}
Todo:
TW: the following 13 lines concern only the sliding spring and could be moved into the if statement above
Todo:
Define the center this way or are radii involved? Or maybe just use middle of overlap region? Thomas: Yes, we should correct that for polydispersed problems; however, then we also have to correct it in StatisticsVector::gatherContactStatistics.

The second particle (i.e. the particle the force acts on) is always a flow particle

Todo:
{Is it the first particle the force acts on?}
245  {
246 <<<<<<< .mine
247  // For particles of the same species, set species vector to Species(PI);
248  // for particles of different species, set species vector to MixedSpecies(PI,PJ)
249  BaseSpecies* pSpecies = (PI->getIndSpecies() == PJ->getIndSpecies()) ? &Species[PI->getIndSpecies()] : speciesHandler.getMixedObject(PI->getIndSpecies(), PJ->getIndSpecies());
250 =======
251  //this is because the rough bottom allows overlapping fixed particles
252  if (P1->isFixed() && P2->isFixed())
253  return;
254 >>>>>>> .r824
255 
260 
261  //Dificult cases for tangential springs in comination with periodic walls are:
262  //
263  // A new contact over a periodic wall:
264  // Starting situation: There are no tangential springs stored at all
265  // Creating periodic particles: There are no tangential springs so nothing happens
266  // Contact detection: There are two contacts detected, one (CA) for P1 with P2per and one (CB) for P2 with P1per
267  // Switch to the 4 cases:
268  // CA: PI=P1, PJ=P2per PJreal=P2
269  // CB: PI=P2, PJ=P1per PJreal=P1
270  // Reconnecting springs:
271  // CA: PI=P1 and PJ=P2per do not have a spring, so a new spring is created in either PI or PJ (could be in periodic particle)
272  // CB: PI=P2 and PJ=P1Per do not have a spring, so a new spring is created in either PI or PJ (could be in periodic particle)
273  // Removing periodic particles: One of the springs will be stored in a periodic particle and thus removed, the other spring is kept and used (this is the real particle with the kighest index)
274 
275  // Reconnect to a contact over a periodic wall:
276  // Starting situation: There is a tangential spring stored in the particle with the highest index, lets assume this is P1
277  // Creating periodic particles: P1 has a tangential spring, so P1per has a reversed copy of this spring
278  // Contact detection: There are two contacts detected, one (CA) for P1 with P2per and one (CB) for P2 with P1per
279  // Switch to the 4 cases:
280  // CA: PI=P1, PJ=P2per PJreal=P2
281  // CB: PI=P2, PJ=P1per PJreal=P1
282  // Reconnecting springs:
283  // CA: PI=P1 and PJ=P2per have a spring in P1, so we reconnect to this spring in the normal way and integrate it.
284  // CB: PI=P2 and PJ=P1Per have a spring in P1per, however this spring has to be a reversed spring since it is in PJ.
285  // Removing periodic particles: The spring in P1per is removed
286 
287  //Cases (start from P1 and P2 and go to PI and PJ
288  //1 Normal-Normal ->PI=max(P1,P2), PJ=min(P1,P2)
289  //2 Periodic-Normal ->(PI=P2 PJ=real(P1))
290  //3 Normal-Periodic ->(PI=P1 PJ=real(P2))
291  //4 Periodic-Periodic ->do nothing
292 
293  //Just some statements to handle the 4 cases
294  BaseParticle *PI, *PJ, *PJreal;
295  bool isPeriodic;
296  BaseParticle *P1Per = P1->getPeriodicFromParticle();
297  BaseParticle *P2Per = P2->getPeriodicFromParticle();
298  if (P1Per == NULL)
299  {
300  if (P2Per == NULL) //N-N
301  if (P1->getIndex() < P2->getIndex())
302  {
303  PI = P2;
304  PJ = P1;
305  PJreal = P1;
306  isPeriodic = false;
307  }
308  else
309  {
310  PI = P1;
311  PJ = P2;
312  PJreal = P2;
313  isPeriodic = false;
314  }
315  else //N-P
316  {
317  PI = P1;
318  PJ = P2;
319  PJreal = P2Per;
320  isPeriodic = true;
321  }
322  }
323  else
324  {
325  if (P2Per == NULL) //P-N
326  {
327  PI = P2;
328  PJ = P1;
329  PJreal = P1Per;
330  isPeriodic = true;
331  }
332  else
333  return;
334  }
335 
336 #ifdef DEBUG_OUTPUT
337  std::cerr << "In computing internal forces between particle "<<PI->getPosition()<<" and "<<PJ->getPosition()<<std::endl;
338 #endif
339 
340  //Get the square of the distance between particle i and particle j
341  Mdouble dist_squared = Vec3D::getDistanceSquared(PI->getPosition(), PJ->getPosition());
342  Mdouble interactionRadii_sum = PI->getInteractionRadius() + PJ->getInteractionRadius();
343 
344 #ifdef DEBUG_OUTPUT_FULL
345  std::cerr << "Square of distance between " << dist_squared << " square sum of radii " << radii_sum*radii_sum <<std::endl;
346 #endif
347 
348  // True if the particles are in contact
349  if (dist_squared < (interactionRadii_sum * interactionRadii_sum))
350  {
351  // For particles of the same species, set species vector to Species(PI);
352  // for particles of different species, set species vector to MixedSpecies(PI,PJ)
353  CSpecies* pSpecies = (PI->getIndSpecies() == PJ->getIndSpecies()) ? &Species[PI->getIndSpecies()] : getMixedSpecies(PI->getIndSpecies(), PJ->getIndSpecies());
354 
355  // Calculate distance between the particles
356  Mdouble dist = sqrt(dist_squared);
357 
358  // Compute normal vector
359 
360  Vec3D normal = (PI->getPosition() - PJ->getPosition()) / dist;
361 
362  // Compute the overlap between the particles
363  Mdouble radii_sum = PI->getRadius() + PJ->getRadius();
364  Mdouble deltan = std::max(0.0, radii_sum - dist);
365 
366  Vec3D force = Vec3D(0, 0, 0);
367  ;
368  Vec3D forcet;
369  forcet.setZero();
370  Vec3D forcerolling;
371  forcerolling.setZero();
372  Vec3D forcetorsion;
373  forcetorsion.setZero();
374  Mdouble fdotn = 0;
375  CTangentialSpring* TangentialSpring = NULL;
376 
377  //evaluate shortrange non-contact forces
378  if (pSpecies->getAdhesionForceType() != AdhesionForceType::NONE)
379  fdotn += computeShortRangeForceWithParticle(PI, PJ, PJreal, pSpecies, dist);
380 
381  if (deltan > 0) //if contact forces
382  {
383 
384  // Compute the relative velocity vector v_ij
385  Vec3D vrel;
386  if (!pSpecies->getSlidingFrictionCoefficient())
387  {
388  vrel = (PI->getVelocity() - PJ->getVelocity());
389  }
390  else
391  {
392  vrel = (PI->getVelocity() - PJ->getVelocity()) + Vec3D::cross(normal, PI->getAngularVelocity() * (PI->getRadius() - .5 * deltan) + PJ->getAngularVelocity() * (PJ->getRadius() - .5 * deltan));
393  }
394 
395  // Compute the projection of vrel onto the normal (can be negative)
396  Mdouble vdotn = -Vec3D::dot(vrel, normal);
397 
398  //update restitution coeff
399  if (pSpecies->getIsRestitutionCoefficientConstant())
400  pSpecies->updateDissipation(PI->getMass(), PJ->getMass());
401  Mdouble a = 0, R = 0;
402 
403  // Compute normal force on particle i due to contact
404  if (pSpecies->getForceType() == ForceType::HERTZ_MINDLIN || pSpecies->getForceType() == ForceType::HERTZ_MINDLIN_DERESIEWICZ)
405  {
406  //R is twice the effective radius
407  R = 2.0 * PI->getRadius() * PJ->getRadius() / (PI->getRadius() + PJ->getRadius());
408  a = sqrt(R * deltan);
409  //pSpecies->getStiffness() stores the elastic modulus
410  Mdouble kn = 4. / 3. * pSpecies->getStiffness() * a;
411  fdotn += kn * deltan + pSpecies->getDissipation() * vdotn;
412  }
413  else
414  {
415  fdotn += pSpecies->getStiffness() * deltan + pSpecies->getDissipation() * vdotn;
416  }
417  force += normal * fdotn;
418 
419  //If tangential forces are present
420  if (pSpecies->getSlidingFrictionCoefficient() || pSpecies->getRollingFrictionCoefficient() || pSpecies->getTorsionFrictionCoefficient())
421  {
422  //call tangential spring
423  if (pSpecies->getSlidingStiffness() || pSpecies->getRollingStiffness() || pSpecies->getTorsionStiffness())
424  TangentialSpring = getTangentialSpring(PI, PJ, PJreal);
425 
426  //Compute norm of normal force
427  Mdouble norm_fn = fabs(fdotn);
428 
429  //calculate sliding friction
430  if (pSpecies->getSlidingFrictionCoefficient())
431  {
432  //Compute the tangential component of vrel
433  Vec3D vrelt = vrel + normal * vdotn;
434  //Compute norm of vrelt
435  Mdouble vdott = vrelt.getLength();
436 
437  if (pSpecies->getSlidingStiffness())
438  {
439  Vec3D* delta = &(TangentialSpring->delta);
440 
441  //Integrate the spring
443  //(*delta) += vrelt * dt - Vec3D::Dot(*delta,normal)*normal;
444  Vec3D ddelta = (vrelt - Vec3D::dot(*delta, PI->getVelocity() - PJ->getVelocity()) * normal / dist) * getTimeStep();
445  (*delta) += ddelta;
446 
447  //Calculate test force including viscous force
448  if (pSpecies->getForceType() == ForceType::HERTZ_MINDLIN)
449  {
450  //pSpecies->getSlidingStiffness() stores the elastic shear modulus
451  Mdouble kt = 8. * pSpecies->getSlidingStiffness() * a;
452  TangentialSpring->SlidingForce += -kt * ddelta;
453  forcet = TangentialSpring->SlidingForce - pSpecies->getSlidingDissipation() * vrelt;
454  }
455  else if (pSpecies->getForceType() == ForceType::HERTZ_MINDLIN_DERESIEWICZ)
456  {
457  //pSpecies->getSlidingStiffness() stores the elastic shear modulus
458  forcet = TangentialSpring->SlidingForce - pSpecies->getSlidingDissipation() * vrelt;
459  Mdouble kt = 8. * pSpecies->getSlidingStiffness() * a * std::pow(1 - Vec3D::getLength(forcet) / pSpecies->getSlidingFrictionCoefficient() / fdotn, 0.33);
460  TangentialSpring->SlidingForce += -kt * ddelta;
461  forcet = TangentialSpring->SlidingForce - pSpecies->getSlidingDissipation() * vrelt;
462  }
463  else
464  {
465  forcet = (-pSpecies->getSlidingDissipation()) * vrelt - pSpecies->getSlidingStiffness() * (*delta);
466  }
467  Mdouble forcet2 = forcet.getLengthSquared();
468 
469  //tangential forces are modelled by a spring-damper of elastisity kt and viscosity getSlidingDissipation() (sticking),
470  //but the force is limited by Coulomb friction (sliding):
471  //f_t = -getSlidingDissipation()*vrelt, if getSlidingDissipation()*vrelt<=mu_s*fdotn, f_t=mu+s*fdotn*t, else
472  double muact = (TangentialSpring->sliding) ? (pSpecies->getSlidingFrictionCoefficient()) : (pSpecies->getSlidingFrictionCoefficientStatic()); // select mu from previous sliding mode
473  if (forcet2 <= mathsFunc::square(muact * norm_fn))
474  {
475  //sticking++;
476  TangentialSpring->sliding = false;
477  }
478  else
479  {
480  //sliding++;
481  TangentialSpring->sliding = true;
483  Mdouble norm_forcet = sqrt(forcet2);
484  forcet *= pSpecies->getSlidingFrictionCoefficient() * norm_fn / norm_forcet;
485  TangentialSpring->SlidingForce = forcet + pSpecies->getSlidingDissipation() * vrelt;
486  (*delta) = TangentialSpring->SlidingForce / (-pSpecies->getSlidingStiffness());
487  }
488 
489  //Add tangential force to total force
490  force += forcet;
491 
492  }
493  else
494  { //if no tangential spring
495  //tangential forces are modelled by a damper of viscosity getSlidingDissipation() (sticking),
496  //but the force is limited by Coulomb friction (sliding):
497  //f_t = -getSlidingDissipation()*vrelt, if getSlidingDissipation()*vrelt<=mu_s*fdotn, f_t=mu+s*fdotn*t, else
498  if (vdott * pSpecies->getSlidingDissipation() <= pSpecies->getSlidingFrictionCoefficientStatic() * norm_fn)
499  { //sticking++;
500  forcet = -pSpecies->getSlidingDissipation() * vrelt;
501  }
502  else
503  { //sliding++;
504  //set force to Coulomb limit
505  forcet = -(pSpecies->getSlidingFrictionCoefficient() * norm_fn / vdott) * vrelt;
506  }
507  //Add tangential force to total force
508  force += forcet;
509  }
510  }
511 
512  //calculate rolling friction
513  if (pSpecies->getRollingFrictionCoefficient())
514  {
515  //From Luding2008, objective rolling velocity (eq 15) w/o 2.0!
516  Mdouble reducedRadiusI = PI->getRadius() - .5 * deltan;
517  Mdouble reducedRadiusJ = PJ->getRadius() - .5 * deltan;
518  Mdouble reducedRadiusIJ = 2.0 * reducedRadiusI * reducedRadiusJ / (reducedRadiusI + reducedRadiusJ);
519  Vec3D vrolling = -reducedRadiusIJ * Vec3D::cross(normal, PI->getAngularVelocity() - PJ->getAngularVelocity());
520  if (pSpecies->getRollingStiffness())
521  {
522  Vec3D* RollingSpring = &(TangentialSpring->RollingSpring);
523 
524  //Integrate the spring
525  (*RollingSpring) += vrolling * getTimeStep(); // - Vec3D::Dot(*RollingSpring,normal)*normal;
526 
527  //Calculate test force including viscous force
528  forcerolling = (-pSpecies->getRollingDissipation()) * vrolling - pSpecies->getRollingStiffness() * (*RollingSpring);
529  Mdouble forcerolling2 = forcerolling.getLengthSquared();
530 
531  //tangential forces are modelled by a spring-damper of elastisity kt and viscosity getSlidingDissipation() (sticking),
532  //but the force is limited by Coulomb friction (sliding):
533  //f_t = -getSlidingDissipation()*vrelt, if getSlidingDissipation()*vrelt<=mu_s*fdotn, f_t=mu+s*fdotn*t, else
534  double muact = (TangentialSpring->slidingRolling) ? (pSpecies->getRollingFrictionCoefficient()) : (pSpecies->getRollingFrictionCoefficientStatic()); // select mu from previous sliding mode
535  if (forcerolling2 <= mathsFunc::square(muact * norm_fn))
536  {
537  //sticking++;
538  TangentialSpring->slidingRolling = false;
539  }
540  else
541  {
542  //sliding++;
543  TangentialSpring->slidingRolling = true;
544  forcerolling *= pSpecies->getRollingFrictionCoefficient() * norm_fn / sqrt(forcerolling2);
545  (*RollingSpring) = (forcerolling + pSpecies->getRollingDissipation() * vrolling) / (-pSpecies->getRollingStiffness());
546  }
547 
548  //Add tangential force to torque
549  Vec3D Torque = reducedRadiusIJ * Vec3D::cross(normal, forcerolling);
550  PI->addTorque(Torque);
551  PJreal->addTorque(-Torque);
552  }
553  }
554 
555  //calculate torsive friction
556  if (pSpecies->getTorsionFrictionCoefficient())
557  {
558  //From Luding2008, spin velocity (eq 16) w/o 2.0!
559  Mdouble RadiusIJ = 2.0 * PI->getRadius() * PJ->getRadius() / (PI->getRadius() + PJ->getRadius());
560  Vec3D vtorsion = RadiusIJ * Vec3D::dot(normal, PI->getAngularVelocity() - PJ->getAngularVelocity()) * normal;
561  if (pSpecies->getTorsionStiffness())
562  {
563  //~ std::cout << "Error; not yet implemented" << std::endl;
564  //~ exit(-1);
565  Vec3D* TorsionSpring = &(TangentialSpring->TorsionSpring);
566 
567  //Integrate the spring
568  (*TorsionSpring) = Vec3D::dot((*TorsionSpring) + vtorsion * getTimeStep(), normal) * normal;
569 
570  //Calculate test force including viscous force
571  forcetorsion = (-pSpecies->getTorsionDissipation()) * vtorsion - pSpecies->getTorsionStiffness() * (*TorsionSpring);
572  Mdouble forcetorsion2 = forcetorsion.getLengthSquared();
573 
574  //tangential forces are modelled by a spring-damper of elastisity kt and viscosity getSlidingDissipation() (sticking),
575  //but the force is limited by Coulomb friction (sliding):
576  //f_t = -getSlidingDissipation()*vrelt, if getSlidingDissipation()*vrelt<=mu_s*fdotn, f_t=mu+s*fdotn*t, else
577  double muact = (TangentialSpring->slidingTorsion) ? (pSpecies->getTorsionFrictionCoefficient()) : (pSpecies->getTorsionFrictionCoefficientStatic()); // select mu from previous sliding mode
578  if (forcetorsion2 <= mathsFunc::square(muact * norm_fn))
579  {
580  //sticking++;
581  TangentialSpring->slidingTorsion = false;
582  }
583  else
584  {
585  //sliding++;
586  TangentialSpring->slidingTorsion = true;
587  //~ std::cout << "sliding " << std::endl;
588  forcetorsion *= pSpecies->getTorsionFrictionCoefficient() * norm_fn / sqrt(forcetorsion2);
589  (*TorsionSpring) = (forcetorsion + pSpecies->getTorsionDissipation() * vtorsion) / (-pSpecies->getTorsionStiffness());
590  }
591 
592  //Add tangential force to torque
593  Vec3D torque = RadiusIJ * forcetorsion;
594  PI->addTorque(torque);
595  PJreal->addTorque(-torque);
596 
597  }
598  }
599  } //end if tangential forces
600 
602  //Make force Hertzian (note: deltan is normalized by the normal distanct of two particles in contact, as in Silbert)
603  if (pSpecies->getForceType() == ForceType::HERTZ)
604  force *= sqrt(deltan / (PI->getRadius() + PJ->getRadius()));
605 
606  }
607  else
608  { //end if contact forces
609  force += fdotn * normal;
610  }
611 
612  //Add forces to total force
613  //PI->addForce(force);
614  //if (!isPeriodic)
615  // PJreal->addForce(-force);
616 
617  // Add torque due to tangential forces: t = Vec3D::Cross(l,f), l=dist*Wall.normal
618  //if (pSpecies->getSlidingFrictionCoefficient())
619  //{
620  // Vec3D Vec3D::Cross = Vec3D::Cross(normal, force);
621  // PI->addTorque(-Vec3D::Cross * (PI->getRadius() - .5 * deltan));
622  // if (!isPeriodic)
623  // PJreal->addTorque(-Vec3D::Cross * (PJ->getRadius() - .5 * deltan));
624  //}
625  //BEGIN: CHANGE
626  if (InPeriodicBox(PI) != InPeriodicBox(PJreal))
627  {
628  if (InPeriodicBox(PI))
629  {
631  {
632  //~ if (PJreal->isFixed()) cout << "PJreal" << endl;
633  PJreal->addForce(-force);
634  }
635  }
636  else
637  {
639  {
640  //~ if (Particles[PI].isFixed()) cout << "PI" << endl;
641  PI->addForce(force);
642  }
643  }
644  }
645  else
646  {
647  PI->addForce(force);
648  PJreal->addForce(-force);
649  }
650  //END: CHANGE
651 
652  // Add torque due to tangential forces: t = Vec3D::Cross(l,f), l=dist*Wall.normal
653  //if (pSpecies->getSlidingFrictionCoefficient()) {
654  // Vec3D Vec3D::Cross = Vec3D::Cross(normal, force);
655  // PI ->add_Torque(-Vec3D::Cross * (PI->getRadius() - .5 * deltan));
656  // PJreal->add_Torque(-Vec3D::Cross * (PJ->getRadius() - .5 * deltan));
657  //}
658  //BEGIN: CHANGE
659  if (pSpecies->getSlidingFrictionCoefficient())
660  {
661  Vec3D cross = Vec3D::cross(normal, force);
662  if (InPeriodicBox(PI) != InPeriodicBox(PJreal))
663  {
664  if (InPeriodicBox(PI))
665  {
667  {
668  PJreal->addTorque(-cross * (PJreal->getRadius() - .5 * deltan));
669  }
670  }
671  else
672  {
674  {
675  PI->addTorque(-cross * (PI->getRadius() - .5 * deltan));
676  }
677  }
678  }
679  else
680  {
681  PI->addTorque(-cross * (PI->getRadius() - .5 * deltan));
682  PJreal->addTorque(-cross * (PJreal->getRadius() - .5 * deltan));
683  }
684  }
685  //END:CHANGE
686 
687  // output for ene and stat files:
688  if (eneFile.getSaveCurrentTimeStep())
689  {
690  if (!isPeriodic)
691  addElasticEnergy(0.5 * (pSpecies->getStiffness() * mathsFunc::square(deltan) + (TangentialSpring ? (pSpecies->getSlidingStiffness() * TangentialSpring->delta.getLengthSquared() + pSpecies->getRollingStiffness() * TangentialSpring->RollingSpring.getLengthSquared() + pSpecies->getTorsionStiffness() * TangentialSpring->TorsionSpring.getLengthSquared()) : 0.0)));
692  else
693  addElasticEnergy(0.25 * (pSpecies->getStiffness() * mathsFunc::square(deltan) + (TangentialSpring ? (pSpecies->getSlidingStiffness() * TangentialSpring->delta.getLengthSquared() + pSpecies->getRollingStiffness() * TangentialSpring->RollingSpring.getLengthSquared() + pSpecies->getTorsionStiffness() * TangentialSpring->TorsionSpring.getLengthSquared()) : 0.0)));
694  }
695  if (fStatFile.getSaveCurrentTimeStep() || statFile.getSaveCurrentTimeStep() || getDoCGAlways())
696  {
697 
698  Mdouble fdott = forcet.getLength();
699  Mdouble deltat_norm = TangentialSpring ? (-TangentialSpring->delta.getLength()) : 0.0;
700 
703  Vec3D centre = 0.5 * (PI->getPosition() + PJ->getPosition());
704 
708 
709  if (!PI->isFixed())
710  {
711  if (statFile.getSaveCurrentTimeStep() || getDoCGAlways())
712  gatherContactStatistics(PJreal->getIndex(), PI->getIndex(), centre, deltan, deltat_norm, fdotn, fdott, -normal, -(fdott ? forcet / fdott : forcet));
713  if (fStatFile.getSaveCurrentTimeStep())
714  fStatFile.getFstream() << getTime() << " " << PJreal->getIndex() << " " << PI->getIndex() << " " << centre << " " << radii_sum - dist << " " << deltat_norm << " " << fdotn << " " << fdott << " " << -normal << " " << -(fdott ? forcet / fdott : forcet) << std::endl;
715  }
716  if (!PJreal->isFixed() && !isPeriodic)
717  {
718  if (statFile.getSaveCurrentTimeStep() || getDoCGAlways())
719  gatherContactStatistics(PI->getIndex(), PJreal->getIndex(), centre, deltan, deltat_norm, fdotn, fdott, normal, (fdott ? forcet / fdott : forcet));
720  if (fStatFile.getSaveCurrentTimeStep())
721  fStatFile.getFstream() << getTime() << " " << PI->getIndex() << " " << PJreal->getIndex() << " " << centre << " " << radii_sum - dist << " " << deltat_norm << " " << fdotn << " " << fdott << " " << normal << " " << (fdott ? forcet / fdott : forcet) << std::endl;
722  }
723  }
724  } // end if particle i and j are overlapping
725  }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
@ R
Definition: StatisticsVector.h:21
virtual const Vec3D & getAngularVelocity() const
Returns the angular velocity of this interactable.
Definition: BaseInteractable.cc:319
void addTorque(const Vec3D &addTorque)
Adds an amount to the torque on this BaseInteractable.
Definition: BaseInteractable.cc:110
virtual const Vec3D & getVelocity() const
Returns the velocity of this interactable.
Definition: BaseInteractable.cc:307
void addForce(const Vec3D &addForce)
Adds an amount to the force on this BaseInteractable.
Definition: BaseInteractable.cc:94
unsigned int getIndSpecies() const
Returns the index of the species associated with the interactable object.
Definition: BaseInteractable.h:67
unsigned int getIndex() const
Returns the index of the object in the handler.
Definition: BaseObject.h:97
Mdouble getRadius() const
Returns the particle's radius.
Definition: BaseParticle.h:331
Mdouble getMass() const
Returns the particle's mass.
Definition: BaseParticle.h:305
BaseSpecies is the class from which all other species are derived.
Definition: BaseSpecies.h:29
double get_PeriodicBoxLength()
Definition: ChuteWithPeriodicInflow.h:850
File eneFile
An instance of class File to handle in- and output into a .ene file.
Definition: DPMBase.h:1494
File fStatFile
An instance of class File to handle in- and output into a .fstat file.
Definition: DPMBase.h:1489
Mdouble getTimeStep() const
Returns the simulation time step.
Definition: DPMBase.cc:1241
Mdouble getTime() const
Returns the current simulation time.
Definition: DPMBase.cc:799
void gatherContactStatistics()
Definition: DPMBase.cc:1887
File statFile
An instance of class File to handle in- and output into a .stat file.
Definition: DPMBase.h:1504
std::fstream & getFstream()
Allows to access the member variable File::fstream_.
Definition: File.cc:131
BaseParticle * getLargestParticle() const
Returns the pointer of the largest particle in the particle handler. When mercury is running in paral...
Definition: ParticleHandler.cc:528
std::enable_if<!std::is_pointer< typename U::MixedSpeciesType >::value, typename U::MixedSpeciesType * >::type getMixedObject(const U *S, const U *T)
Definition: SpeciesHandler.h:52
Contains material and contact force properties.
Definition: Species.h:14
static Mdouble getLengthSquared(const Vec3D &a)
Calculates the squared length of a Vec3D: .
Definition: Kernel/Math/Vector.h:324
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:143
static Mdouble getLength(const Vec3D &a)
Calculates the length of a Vec3D: .
Definition: Vector.cc:350
static Mdouble getDistanceSquared(const Vec3D &a, const Vec3D &b)
Calculates the squared distance between two Vec3D: .
Definition: Kernel/Math/Vector.h:303
void setZero()
Sets all elements to zero.
Definition: Vector.cc:23
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:56
Mdouble getLength() const
Calculates the length of this Vec3D: .
Definition: Vector.cc:339
#define max(a, b)
Definition: datatypes.h:23
const Scalar * a
Definition: level2_cplx_impl.h:32
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
int delta
Definition: MultiOpt.py:96
void normal(const Vector< double > &x, Vector< double > &normal)
Definition: free_surface_rotation.cc:65
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
T square(const T val)
squares a number
Definition: ExtendedMath.h:86
void cross(const Vector< double > &A, const Vector< double > &B, Vector< double > &C)
Definition: oomph-lib/src/generic/Vector.h:319

References a, BaseInteractable::addForce(), BaseInteractable::addTorque(), Vec3D::cross(), oomph::VectorHelpers::cross(), MultiOpt::delta, Vec3D::dot(), DPMBase::eneFile, boost::multiprecision::fabs(), DPMBase::fStatFile, DPMBase::gatherContactStatistics(), get_PeriodicBoxLength(), BaseInteractable::getAngularVelocity(), Vec3D::getDistanceSquared(), File::getFstream(), BaseObject::getIndex(), BaseInteractable::getIndSpecies(), ParticleHandler::getLargestParticle(), Vec3D::getLength(), Vec3D::getLengthSquared(), BaseParticle::getMass(), SpeciesHandler::getMixedObject(), BaseParticle::getPeriodicFromParticle(), BaseInteractable::getPosition(), BaseParticle::getRadius(), DPMBase::getTime(), DPMBase::getTimeStep(), BaseInteractable::getVelocity(), InPeriodicBox(), BaseParticle::isFixed(), max, WallFunction::normal(), DPMBase::particleHandler, Eigen::bfloat16_impl::pow(), R, Vec3D::setZero(), DPMBase::speciesHandler, sqrt(), mathsFunc::square(), DPMBase::statFile, and Vec3D::X.

◆ ExtendInWidth()

void ChuteWithPeriodicInflow::ExtendInWidth ( int  numRepetitionsInWidth)
inline
97  {
99  // creates bottom outside periodic chute of species 1
100  double widthPeriodicChute = getYMax();
102  particleHandler.setStorageCapacity(N * numRepetitionsInWidth);
103  for (int j = 1; j < numRepetitionsInWidth; j++)
104  {
105  for (int i = 0; i < N; i++)
106  {
107  //if (particleHandler.getObject(i)->isFixed()) {
109  //Particles.back()->get_IndSpecies()=1;
110  particleHandler.getLastObject()->move(Vec3D(0.0, widthPeriodicChute * j, 0.0));
111  //}
112  }
113  }
114  setYMax(numRepetitionsInWidth * widthPeriodicChute);
115  perw->set(Vec3D(0, 1, 0), getYMin(), getYMax());
116  //setHGridNumberOfBucketsToPower(particleHandler.getNumberOfObjects());
117  }
Mdouble getYMin() const
If the length of the problem domain in y-direction is YMax - YMin, then getYMin() returns YMin.
Definition: DPMBase.h:616
void setYMax(Mdouble newYMax)
Sets the value of YMax, the upper bound of the problem domain in the y-direction.
Definition: DPMBase.cc:1182
Mdouble getYMax() const
If the length of the problem domain in y-direction is YMax - YMin, then getYMax() returns XMax.
Definition: DPMBase.h:622
void set(Vec3D normal, Mdouble distanceLeft, Mdouble distanceRight)
Defines a PeriodicBoundary by its normal and positions.
Definition: PeriodicBoundary.cc:63

References ParticleHandler::addObject(), DPMBase::boundaryHandler, BaseParticle::copy(), BaseHandler< T >::getLastObject(), ParticleHandler::getNumberOfObjects(), BaseHandler< T >::getObject(), DPMBase::getYMax(), DPMBase::getYMin(), i, j, BaseInteractable::move(), N, DPMBase::particleHandler, PeriodicBoundary::set(), BaseHandler< T >::setStorageCapacity(), and DPMBase::setYMax().

Referenced by ChuteWithPeriodicInflow().

◆ get_PeriodicBoxLength()

double ChuteWithPeriodicInflow::get_PeriodicBoxLength ( )
inline
851  {
852  return PeriodicBoxLength;
853  }
double PeriodicBoxLength
stores the length of the periodic box
Definition: ChuteWithPeriodicInflow.h:876

References PeriodicBoxLength.

Referenced by computeInternalForces().

◆ get_PeriodicBoxNSpecies()

int ChuteWithPeriodicInflow::get_PeriodicBoxNSpecies ( )
inline
859  {
860  return PeriodicBoxNSpecies;
861  }
int PeriodicBoxNSpecies
stores the number of species in the periodic box
Definition: ChuteWithPeriodicInflow.h:878

References PeriodicBoxNSpecies.

Referenced by integrateBeforeForceComputation(), and loadPeriodicBox().

◆ getInfo()

double ChuteWithPeriodicInflow::getInfo ( const BaseParticle p) const
inlineoverridevirtual

A virtual function that returns some user-specified information about a particle.

Returns
double

Reimplemented from DPMBase.

38  {
39  return P.getIndSpecies();
40  }
double P
Uniform pressure.
Definition: TwenteMeshGluing.cpp:77

References Global_Physical_Variables::P.

◆ InPeriodicBox()

bool ChuteWithPeriodicInflow::InPeriodicBox ( BaseParticle P)
inline
868  {
869  return P->getIndSpecies() < PeriodicBoxNSpecies;
870  }

References Global_Physical_Variables::P, and PeriodicBoxNSpecies.

Referenced by Check_and_Duplicate_Periodic_Particle(), and computeInternalForces().

◆ integrateBeforeForceComputation()

void ChuteWithPeriodicInflow::integrateBeforeForceComputation ( )
inlinevirtual

Update particles' and walls' positions and velocities before force computation.

This is where the integration is done, at the moment it is velocity Verlet integration and is done before the forces are computed. See http://en.wikipedia.org/wiki/Verlet_integration#Velocity_Verlet

Performs integration - i.e. updating particle's positions, velocities and accelerations - for all particles and walls within the system (i.e. in the particleHandler and wallHandler). Integration is performed using the BaseParticle::integrateBeforeForceComputation() function.

The velocity Verlet algorithm requires us to integrate twice each time step: both before and after the force computation. This method is therefore used in conjunction with DPMBase::integrateAfterForceComputation(). See http://en.wikipedia.org/wiki/Verlet_integration#Velocity_Verlet for details.

Reimplemented from DPMBase.

156  {
157  for (std::vector<BaseParticle*>::iterator it = particleHandler.begin(); it != particleHandler.end(); ++it)
158  {
159  (*it)->integrateBeforeForceComputation(getTimeStep());
160 
161  // This shifts particles that moved through periodic walls
163  PeriodicBoundary* perw1 = static_cast<PeriodicBoundary*>(boundaryHandler.getObject(1));
164  if ((*it)->getIndSpecies() == 0 && perw->getDistance(**it) < 0)
165  {
166  //if (iP->getIndSpecies()==0 && static_cast<PeriodicWall*>(boundaryHandler.getObject(0))->getDistance(*iP)<0) {
167  if (!perw->isClosestToLeftBoundary())
168  {
171  particleHandler.addObject((*it)->copy());
174  perw->shiftPosition((*it));
176  {
177  hGridRemoveParticle((*it));
178  hGridUpdateParticle((*it));
179  }
180  }
181  else
182  {
183  static int count = -1;
184  count++;
185  if (!(count % dataFile.getSaveCount()))
186 
187  std::cout << "Warning: Particle " << (*it)->getIndex() << " is left of xmin: x="
188  << (*it)->getPosition().X << ", v_x=" << (*it)->getVelocity().X << std::endl;
189  }
190  }
191  if (getIsPeriodic() > 1 && perw1->getDistance(**it) < 0)
192  {
193  perw1->shiftPosition(*it);
195  {
196  hGridRemoveParticle(*it);
197  hGridUpdateParticle(*it);
198  }
199  }
200  }
201 
202  }
const std::vector< T * >::const_iterator begin() const
Gets the begin of the const_iterator over all Object in this BaseHandler.
Definition: BaseHandler.h:698
unsigned int getStorageCapacity() const
Gets the storage capacity of this BaseHandler.
Definition: BaseHandler.h:670
const std::vector< T * >::const_iterator end() const
Gets the end of the const_iterator over all BaseBoundary in this BaseHandler.
Definition: BaseHandler.h:712
int get_PeriodicBoxNSpecies()
Definition: ChuteWithPeriodicInflow.h:858
File dataFile
An instance of class File to handle in- and output into a .data file.
Definition: DPMBase.h:1484
unsigned int getSaveCount() const
Gets File::saveCount_.
Definition: File.cc:233
void hGridUpdateParticle(BaseParticle *obj) override
Updates the cell (not the level) of a BaseParticle.
Definition: Mercury3D.cc:340
void hGridRemoveParticle(BaseParticle *obj) override
Removes a BaseParticle from the HGrid.
Definition: Mercury3D.cc:401
bool getHGridUpdateEachTimeStep() const final
Gets whether or not the HGrid is updated every time step.
Definition: MercuryBase.cc:163
Definition: ParticleSpecies.h:16
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:254

References ParticleHandler::addObject(), BaseHandler< T >::begin(), DPMBase::boundaryHandler, DPMBase::dataFile, BaseHandler< T >::end(), get_PeriodicBoxNSpecies(), PeriodicBoundary::getDistance(), MercuryBase::getHGridUpdateEachTimeStep(), BaseInteractable::getIndSpecies(), Chute::getIsPeriodic(), BaseHandler< T >::getLastObject(), ParticleHandler::getNumberOfObjects(), BaseHandler< T >::getObject(), File::getSaveCount(), BaseHandler< T >::getStorageCapacity(), DPMBase::getTimeStep(), Mercury3D::hGridRemoveParticle(), Mercury3D::hGridUpdateParticle(), PeriodicBoundary::isClosestToLeftBoundary(), DPMBase::particleHandler, BaseParticle::setSpecies(), BaseHandler< T >::setStorageCapacity(), PeriodicBoundary::shiftPosition(), and DPMBase::speciesHandler.

◆ loadPeriodicBox()

void ChuteWithPeriodicInflow::loadPeriodicBox ( std::string const  restart_file)
inline

loads periodic chute data from restart file

44  {
45  FillChute = false;
46 
47  std::cout << "constructor" << std::endl;
48 
49  //load particles and chute setup into periodic chute
50  setName(restart_file.c_str());
52 
55 
56  //we don't want to treat the data as restarted, i.e. we start the program at time=0 and create new output files
57  setRestarted(false);
58 
59  //keep file name but create files in the local directory, i.e. remove folder
60  size_t found = restart_file.find_last_of("/\\");
61  setName(restart_file.substr(found + 1).c_str());
62 
63  // adds new species for the flow particles (0 for the periodic inflow particles)
64  for (int i = 0; i < get_PeriodicBoxNSpecies(); i++)
65  addSpecies(Species[i]);
66 
67  setName("ChuteWithPeriodicInflow");
68  }
virtual unsigned int getNumberOfObjects() const
Gets the number of real Object in this BaseHandler. (i.e. no mpi or periodic particles)
Definition: BaseHandler.h:656
void set_PeriodicBoxNSpecies(int new_)
Definition: ChuteWithPeriodicInflow.h:862
void set_PeriodicBoxLength(double new_)
Definition: ChuteWithPeriodicInflow.h:854
void setName(const std::string &name)
Allows to set the name of all the files (ene, data, fstat, restart, stat)
Definition: DPMBase.cc:400
void setRestarted(bool newRestartedFlag)
Allows to set the flag stating if the simulation is to be restarted or not.
Definition: DPMBase.cc:1492
bool readRestartFile(ReadOptions opt=ReadOptions::ReadAll)
Reads all the particle data corresponding to a given, existing . restart file (for more details regar...
Definition: DPMBase.cc:3043
bool found
Definition: MergeRestartFiles.py:24

References FillChute, MergeRestartFiles::found, get_PeriodicBoxNSpecies(), BaseHandler< T >::getNumberOfObjects(), DPMBase::getXMax(), i, DPMBase::readRestartFile(), set_PeriodicBoxLength(), set_PeriodicBoxNSpecies(), DPMBase::setName(), DPMBase::setRestarted(), and DPMBase::speciesHandler.

Referenced by ChuteWithPeriodicInflow().

◆ outputXBallsDataParticlee()

void ChuteWithPeriodicInflow::outputXBallsDataParticlee ( const unsigned int  i,
const unsigned int  format,
std::ostream &  os 
) const
inline
829  {
830  os
831  << particleHandler.getObject(i)->getPosition().X << " "
832  << particleHandler.getObject(i)->getPosition().Y << " "
833  << particleHandler.getObject(i)->getPosition().Z << " "
834  << particleHandler.getObject(i)->getVelocity().X << " "
835  << particleHandler.getObject(i)->getVelocity().Y << " "
836  << particleHandler.getObject(i)->getVelocity().Z << " "
837  //~ << (particleHandler.getObject(i)->isFixed()?0.0:particleHandler.getObject(i)->Force.X) << " "
838  //~ << (particleHandler.getObject(i)->isFixed()?0.0:particleHandler.getObject(i)->Force.Y) << " "
839  //~ << (particleHandler.getObject(i)->isFixed()?0.0:particleHandler.getObject(i)->Force.Z) << " "
840  << particleHandler.getObject(i)->getRadius() << " "
842  << particleHandler.getObject(i)->getOrientation().Y << " "
843  << particleHandler.getObject(i)->getOrientation().Z << " "
847  << particleHandler.getObject(i)->getIndSpecies() << std::endl;
848  }
const Quaternion & getOrientation() const
Returns the orientation of this BaseInteractable.
Definition: BaseInteractable.h:209
Mdouble Y
Definition: Kernel/Math/Vector.h:45
Mdouble Z
Definition: Kernel/Math/Vector.h:45

References BaseInteractable::getAngularVelocity(), BaseInteractable::getIndSpecies(), BaseHandler< T >::getObject(), BaseInteractable::getOrientation(), BaseInteractable::getPosition(), BaseParticle::getRadius(), BaseInteractable::getVelocity(), i, DPMBase::particleHandler, Vec3D::X, Vec3D::Y, and Vec3D::Z.

◆ printTime()

void ChuteWithPeriodicInflow::printTime ( ) const
inlineoverridevirtual

add some particular output

Reimplemented from DPMBase.

206  {
207  int counter = 0;
208 
209  double speed1 = 0;
210  double speed0 = 0;
211  double n1 = 0;
212  double n0 = 0;
213  for (int i = 0; i < particleHandler.getNumberOfObjects(); i++)
214  {
216  counter++;
218  {
220  {
221  speed0 += particleHandler.getObject(i)->getVelocity().X;
222  n0++;
223  }
224  else
225  {
226  speed1 += particleHandler.getObject(i)->getVelocity().X;
227  n1++;
228  }
229  }
230  }
231  speed0 /= n0;
232  speed1 /= n1;
233 
234  std::cout << "t=" << std::setprecision(5) << std::left << std::setw(6) << getTime()
235  << ", n=" << n1 << "(inflow:" << n0 << ")"
236  << ", speed= " << speed1 << "(" << speed0 << ")"
237  << std::endl;
238  //~ cout
239  //~ << " A " << PeriodicBoxLength
240  //~ << " A " << PeriodicBoxNSpecies
241  //~ << " A " << FillChute << endl;
242  }

References BaseInteractable::getIndSpecies(), ParticleHandler::getNumberOfObjects(), BaseHandler< T >::getObject(), DPMBase::getTime(), BaseInteractable::getVelocity(), i, BaseParticle::isFixed(), DPMBase::particleHandler, and Vec3D::X.

◆ set_PeriodicBoxLength()

void ChuteWithPeriodicInflow::set_PeriodicBoxLength ( double  new_)
inline
855  {
856  PeriodicBoxLength = new_;
857  }

References PeriodicBoxLength.

Referenced by loadPeriodicBox().

◆ set_PeriodicBoxNSpecies()

void ChuteWithPeriodicInflow::set_PeriodicBoxNSpecies ( int  new_)
inline
863  {
864  PeriodicBoxNSpecies = new_;
865  }

References PeriodicBoxNSpecies.

Referenced by loadPeriodicBox().

◆ setupInitialConditions()

void ChuteWithPeriodicInflow::setupInitialConditions ( )
inlineoverridevirtual

Do not add bottom.

Reimplemented from DPMBase.

152  {
153  }

◆ writeXBallsScript()

void ChuteWithPeriodicInflow::writeXBallsScript ( ) const
inlinevirtual

This writes a script which can be used to load the xballs problem to display the data just generated.

Todo:
Implement or make pure virtual

First put in all the script lines. All these lines do is move you to the correct directory from any location

Reimplemented from DPMBase.

756  {
757 
758  std::stringstream file_name;
759  std::ofstream script_file;
760  file_name << getName() << ".disp";
761  script_file.open((file_name.str()).c_str());
762 
764  script_file << "#!/bin/bash" << std::endl;
765  script_file << "x=$(echo $0 | cut -c2-)" << std::endl;
766  script_file << "file=$PWD$x" << std::endl;
767  script_file << "dirname=`dirname \"$file\"`" << std::endl;
768  script_file << "cd $dirname" << std::endl;
769 
770  double scale;
771  int format;
772 
773  int width = 1570, height = 860;
774  double ratio = getSystemDimensions() < 3 ? (getXMax() - getXMin()) / (getYMax() - getYMin()) : (getXMax() - getXMin()) / (getZMax() - getZMin());
775  if (ratio > width / height)
776  height = width / ratio;
777  else
778  width = height * ratio;
779  //~ cout << ratio << "r" << endl;
780  if (getSystemDimensions() < 3)
781  { // dim = 1 or 2
782  format = 8;
783  if (getXBallsScale() < 0)
784  {
785  scale = 1.0 / std::max(getYMax() - getYMin(), getXMax() - getXMin());
786  }
787  else
788  {
789  scale = getXBallsScale();
790  }
791  }
792  else
793  { //dim==3
794  format = 14;
795  if (getXBallsScale() < 0)
796  {
797  scale = width / 480 / (getXMax() - getXMin());
798  }
799 
800  else
801  {
802  scale = getXBallsScale();
803  }
804 
805  }
806 
807  script_file << "../xballs -format " << format
808  << " -f " << dataFile.getFullName()
809  << " -s " << scale
810  << " -w " << width + 140
811  << " -h " << height + 140
812  << " -cmode " << getXBallsColourMode()
813  << " -cmax -scala 4 -sort -oh 50 "
815  << " $*";
816  if (getXBallsVectorScale() > -1)
817  {
818  script_file << " -vscale " << getXBallsVectorScale();
819  }
820  script_file.close();
821 
822  //This line changes the file permision and gives the owner (i.e. you) read, write and execute permission to the file
823 #ifdef UNIX
824  chmod((file_name.str().c_str()),S_IRWXU);
825 #endif
826  }
int getXBallsColourMode() const
Get the xballs colour mode (CMode).
Definition: DPMBase.cc:1301
Mdouble getXMin() const
If the length of the problem domain in x-direction is XMax - XMin, then getXMin() returns XMin.
Definition: DPMBase.h:603
double getXBallsVectorScale() const
Returns the scale of vectors used in xballs.
Definition: DPMBase.cc:1321
const std::string & getName() const
Returns the name of the file. Does not allow to change it though.
Definition: DPMBase.cc:377
std::string getXBallsAdditionalArguments() const
Returns the additional arguments for xballs.
Definition: DPMBase.cc:1346
double getXBallsScale() const
Returns the scale of the view in xballs.
Definition: DPMBase.cc:1363
unsigned int getSystemDimensions() const
Returns the system dimensionality.
Definition: DPMBase.cc:1421
Mdouble getZMax() const
If the length of the problem domain in z-direction is ZMax - ZMin, then getZMax() returns ZMax.
Definition: DPMBase.h:634
Mdouble getZMin() const
If the length of the problem domain in z-direction is ZMax - ZMin, then getZMin() returns ZMin.
Definition: DPMBase.h:628
const std::string getFullName() const
Also allows to access the file name, however with additional information which is the file counter,...
Definition: File.cc:148
double height(const double &x)
Height of domain.
Definition: simple_spine_channel.cc:429
string file_name
Definition: Particles2023AnalysisHung.py:321
std::string format(const std::string &str, const std::vector< std::string > &find, const std::vector< std::string > &replace)
Definition: openglsupport.cpp:217

References DPMBase::dataFile, Particles2023AnalysisHung::file_name, format(), File::getFullName(), DPMBase::getName(), DPMBase::getSystemDimensions(), DPMBase::getXBallsAdditionalArguments(), DPMBase::getXBallsColourMode(), DPMBase::getXBallsScale(), DPMBase::getXBallsVectorScale(), DPMBase::getXMax(), DPMBase::getXMin(), DPMBase::getYMax(), DPMBase::getYMin(), DPMBase::getZMax(), DPMBase::getZMin(), Global_Physical_Variables::height(), and max.

Member Data Documentation

◆ FillChute

bool ChuteWithPeriodicInflow::FillChute
private

◆ inflowParticle_

◆ PeriodicBoxLength

double ChuteWithPeriodicInflow::PeriodicBoxLength
private

stores the length of the periodic box

Referenced by get_PeriodicBoxLength(), and set_PeriodicBoxLength().

◆ PeriodicBoxNSpecies

int ChuteWithPeriodicInflow::PeriodicBoxNSpecies
private

stores the number of species in the periodic box

Referenced by get_PeriodicBoxNSpecies(), InPeriodicBox(), and set_PeriodicBoxNSpecies().


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