ls1-MarDyn
ls1-MarDyn molecular dynamics code
Simulation.h
1#ifndef SIMULATION_H_
2#define SIMULATION_H_
3
4#include <memory>
5#include <any>
6
7#include "ensemble/CavityEnsemble.h"
8#include "io/TimerProfiler.h"
9#include "thermostats/VelocityScalingThermostat.h"
10#include "utils/FixedSizeQueue.h"
11#include "utils/FunctionWrapper.h"
12#include "utils/SysMon.h"
13
14// plugins
15#include "plugins/PluginFactory.h"
16
17#if !defined (SIMULATION_SRC) or defined (IN_IDE_PARSER)
18class Simulation;
20extern Simulation* global_simulation;
21#endif
22
23class ParticleInsertion;
24class Ensemble;
25
26
27
29#define _simulation (*global_simulation)
30
31#include <list>
32#include <vector>
33#include <string>
34#include <io/TaskTimingProfiler.h>
35
36#ifdef STEEREO
37class SteereoSimSteering;
38class SteereoCouplingSim;
39#endif
40
41class Domain;
44class CellProcessor;
45class Integrator;
46class PluginBase;
48class InputBase;
49class RDF;
50class FlopCounter;
52class Homogeneous;
53class Planar;
55class MemoryProfiler;
56
57// by Stefan Becker
58const int VELSCALE_THERMOSTAT = 1;
59
60namespace bhfmm {
61class FastMultipoleMethod;
62} // bhfmm
63
64
71private:
72 Simulation(Simulation &simulation);
73 Simulation& operator=(Simulation &simulation);
74 void updateForces();
75
76public:
78 Simulation();
79
82
124 void readXML(XMLfileUnits& xmlconfig);
125
132 static void exit(int exitcode);
133
139 void readConfigFile(std::string filename);
140
144 void initConfigXML(const std::string& inputfilename);
145
155 void prepare_start();
156
180 void simulate();
181
182
191 void pluginEndStepCall(unsigned long simstep);
192
194 void finalize();
195
201 void updateParticleContainerAndDecomposition(double lastTraversalTime, bool useTimers);
202
209
215
217 DomainDecompBase& domainDecomposition() { return *_domainDecomposition; }
218
220 Domain* getDomain() { return _domain; }
221
223 Integrator* getIntegrator() {return _integrator; }
224
226 ParticleContainer* getMoleculeContainer() { return _moleculeContainer; }
227
229 void setNumTimesteps( unsigned long steps ) { _numberOfTimesteps = steps; }
231 unsigned long getNumTimesteps() { return _numberOfTimesteps; }
233 unsigned long getNumInitTimesteps() { return _initSimulation; }
235 unsigned long getSimulationStep() { return _simstep; }
237 void setLoopAbortTime(double time) {
238 global_log->info() << "Max loop-abort-time set: " << time << "\n";
239 _wallTimeEnabled = true;
240 _maxWallTime = time;
241 }
242
243 double getcutoffRadius() const { return _cutoffRadius; }
244 void setcutoffRadius(double cutoffRadius) { _cutoffRadius = cutoffRadius; }
245 double getLJCutoff() const { return _LJCutoffRadius; }
246 void setLJCutoff(double LJCutoffRadius) { _LJCutoffRadius = LJCutoffRadius; }
247 unsigned long getTotalNumberOfMolecules() const;
248
275 double Tfactor(unsigned long simstep);
276
277 void initCanonical(unsigned long t) { this->_initCanonical = t; }
278 void initGrandCanonical(unsigned long t) { this->_initGrandCanonical = t; }
279 void initStatistics(unsigned long t) { this->_initStatistics = t; }
280 unsigned long getInitStatistics() const { return this->_initStatistics; }
281
282 void setSimulationTime(double curtime) { _simulationTime = curtime; }
283 void advanceSimulationTime(double timestep) { _simulationTime += timestep; }
284 double getSimulationTime() { return _simulationTime; }
285
286 void setEnsemble(Ensemble *ensemble) { _ensemble = ensemble; }
287 Ensemble* getEnsemble() { return _ensemble; }
288
289 std::shared_ptr<MemoryProfiler> getMemoryProfiler() {
290 return _memoryProfiler;
291 }
292
293 TimerProfiler* timers() {
294 return &_timerProfiler;
295 }
296
298 double getH() {return h;}
300 void setH(double h_extern) {h = h_extern;}
301
302private:
303
304
305 double _simulationTime;
310 unsigned long _maxMoleculeId;
311
313 double _cutoffRadius;
314
316 double _LJCutoffRadius;
317
331 unsigned _collectThermostatDirectedVelocity;
332
335 int _thermostatType;
336
337 unsigned long _numberOfTimesteps;
338public:
339 unsigned long getNumberOfTimesteps() const;
340
341private:
344 unsigned long _simstep;
347 unsigned long _initSimulation;
349 unsigned long _initCanonical;
351 unsigned long _initGrandCanonical;
353 unsigned long _initStatistics;
354
355 Ensemble* _ensemble;
356
358 ParticleContainer* _moleculeContainer;
359
361 ParticlePairsHandler* _particlePairsHandler;
362
364 CellProcessor* _cellProcessor;
365
367 DomainDecompBase* _domainDecomposition;
368
370 Integrator* _integrator;
371
373 Domain* _domain;
374
376 InputBase* _inputReader;
377
379 std::string _outputPrefix;
380
382 unsigned _momentumInterval;
383
385 Random _rand;
386
388 LongRangeCorrection* _longRangeCorrection;
389
391 TemperatureControl* _temperatureControl;
392
395
397 TimerProfiler _timerProfiler;
398
400 std::shared_ptr<MemoryProfiler> _memoryProfiler;
401
402#ifdef TASKTIMINGPROFILE
404 TaskTimingProfiler* _taskTimingProfiler;
405#endif
406public:
407#ifdef TASKTIMINGPROFILE
408 TaskTimingProfiler* getTaskTimingProfiler(){
409 return _taskTimingProfiler;
410 }
411#endif
414 if(_loopCompTimeSteps==0){
415 return 1.;
416 }
417 double t = _loopCompTime/_loopCompTimeSteps;
418 _loopCompTime = 0.;
419 _loopCompTimeSteps = 0;
420 return t;
421 }
422 void setOutputPrefix( std::string prefix ) { _outputPrefix = prefix; }
423 void setOutputPrefix( char *prefix ) { _outputPrefix = std::string( prefix ); }
424 std::string getOutputPrefix() { return _outputPrefix; }
425
426 void enableFinalCheckpoint() { _finalCheckpoint = true; }
427 void disableFinalCheckpoint() { _finalCheckpoint = false; }
428
429 void useLegacyCellProcessor() { _legacyCellProcessor = true; }
430
431 void enableMemoryProfiler() {
432 _memoryProfiler = std::make_shared<MemoryProfiler>();
433 _memoryProfiler->registerObject(reinterpret_cast<MemoryProfilable**>(&_moleculeContainer));
434 _memoryProfiler->registerObject(reinterpret_cast<MemoryProfilable**>(&_domainDecomposition));
435 }
436
437 void setForcedCheckpointTime(double time) { _forced_checkpoint_time = time; }
438
440 void initialize();
441
446
447 std::list<PluginBase*>* getPluginList(){
448 return &_plugins;
449 }
450
453 void writeGlobalEnergyLog(const double& globalUpot, const double& globalT, const double& globalPressure);
454
455 CellProcessor *getCellProcessor() const;
456
458 void refreshParticleIDs();
459
461 bool keepRunning();
462
463private:
464
466 Timer* _timerForLoad{nullptr};
467
468 Timer _timeFromStart;
469 double _maxWallTime = -1;
470 bool _wallTimeEnabled = false;
471
473 bool _finalCheckpoint;
474
476 bool _legacyCellProcessor = false;
477
483 bool _overlappingP2P {false};
484
486 std::list<PluginBase*> _plugins;
487
493 std::map<std::string, FunctionWrapper> _callbacks;
494
495 VelocityScalingThermostat _velocityScalingThermostat;
496
504 double h;
505
507 double _forced_checkpoint_time;
508
510 double _loopCompTime;
511
512 int _loopCompTimeSteps;
513
515 unsigned long _nWriteFreqGlobalEnergy;
516 std::string _globalEnergyLogFilename;
517
524 struct PrepareStartOptions {
525 bool refreshIDs;
526 } _prepare_start_opt;
527
528 FixedSizeQueue<double> _lastTraversalTimeHistory;
529};
530#endif /*SIMULATION_H_*/
System monitoring class.
Definition: CellProcessor.h:29
handle boundary region and multiple processes
Definition: DomainDecompBase.h:51
This class is used to read in the phasespace and to handle macroscopic values.
Definition: Domain.h:47
Base class for ensembles.
Definition: EnsembleBase.h:47
Definition: FixedSizeQueue.h:15
Definition: FlopCounter.h:60
Definition: Homogeneous.h:15
interface for any kind of input class
Definition: InputBase.h:18
Update velocities and positions.
Definition: Integrator.h:32
Definition: LongRangeCorrection.h:12
Definition: MemoryProfiler.h:13
Definition: MemoryProfiler.h:22
This Interface is used to get access to particles and pairs of particles.
Definition: ParticleContainer.h:69
interface for defining the action performed when processing a pair
Definition: ParticlePairsHandler.h:38
Definition: Planar.h:24
The PluginBase class provides the interface for any kind of output/plugin classes - called "(output) ...
Definition: PluginBase.h:47
This class calculates the Radial Distribution Function (RDF).
Definition: RDF.h:37
Definition: Random.h:12
Controls the simulation process.
Definition: Simulation.h:70
void pluginEndStepCall(unsigned long simstep)
call plugins every nth-simstep
Definition: Simulation.cpp:1216
void setDomainDecomposition(DomainDecompBase *domainDecomposition)
Definition: Simulation.cpp:1321
Domain * getDomain()
Definition: Simulation.h:220
Integrator * getIntegrator()
Definition: Simulation.h:223
double getAndResetOneLoopCompTime()
computational time for one execution of traverseCell
Definition: Simulation.h:413
void performOverlappingDecompositionAndCellTraversalStep(double etime)
Definition: Simulation.cpp:1301
void readXML(XMLfileUnits &xmlconfig)
Read in XML configuration for simulation and all its included objects.
Definition: Simulation.cpp:160
void simulate()
Controls the main loop of the simulation.
Definition: Simulation.cpp:920
void initConfigXML(const std::string &inputfilename)
Opens given XML file and reads in parameters for the simulaion.
Definition: Simulation.cpp:687
double getH()
get Planck constant
Definition: Simulation.h:298
void refreshParticleIDs()
Refresh particle IDs to continuous numbering.
Definition: Simulation.cpp:1408
void readConfigFile(std::string filename)
process configuration file
Definition: Simulation.cpp:675
Simulation()
Definition: Simulation.cpp:85
PluginBase * getPlugin(const std::string &name)
get plugin
Definition: Simulation.cpp:1391
static void exit(int exitcode)
Terminate simulation with given exit code.
Definition: Simulation.cpp:155
void updateParticleContainerAndDecomposition(double lastTraversalTime, bool useTimers)
Definition: Simulation.cpp:1265
void initialize()
Definition: Simulation.cpp:1345
void setLoopAbortTime(double time)
Definition: Simulation.h:237
DomainDecompBase & domainDecomposition()
Definition: Simulation.h:217
unsigned long getSimulationStep()
Definition: Simulation.h:235
void setNumTimesteps(unsigned long steps)
Definition: Simulation.h:229
~Simulation()
Definition: Simulation.cpp:127
unsigned long getNumInitTimesteps()
Definition: Simulation.h:233
void setH(double h_extern)
set Planck constant
Definition: Simulation.h:300
void prepare_start()
calculate all values for the starting timepoint
Definition: Simulation.cpp:764
bool keepRunning()
Checks if Simsteps or MaxWallTime are reached.
Definition: Simulation.cpp:1372
void finalize()
clean up simulation
Definition: Simulation.cpp:1242
double Tfactor(unsigned long simstep)
Temperature increase factor function during automatic equilibration.
Definition: Simulation.cpp:1331
ParticleContainer * getMoleculeContainer()
Definition: Simulation.h:226
void initGlobalEnergyLog()
unsigned long getNumTimesteps()
Definition: Simulation.h:231
Definition: TemperatureControl.h:166
Class for managing timers across the simulation.
Definition: TimerProfiler.h:28
This class is used to measure times in sequential and parallel versions.
Definition: Timer.h:59
Definition: VelocityScalingThermostat.h:9
XML file with unit attributes abstraction.
Definition: xmlfileUnits.h:25
Definition: FastMultipoleMethod.h:45
Definition: L2PCellProcessor.cpp:15
::xsd::cxx::tree::name< char, token > name
C++ type corresponding to the Name XML Schema built-in type.
Definition: vtk-punstructured.h:288
::xsd::cxx::tree::time< char, simple_type > time
C++ type corresponding to the time XML Schema built-in type.
Definition: vtk-punstructured.h:438
::xsd::cxx::tree::string< char, simple_type > string
C++ type corresponding to the string XML Schema built-in type.
Definition: vtk-punstructured.h:270