ls1-MarDyn
ls1-MarDyn molecular dynamics code
|
Controls the simulation process. More...
#include <Simulation.h>
Public Member Functions | |
Simulation () | |
~Simulation () | |
void | readXML (XMLfileUnits &xmlconfig) |
Read in XML configuration for simulation and all its included objects. More... | |
void | readConfigFile (std::string filename) |
process configuration file More... | |
void | initConfigXML (const std::string &inputfilename) |
Opens given XML file and reads in parameters for the simulaion. More... | |
void | prepare_start () |
calculate all values for the starting timepoint More... | |
void | simulate () |
Controls the main loop of the simulation. More... | |
void | pluginEndStepCall (unsigned long simstep) |
call plugins every nth-simstep More... | |
void | finalize () |
clean up simulation | |
void | updateParticleContainerAndDecomposition (double lastTraversalTime, bool useTimers) |
void | performOverlappingDecompositionAndCellTraversalStep (double etime) |
void | setDomainDecomposition (DomainDecompBase *domainDecomposition) |
DomainDecompBase & | domainDecomposition () |
Domain * | getDomain () |
Integrator * | getIntegrator () |
ParticleContainer * | getMoleculeContainer () |
void | setNumTimesteps (unsigned long steps) |
unsigned long | getNumTimesteps () |
unsigned long | getNumInitTimesteps () |
unsigned long | getSimulationStep () |
void | setLoopAbortTime (double time) |
double | getcutoffRadius () const |
void | setcutoffRadius (double cutoffRadius) |
double | getLJCutoff () const |
void | setLJCutoff (double LJCutoffRadius) |
unsigned long | getTotalNumberOfMolecules () const |
double | Tfactor (unsigned long simstep) |
Temperature increase factor function during automatic equilibration. More... | |
void | initCanonical (unsigned long t) |
void | initGrandCanonical (unsigned long t) |
void | initStatistics (unsigned long t) |
unsigned long | getInitStatistics () const |
void | setSimulationTime (double curtime) |
void | advanceSimulationTime (double timestep) |
double | getSimulationTime () |
void | setEnsemble (Ensemble *ensemble) |
Ensemble * | getEnsemble () |
std::shared_ptr< MemoryProfiler > | getMemoryProfiler () |
TimerProfiler * | timers () |
double | getH () |
get Planck constant | |
void | setH (double h_extern) |
set Planck constant | |
unsigned long | getNumberOfTimesteps () const |
double | getAndResetOneLoopCompTime () |
computational time for one execution of traverseCell | |
void | setOutputPrefix (std::string prefix) |
void | setOutputPrefix (char *prefix) |
std::string | getOutputPrefix () |
void | enableFinalCheckpoint () |
void | disableFinalCheckpoint () |
void | useLegacyCellProcessor () |
void | enableMemoryProfiler () |
void | setForcedCheckpointTime (double time) |
void | initialize () |
PluginBase * | getPlugin (const std::string &name) |
get plugin More... | |
std::list< PluginBase * > * | getPluginList () |
void | initGlobalEnergyLog () |
void | writeGlobalEnergyLog (const double &globalUpot, const double &globalT, const double &globalPressure) |
CellProcessor * | getCellProcessor () const |
void | refreshParticleIDs () |
Refresh particle IDs to continuous numbering. | |
bool | keepRunning () |
Checks if Simsteps or MaxWallTime are reached. | |
Static Public Member Functions | |
static void | exit (int exitcode) |
Terminate simulation with given exit code. More... | |
Controls the simulation process.
Simulation parameters are provided via a xml config file or can be set directly via the corresponding methods.
Simulation::Simulation | ( | ) |
Instantiate simulation object
Simulation::~Simulation | ( | ) |
destruct simulation object
|
inline |
Return a reference to the domain decomposition used in the simulation
|
static |
Terminate simulation with given exit code.
The exit method takes care over the right way to terminate the application in a correct way for the different parallelization schemes. e.g. terminating other processes in MPI parallel execution mode.
|
inline |
Get pointer to the domain
|
inline |
Get pointer to the integrator
|
inline |
Get pointer to the molecule container
|
inline |
Get initial number of steps
|
inline |
Get the number of time steps to be performed in the simulation
PluginBase * Simulation::getPlugin | ( | const std::string & | name | ) |
get plugin
|
inline |
Get the number of the actual time step currently processed in the simulation.
void Simulation::initConfigXML | ( | const std::string & | inputfilename | ) |
Opens given XML file and reads in parameters for the simulaion.
[in] | inputfilename | filename of the XML input file |
void Simulation::initGlobalEnergyLog | ( | ) |
Global energy log
void Simulation::initialize | ( | ) |
initialize all member variables with a suitable value
void Simulation::performOverlappingDecompositionAndCellTraversalStep | ( | double | etime | ) |
Performs both the decomposition and the cell traversal in an overlapping way. The overlapping is needed to speed up the overall computation. The order of cells traversed will be different, than for the non-overlapping case, slightly different results are possible.
void Simulation::pluginEndStepCall | ( | unsigned long | simstep | ) |
call plugins every nth-simstep
The present method serves as a redirection to the actual plugins. That includes a) particular plugin objects included in _plugins,
[in] | simstep | timestep of the plugins |
void Simulation::prepare_start | ( | ) |
calculate all values for the starting timepoint
After the input file has been read in, only the information which is directly in the inputfile is known at time step zero. This includes e.g. Molecule position and velocities, but not forces which act on the molecules or macroscopic values at the potential. This method has to ensure, that all values are available at time step zero
Init TemperatureControl beta_trans, beta_rot log-files, register as observer if plugin DistControl is in use.
refresh particle IDs
void Simulation::readConfigFile | ( | std::string | filename | ) |
process configuration file
calls initConfigXML
[in] | filename | filename of the input file |
void Simulation::readXML | ( | XMLfileUnits & | xmlconfig | ) |
Read in XML configuration for simulation and all its included objects.
The following xml object structure is handled by this method:
Prepare start options, affecting behavior of method prepare_start()
void Simulation::setDomainDecomposition | ( | DomainDecompBase * | domainDecomposition | ) |
Set the private _domainDecomposition variable to a new pointer.
domainDecomposition | the new va |
|
inline |
Set Loop Time Limit in seconds
|
inline |
Set the number of time steps to be performed in the simulation
void Simulation::simulate | ( | ) |
Controls the main loop of the simulation.
precondition for this method is that initialize() has been called. The main loop calls all methods which have to be called in each iteration step, e.g. initializes the molecule Container, calculates the forces, ...
For the integration, a seperate integration object is used. This method is written in a way that should make it possible to use different integrators without having to change anything. Whenever something is done which might make it necessary for an integrator to do something, the integrator is informed about that using the integrators corresponding method (see integrator documentation) It follows a coarse outline of what has to be done:
< timer for the entire simulation loop (synced)
< timer for decomposition: sub-timer of loopTimer
< timer for computation: sub-timer of loopTimer
< timer for io in simulation loop: sub-timer of loopTimer
< timer for force calculation: sub-timer of computationTimer
< timer for measuring MPI-OMP communication time: sub-timer of decompositionTimer
double Simulation::Tfactor | ( | unsigned long | simstep | ) |
Temperature increase factor function during automatic equilibration.
[in] | current | simulation time step |
The present version of Mardyn provides the option of equilibrating the system at an increased temperature. This has, of course, the advantage that the equilibration is accelerated. In a special case that is of particular interest to a relevant subset of the developers, namely MD simulation of nucleation, it also reduces the amount of clusters formed during equilibration. Then, the actual nucleation process can be investigated for the equilibrated system.
Equilibration at an increased temperature is nothing new, in principle it could be considered a variant of simulated annealing. For the purpose of accelerating the relaxation, starting from an unrealistic initial state, it was e.g. proposed in the thesis of M. Kreitmeir.
The key parameter to regulating the equilibration is this->_initCanonical. The system is gradually heated and cooled again, and from this->_initCanonical on the temperature is constant (i.e. this method returns 1.0).
Thus, if temperature increase is to be avoided, the user should set the "initCanonical" parameter to "0".
void Simulation::updateParticleContainerAndDecomposition | ( | double | lastTraversalTime, |
bool | useTimers | ||
) |
The following things have to be done here: