DrumRot Class Reference
+ Inheritance diagram for DrumRot:

Public Member Functions

void setupInitialConditions () override
 This function allows to set the initial conditions for our problem to be solved, by default particle locations are randomly set. Remember particle properties must also be defined here. More...
 
void actionsBeforeTimeStep () override
 A virtual function which allows to define operations to be executed before the new time step. More...
 
void actionsOnRestart () override
 A virtual function where the users can add extra code which is executed only when the code is restarted. More...
 
void setDrumRadius (double radius)
 
void setRevolutionSpeed (double rpm)
 
void setFractionalPolydispersity (double fpd)
 
void setDrumFillFraction (double dff)
 
void setFrictionCoeff (double pwf, double ppf)
 
void setCOR (double drumCOR, double COR1)
 
void setSlidingFriction (double drum, double f1)
 
void setRollingFriction (double drum, double f1)
 
void setTorsionFriction (double drum, double f1)
 
double getDrumRadius ()
 
- 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< int > get1DParametersFromRunNumber (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< int > get2DParametersFromRunNumber (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< int > get3DParametersFromRunNumber (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...
 
virtual void writeXBallsScript () const
 This writes a script which can be used to load the xballs problem to display the data just generated. More...
 
virtual Mdouble getInfo (const BaseParticle &P) const
 A virtual function that returns some user-specified information about a particle. 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
 
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 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 getMin () const
 
Vec3D getMax () const
 
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)
 
int getNumberOfOMPThreads () const
 
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 allowd to be negaitve 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 () const
 Returns the global gravitational potential energy stored in the system. More...
 
Mdouble getRotationalEnergy () const
 JMFT 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< unsigned > getNumberOfDomains ()
 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)
 
void writePythonFileForVTKVisualisation () const
 
void setWritePythonFileForVTKVisualisation (bool forceWritePythonFileForVTKVisualisation)
 
bool getWritePythonFileForVTKVisualisation () const
 
WallVTKWritergetWallVTKWriter ()
 

Private Attributes

double radiusS1
 
double rhoS1
 
double massS1
 
double CORDrum
 
double CORS1
 
double tc
 
double sizeRatio
 
double densityRatio
 
double drumFillFraction
 
double volumeFraction
 
double particleWallFriction
 
double particleParticleFriction
 
int numS1
 
int numS1ToBeInserted
 
double drumRadius
 
double revolutionsPerSecond
 
double fractionalPolydispersity
 
double slidingFrictionDrum
 
double slidingFriction1
 
double rollingFrictionDrum
 
double rollingFriction1
 
double torsionFrictionDrum
 
double torsionFriction1
 
int step
 
double checkTime
 

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...
 
- 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...
 
- 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 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 integrateBeforeForceComputation ()
 Update particles' and walls' positions and velocities before force computation. More...
 
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 void printTime () const
 Displays the current simulation time and the maximum simulation duration. 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 setSoftStop ()
 function for setting sigaction constructor. More...
 
- Static Protected Member Functions inherited from DPMBase
static void signalHandler (int signal)
 signal handler function. More...
 

Member Function Documentation

◆ actionsBeforeTimeStep()

void DrumRot::actionsBeforeTimeStep ( )
inlineoverridevirtual

A virtual function which allows to define operations to be executed before the new time step.

no implementation but can be overidden in its derived classes.

Reimplemented from DPMBase.

174  {
175  //wallHandler.getObject(0)->setOrientation(Vec3D(0.0,1.0,0.0));
176  //wallHandler.getObject(1)->setOrientation(Vec3D(0.0,-1.0,0.0));
177  //wallHandler.getObject(2)->setOrientation(Vec3D(0.0,1.0,0.0));
178 
179  if (step==2)
180  {
181  if (getTime() > checkTime)
182  {
183  logger(INFO, "Current KE %", getKineticEnergy());
184  if (getKineticEnergy() < (10.0))
185  {
186  step = 3;
187  setTime(0.0);
188  setSaveCount(2000);
189  double drumStartTime = getTime();
190  logger(INFO, "\n \n \n"
191  "Start the drum rotation, rpm = %"
192  "--------------------------"
193  "\n\n\n", revolutionsPerSecond * 60.0);
194  // rotate the drum
196  Vec3D(0.0, revolutionsPerSecond * constants::pi * 2.0, 0.0));
198  Vec3D(0.0, revolutionsPerSecond * constants::pi * 2.0, 0.0));
200  Vec3D(0.0, revolutionsPerSecond * constants::pi * 2.0, 0.0));
201 
202  //wallHandler.getObject(0)->setOrientation(Vec3D(0.0,1.0,0.0));
203  //wallHandler.getObject(1)->setOrientation(Vec3D(0.0,-1.0,0.0));
204  //wallHandler.getObject(2)->setOrientation(Vec3D(0.0,1.0,0.0));
205  }
206  else
207  {
208  checkTime = getTime() + 0.1;
209  }
210  }
211  }
212  }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ INFO
T * getObject(const unsigned int id)
Gets a pointer to the Object at the specified index in the BaseHandler.
Definition: BaseHandler.h:613
void setAngularVelocity(const Vec3D &angularVelocity)
set the angular velocity of the BaseInteractble.
Definition: BaseInteractable.cc:360
void setSaveCount(unsigned int saveCount)
Sets File::saveCount_ for all files (ene, data, fstat, restart, stat)
Definition: DPMBase.cc:408
Mdouble getTime() const
Returns the current simulation time.
Definition: DPMBase.cc:808
Mdouble getKineticEnergy() const
Returns the global kinetic energy stored in the system.
Definition: DPMBase.cc:1544
WallHandler wallHandler
An object of the class WallHandler. Contains pointers to all the walls created.
Definition: DPMBase.h:1447
void setTime(Mdouble time)
Sets a new value for the current simulation time.
Definition: DPMBase.cc:836
double revolutionsPerSecond
Definition: MonodispersedDrum.cpp:309
int step
Definition: MonodispersedDrum.cpp:321
double checkTime
Definition: MonodispersedDrum.cpp:322
Definition: Vector.h:51
const Mdouble pi
Definition: ExtendedMath.h:45

References checkTime, DPMBase::getKineticEnergy(), BaseHandler< T >::getObject(), DPMBase::getTime(), INFO, logger, constants::pi, revolutionsPerSecond, BaseInteractable::setAngularVelocity(), DPMBase::setSaveCount(), DPMBase::setTime(), step, and DPMBase::wallHandler.

◆ actionsOnRestart()

void DrumRot::actionsOnRestart ( )
inlineoverridevirtual

A virtual function where the users can add extra code which is executed only when the code is restarted.

no implementation but can be overidden in its derived classes.

Reimplemented from DPMBase.

215  {
216 
217  //setTimeMax(200.0);
218  setSaveCount(2000);
219 
220  logger(INFO, "\n \n \n"
221  "set drum rotation speed (rpm) to %"
222  "--------------------------"
223  "\n\n\n", revolutionsPerSecond * 60.0);
224 
226  Vec3D(0.0, revolutionsPerSecond * constants::pi * 2.0, 0.0)); //rad/s
229  }

References BaseHandler< T >::getObject(), INFO, logger, constants::pi, revolutionsPerSecond, BaseInteractable::setAngularVelocity(), DPMBase::setSaveCount(), and DPMBase::wallHandler.

◆ getDrumRadius()

double DrumRot::getDrumRadius ( )
inline
283  {
284  return drumRadius;
285  }
double drumRadius
Definition: MonodispersedDrum.cpp:308

References drumRadius.

Referenced by main().

◆ setCOR()

void DrumRot::setCOR ( double  drumCOR,
double  COR1 
)
inline
258  {
259  CORDrum = drumCOR;
260  CORS1 = COR1;
261  }
double CORS1
Definition: MonodispersedDrum.cpp:296
double CORDrum
Definition: MonodispersedDrum.cpp:296

References CORDrum, and CORS1.

Referenced by main().

◆ setDrumFillFraction()

void DrumRot::setDrumFillFraction ( double  dff)
inline
247  {
248  drumFillFraction = dff;
249  }
double drumFillFraction
Definition: MonodispersedDrum.cpp:300

References drumFillFraction.

Referenced by main().

◆ setDrumRadius()

void DrumRot::setDrumRadius ( double  radius)
inline
232  {
233  drumRadius = radius;
234  }

References drumRadius.

Referenced by main().

◆ setFractionalPolydispersity()

void DrumRot::setFractionalPolydispersity ( double  fpd)
inline
242  {
244  }
double fractionalPolydispersity
Definition: MonodispersedDrum.cpp:310

References fractionalPolydispersity.

Referenced by main().

◆ setFrictionCoeff()

void DrumRot::setFrictionCoeff ( double  pwf,
double  ppf 
)
inline
252  {
253  particleWallFriction = pwf;
255  }
double particleParticleFriction
Definition: MonodispersedDrum.cpp:303
double particleWallFriction
Definition: MonodispersedDrum.cpp:303

References particleParticleFriction, and particleWallFriction.

◆ setRevolutionSpeed()

void DrumRot::setRevolutionSpeed ( double  rpm)
inline
237  {
238  revolutionsPerSecond = rpm/60.0;// non-dimensionalised based on 3 mm particles and g=9.81
239  }

References revolutionsPerSecond.

Referenced by main().

◆ setRollingFriction()

void DrumRot::setRollingFriction ( double  drum,
double  f1 
)
inline
271  {
272  rollingFrictionDrum = drum;
273  rollingFriction1 = f1;
274  }
double rollingFrictionDrum
Definition: MonodispersedDrum.cpp:315
double rollingFriction1
Definition: MonodispersedDrum.cpp:316

References rollingFriction1, and rollingFrictionDrum.

Referenced by main().

◆ setSlidingFriction()

void DrumRot::setSlidingFriction ( double  drum,
double  f1 
)
inline
265  {
266  slidingFrictionDrum = drum;
267  slidingFriction1 = f1;
268  }
double slidingFriction1
Definition: MonodispersedDrum.cpp:313
double slidingFrictionDrum
Definition: MonodispersedDrum.cpp:312

References slidingFriction1, and slidingFrictionDrum.

Referenced by main().

◆ setTorsionFriction()

void DrumRot::setTorsionFriction ( double  drum,
double  f1 
)
inline
277  {
278  torsionFrictionDrum = drum;
279  torsionFriction1 = f1;
280  }
double torsionFrictionDrum
Definition: MonodispersedDrum.cpp:318
double torsionFriction1
Definition: MonodispersedDrum.cpp:319

References torsionFriction1, and torsionFrictionDrum.

Referenced by main().

◆ setupInitialConditions()

void DrumRot::setupInitialConditions ( )
inlineoverridevirtual

This function allows to set the initial conditions for our problem to be solved, by default particle locations are randomly set. Remember particle properties must also be defined here.

A virtual function with no implementation but can be overriden.

Todo:
I (Anthony) wants to change this to be an external function. This has a lot of advantages especially when using copy-constructors. This is a major change and will break other codes, so therefore has to be done carefully.

This sets up the particles initial conditions it is as you expect the user to override this. By default the particles are randomly distributed

Reimplemented from DPMBase.

42  {
43  radiusS1 = 0.0005; // 1mm diameter
44 
45  rhoS1 = 2500.0;
46 
47  massS1 = 4 / 3 * constants::pi * pow(radiusS1, 3.0) * rhoS1;
48 
49  double fillVolume = drumFillFraction * constants::pi * pow(drumRadius, 2.0) * (std::abs(getYMax() - getYMin()));
50  numS1 = (fillVolume) / (4. / 3. * constants::pi * pow(radiusS1, 3.0)); // Drum volume - volume large particles
51  //std::cout << "fillVolume" << fillVolume << "total particle volume" << numS1*4 / 3 * constants::pi * pow(radiusS1, 3.0) + numS2*4 / 3 * constants::pi * pow(radiusS2, 3.0) << std::endl;
52  //std::cout<< "fillVolume = " << fillVolume << ", drumFillFraction = "<< drumFillFraction << ", S2Volume = " << S2Volume << std::endl;
53  logger(INFO, "Expected numS1 = %", numS1);
55 
56  //std::cout << " mass " << massS1 << " " << massS2 << std::endl;
57 
58  tc = 1 / 4000.0;
59  //original value
60  //tc = 0.005;
61 
63 
66 
67  double RPSInitial = 0.0;
68 
69  speciesDrum->setDensity(rhoS1);
70  speciesDrum->setCollisionTimeAndRestitutionCoefficient(tc, CORDrum, massS1);
71  speciesDrum->setSlidingDissipation(speciesDrum->getDissipation()*2./7.);
72  speciesDrum->setSlidingStiffness(speciesDrum->getStiffness()*2./7.);
73  speciesDrum->setSlidingFrictionCoefficient(slidingFrictionDrum);
74  speciesDrum->setRollingStiffness(speciesDrum->getStiffness()*2.0/7.0);
75  speciesDrum->setRollingFrictionCoefficient(rollingFrictionDrum);
76  speciesDrum->setRollingDissipation(speciesDrum->getDissipation()*2./7.);
77  speciesDrum->setTorsionStiffness(speciesDrum->getStiffness()*2.0/7.0);
78  speciesDrum->setTorsionFrictionCoefficient(torsionFrictionDrum);
79  speciesDrum->setTorsionDissipation(speciesDrum->getDissipation()*2./7.);
80  //
81 
82  //
83  speciesS1->setDensity(rhoS1);
84  speciesS1->setCollisionTimeAndRestitutionCoefficient(tc, CORS1, massS1);
85  speciesS1->setSlidingDissipation(speciesS1->getDissipation()*2./7.);
86  speciesS1->setSlidingStiffness(speciesS1->getStiffness()*2./7.);
87  speciesS1->setSlidingFrictionCoefficient(slidingFriction1);
88  speciesS1->setRollingStiffness(speciesS1->getStiffness()*2.0/7.0);
89  speciesS1->setRollingFrictionCoefficient(rollingFriction1);
90  speciesS1->setRollingDissipation(speciesS1->getDissipation()*2./7.);
91  speciesS1->setTorsionStiffness(speciesS1->getStiffness()*2.0/7.0);
92  speciesS1->setTorsionFrictionCoefficient(torsionFriction1);
93  speciesS1->setTorsionDissipation(speciesS1->getDissipation()*2./7.);
94  //
95 
96  auto speciesDrumAndS1 = speciesHandler.getMixedObject(speciesDrum,speciesS1);
97  speciesDrumAndS1->setCollisionTimeAndRestitutionCoefficient(tc, ((CORS1 + CORDrum) / 2) , massS1, massS1);
98  speciesDrumAndS1->setSlidingDissipation(speciesDrumAndS1->getDissipation()*2./7.);
99  speciesDrumAndS1->setSlidingFrictionCoefficient( ((slidingFrictionDrum + slidingFriction1)/2));
100  speciesDrumAndS1->setSlidingStiffness(speciesDrumAndS1->getStiffness()*2.0/7.0);
101  speciesDrumAndS1->setRollingStiffness(speciesDrumAndS1->getStiffness()*2.0/7.0);
102  speciesDrumAndS1->setRollingFrictionCoefficient(((rollingFrictionDrum + rollingFriction1)/2));
103  speciesDrumAndS1->setRollingDissipation(speciesDrumAndS1->getDissipation()*2./7.);
104  speciesDrumAndS1->setTorsionStiffness(speciesDrumAndS1->getStiffness()*2.0/7.0);
105  speciesDrumAndS1->setTorsionFrictionCoefficient(((torsionFrictionDrum + torsionFriction1)/2));
106  speciesDrumAndS1->setTorsionDissipation(speciesDrumAndS1->getDissipation()*2./7.);
107  //
108 
109 
110  //
111  Vec3D drumCenter = {0.5*(getXMin() + getXMax()),
112  0.5*(getYMin() + getYMax()),
113  0.5*(getZMin() + getZMax())};
114 
115  wallHandler.clear();
116 
118  drumWall->setSpecies(speciesDrum);
119  drumWall->setPosition(drumCenter);
120  drumWall->setOrientation(Vec3D(0.0,1.0,0.0));
121  drumWall->addObject(Vec3D(1,0,0), Vec3D(drumRadius,0.0,0.0));
122  drumWall->setAngularVelocity(Vec3D(0.0,RPSInitial * 2.0 * constants::pi,0.0));
123 
124  InfiniteWall w0;
125  w0.setSpecies(speciesDrum);
126  w0.set(Vec3D(0.,-1.,0.),Vec3D(drumCenter.X,getYMin(),drumCenter.Z));
128  w0.set(Vec3D(0.,1.,0.),Vec3D(drumCenter.X,getYMax(),drumCenter.Z));
130 
132  double radius = 0.0;
133  int numS1Inserted=0;
134  Vec3D pos;
135  double R, C;
136  int failCounter = 0;
137 
138  while (numS1Inserted < numS1)
139  {
140  failCounter = 0;
141  do
142  {
143  p0.setSpecies(speciesS1);
144  p0.setRadius(radiusS1);
145  R = random.getRandomNumber(2*radiusS1,drumRadius-2.0*p0.getRadius()); // Find a random radial position within the cylinder
146  C = random.getRandomNumber(0.0*constants::pi,2.0*constants::pi); // Find a random angular position within the cylinder
147  // Set the position in the cylinder with a random z-coordinate
149  p0.setVelocity(Vec3D(0.0,0.0,0.0)); // Set initial particle velocity
150  failCounter++; // Adding 1 to failcounter
151  if (failCounter==1000)
152  {
153  break; // If we failed 1000 times to find a non-contact position we break the placement loop (makes sure that no infinite loops are created in combination with line 130)
154  }
155  } while(!checkParticleForInteraction(p0));
157  numS1Inserted++;
158  hGridRebuild();
159  }
160  logger(INFO, "Number of S1 particles inserted %", numS1Inserted);
161 
162  if (numS1Inserted == numS1)
163  {
164  step = 2;
165  checkTime = getTime() + 1.0;
166  logger(INFO, "\n \n \n"
167  "Particles settling down, checkTime = %"
168  "--------------------------"
169  "\n\n\n", checkTime);
170  }
171  }
Species< LinearViscoelasticNormalSpecies, FrictionSpecies > LinearViscoelasticFrictionSpecies
Definition: LinearViscoelasticFrictionSpecies.h:34
@ R
Definition: StatisticsVector.h:42
Use AxisymmetricIntersectionOfWalls to Screw Screw::read Screw::read Screw::read define axisymmetric ...
Definition: AxisymmetricIntersectionOfWalls.h:126
virtual void clear()
Empties the whole BaseHandler by removing all Objects and setting all other variables to 0.
Definition: BaseHandler.h:528
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:379
void setVelocity(const Vec3D &velocity)
set the velocity of the BaseInteractable.
Definition: BaseInteractable.cc:350
virtual void setPosition(const Vec3D &position)
Sets the position of this BaseInteractable.
Definition: BaseInteractable.h:239
Mdouble getRadius() const
Returns the particle's radius.
Definition: BaseParticle.h:348
virtual void setRadius(Mdouble radius)
Sets the particle's radius_ (and adjusts the mass_ accordingly, based on the particle's species)
Definition: BaseParticle.cc:553
void setSpecies(const ParticleSpecies *species)
Definition: BaseParticle.cc:818
void setSpecies(const ParticleSpecies *species)
Defines the species of the current wall.
Definition: BaseWall.cc:169
Mdouble getXMin() const
If the length of the problem domain in x-direction is XMax - XMin, then getXMin() returns XMin.
Definition: DPMBase.h:619
Mdouble getXMax() const
If the length of the problem domain in x-direction is XMax - XMin, then getXMax() returns XMax.
Definition: DPMBase.h:626
SpeciesHandler speciesHandler
A handler to that stores the species type i.e. LinearViscoelasticSpecies, etc.
Definition: DPMBase.h:1427
Mdouble getYMin() const
If the length of the problem domain in y-direction is YMax - YMin, then getYMin() returns YMin.
Definition: DPMBase.h:632
ParticleHandler particleHandler
An object of the class ParticleHandler, contains the pointers to all the particles created.
Definition: DPMBase.h:1437
RNG random
This is a random generator, often used for setting up the initial conditions etc.....
Definition: DPMBase.h:1432
Mdouble getYMax() const
If the length of the problem domain in y-direction is YMax - YMin, then getYMax() returns XMax.
Definition: DPMBase.h:638
Mdouble getZMax() const
If the length of the problem domain in z-direction is ZMax - ZMin, then getZMax() returns ZMax.
Definition: DPMBase.h:650
Mdouble getZMin() const
If the length of the problem domain in z-direction is ZMax - ZMin, then getZMin() returns ZMin.
Definition: DPMBase.h:644
int numS1
Definition: MonodispersedDrum.cpp:305
double massS1
Definition: MonodispersedDrum.cpp:294
double radiusS1
Definition: MonodispersedDrum.cpp:292
int numS1ToBeInserted
Definition: MonodispersedDrum.cpp:306
double tc
Definition: MonodispersedDrum.cpp:296
double rhoS1
Definition: MonodispersedDrum.cpp:293
A infinite wall fills the half-space {point: (position_-point)*normal_<=0}.
Definition: InfiniteWall.h:48
void set(Vec3D normal, Vec3D point)
Defines a standard wall, given an outward normal vector s.t. normal*x=normal*point for all x of the w...
Definition: InfiniteWall.cc:118
bool checkParticleForInteraction(const BaseParticle &P) final
Checks if given BaseParticle has an interaction with a BaseWall or other BaseParticle.
Definition: MercuryBase.cc:594
void hGridRebuild()
This sets up the parameters required for the contact model.
Definition: MercuryBase.cc:204
Mdouble getRandomNumber()
This is a random generating routine can be used for initial positions.
Definition: RNG.cc:143
void clear() override
Empties the whole BaseHandler by removing all Objects and setting all other variables to 0.
Definition: SpeciesHandler.h:54
std::enable_if<!std::is_pointer< typename U::MixedSpeciesType >::value, typename U::MixedSpeciesType * >::type getMixedObject(const U *S, const U *T)
Definition: SpeciesHandler.h:74
A spherical particle is the most simple particle used in MercuryDPM.
Definition: SphericalParticle.h:37
Mdouble Z
Definition: Vector.h:66
Mdouble X
the vector components
Definition: Vector.h:66
Mdouble cos(Mdouble x)
Definition: ExtendedMath.cc:64
Mdouble sin(Mdouble x)
Definition: ExtendedMath.cc:44

References MercuryBase::checkParticleForInteraction(), checkTime, BaseHandler< T >::clear(), SpeciesHandler::clear(), BaseHandler< T >::copyAndAddObject(), CORDrum, CORS1, mathsFunc::cos(), drumFillFraction, drumRadius, SpeciesHandler::getMixedObject(), BaseParticle::getRadius(), RNG::getRandomNumber(), DPMBase::getTime(), DPMBase::getXMax(), DPMBase::getXMin(), DPMBase::getYMax(), DPMBase::getYMin(), DPMBase::getZMax(), DPMBase::getZMin(), MercuryBase::hGridRebuild(), INFO, logger, massS1, numS1, numS1ToBeInserted, DPMBase::particleHandler, constants::pi, R, radiusS1, DPMBase::random, rhoS1, rollingFriction1, rollingFrictionDrum, InfiniteWall::set(), BaseInteractable::setPosition(), BaseParticle::setRadius(), BaseParticle::setSpecies(), BaseWall::setSpecies(), BaseInteractable::setVelocity(), mathsFunc::sin(), slidingFriction1, slidingFrictionDrum, DPMBase::speciesHandler, step, tc, torsionFriction1, torsionFrictionDrum, DPMBase::wallHandler, Vec3D::X, and Vec3D::Z.

Member Data Documentation

◆ checkTime

double DrumRot::checkTime
private

◆ CORDrum

double DrumRot::CORDrum
private

Referenced by setCOR(), and setupInitialConditions().

◆ CORS1

double DrumRot::CORS1
private

Referenced by setCOR(), and setupInitialConditions().

◆ densityRatio

double DrumRot::densityRatio
private

◆ drumFillFraction

double DrumRot::drumFillFraction
private

◆ drumRadius

double DrumRot::drumRadius
private

◆ fractionalPolydispersity

double DrumRot::fractionalPolydispersity
private

◆ massS1

double DrumRot::massS1
private

Referenced by setupInitialConditions().

◆ numS1

int DrumRot::numS1
private

Referenced by setupInitialConditions().

◆ numS1ToBeInserted

int DrumRot::numS1ToBeInserted
private

Referenced by setupInitialConditions().

◆ particleParticleFriction

double DrumRot::particleParticleFriction
private

Referenced by setFrictionCoeff().

◆ particleWallFriction

double DrumRot::particleWallFriction
private

Referenced by setFrictionCoeff().

◆ radiusS1

double DrumRot::radiusS1
private

Referenced by setupInitialConditions().

◆ revolutionsPerSecond

double DrumRot::revolutionsPerSecond
private

◆ rhoS1

double DrumRot::rhoS1
private

Referenced by setupInitialConditions().

◆ rollingFriction1

double DrumRot::rollingFriction1
private

◆ rollingFrictionDrum

double DrumRot::rollingFrictionDrum
private

◆ sizeRatio

double DrumRot::sizeRatio
private

◆ slidingFriction1

double DrumRot::slidingFriction1
private

◆ slidingFrictionDrum

double DrumRot::slidingFrictionDrum
private

◆ step

int DrumRot::step
private

◆ tc

double DrumRot::tc
private

Referenced by setupInitialConditions().

◆ torsionFriction1

double DrumRot::torsionFriction1
private

◆ torsionFrictionDrum

double DrumRot::torsionFrictionDrum
private

◆ volumeFraction

double DrumRot::volumeFraction
private

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