Molecular Dynamics Simulation  1.0
Public Member Functions | Private Attributes | Friends | List of all members
Simulation Class Reference

Class to run a simulation. More...

#include <Simulation.h>

Collaboration diagram for Simulation:
Collaboration graph

Public Member Functions

 Simulation (const std::vector< Particle > &particles, const SimulationParams &params, IntegrationMethod integration_method=IntegrationMethod::VERLET)
 Construct a new Simulation object and initialize all the necessary components. More...
 
 ~Simulation ()
 
SimulationOverview runSimulation ()
 Runs the simulation, using the parameters given at construction and returns a SimulationOverview object containing some data. More...
 

Private Attributes

const SimulationParamsparams
 Reference to the simulation parameters object. More...
 
std::unique_ptr< ParticleContainerparticle_container
 Reference to the ParticleContainer on whose content the simulation is performed. More...
 
std::unique_ptr< IntegrationFunctorintegration_functor
 Functor used to integrate the particles. More...
 

Friends

class ProgressBarInterceptor
 
class FrameWriterInterceptor
 
class TemperatureSensorInterceptor
 
class ThermostatInterceptor
 
class ParticleUpdateCounterInterceptor
 
class RadialDistributionFunctionInterceptor
 
class DiffusionFunctionInterceptor
 
class VelocityProfileInterceptor
 

Detailed Description

Class to run a simulation.

This class collects all the components needed to run a simulation, and provides a method to run it.

Definition at line 20 of file Simulation.h.

Constructor & Destructor Documentation

◆ Simulation()

Simulation::Simulation ( const std::vector< Particle > &  particles,
const SimulationParams params,
IntegrationMethod  integration_method = IntegrationMethod::VERLET 
)

Construct a new Simulation object and initialize all the necessary components.

Parameters
particlesReference to the ParticleContainer on whose content the simulation is performed
paramsParameters for the simulation. See the class SimulationParams for more information
integration_methodThe integration method to use for the simulation (Default: IntegrationMethod::VERLET)

Definition at line 21 of file Simulation.cpp.

22  : params(params), integration_functor(get_integration_functor(integration_method)) {
23  // Create particle container
24  if (std::holds_alternative<SimulationParams::LinkedCellsType>(params.container_type)) {
25  auto lc_type = std::get<SimulationParams::LinkedCellsType>(params.container_type);
27  std::make_unique<LinkedCellsContainer>(lc_type.domain_size, lc_type.cutoff_radius, lc_type.boundary_conditions);
28  } else if (std::holds_alternative<SimulationParams::DirectSumType>(params.container_type)) {
29  particle_container = std::make_unique<DirectSumContainer>();
30  } else {
31  throw std::runtime_error("Unknown container type");
32  }
33 
34  // Add particles to container
35  particle_container->reserve(initial_particles.size());
36  for (auto& particle : initial_particles) {
37  particle_container->addParticle(particle);
38  }
39 }
std::unique_ptr< IntegrationFunctor > get_integration_functor(IntegrationMethod method)
Returns the corresponding integration functor for the given integration method.
std::variant< DirectSumType, LinkedCellsType > container_type
Type of the particle container.
std::unique_ptr< ParticleContainer > particle_container
Reference to the ParticleContainer on whose content the simulation is performed.
Definition: Simulation.h:50
std::unique_ptr< IntegrationFunctor > integration_functor
Functor used to integrate the particles.
Definition: Simulation.h:55
const SimulationParams & params
Reference to the simulation parameters object.
Definition: Simulation.h:45

◆ ~Simulation()

Simulation::~Simulation ( )
default

Member Function Documentation

◆ runSimulation()

SimulationOverview Simulation::runSimulation ( )

Runs the simulation, using the parameters given at construction and returns a SimulationOverview object containing some data.

Returns
SimulationOverview object containing some data about the simulation performed

Definition at line 43 of file Simulation.cpp.

43  {
44  size_t iteration = params.start_iteration;
45  double simulated_time = params.start_iteration * params.delta_t;
46 
47  // Calculate initial forces
48  particle_container->prepareForceCalculation();
49  particle_container->applySimpleForces(params.simple_forces);
50  particle_container->applyPairwiseForces(params.pairwise_forces);
51  particle_container->applyTargettedForces(params.targetted_forces, simulated_time);
52 
53  Logger::logger->info("Simulation started...");
54 
55  std::time_t t_start_helper = std::chrono::system_clock::to_time_t(std::chrono::system_clock::now());
56  Logger::logger->info("Start time: {}", fmt::format("{:%A %Y-%m-%d %H:%M:%S}", fmt::localtime(t_start_helper)));
57 
58  // Notify interceptors that the simulation is about to start
59  for (auto& interceptor : params.interceptors) {
60  (*interceptor).onSimulationStart(*this);
61  }
62 
63  auto t_start = std::chrono::high_resolution_clock::now();
64 
65  while (simulated_time < params.end_time) {
67  simulated_time);
68 
69  ++iteration;
70  simulated_time += params.delta_t;
71 
72  // Notify interceptors of the current iteration
73  for (auto& interceptor : params.interceptors) {
74  (*interceptor).notify(iteration, *this);
75  }
76  }
77 
78  auto t_end = std::chrono::high_resolution_clock::now();
79 
80  // Notify interceptors that the simulation has ended
81  for (auto& interceptor : params.interceptors) {
82  (*interceptor).onSimulationEnd(iteration, *this);
83  }
84 
85  Logger::logger->info("Simulation finished.");
86  Logger::logger->info("End time: {}", fmt::format("{:%A %Y-%m-%d %H:%M:%S}", fmt::localtime(t_end)));
87 
88  std::vector<std::string> interceptor_summaries;
89  for (auto& interceptor : params.interceptors) {
90  auto summary = std::string(*interceptor);
91  if (!summary.empty()) {
92  interceptor_summaries.push_back(summary);
93  }
94  }
95 
96  auto total_time_ms = std::chrono::duration_cast<std::chrono::milliseconds>(t_end - t_start).count();
97 
98  SimulationOverview overview{params, total_time_ms / 1000.0, iteration, interceptor_summaries,
99  std::vector<Particle>(particle_container->begin(), particle_container->end())};
100 
101  return overview;
102 }
static std::shared_ptr< spdlog::logger > logger
Publically accessible shared pointer to the logger.
Definition: Logger.h:35
Class to store some overview data of an executed simulation.
std::vector< std::shared_ptr< PairwiseForceSource > > pairwise_forces
Pairwise Forces to be applied to the particles.
double delta_t
Time step of a single simulation iteration.
std::vector< std::shared_ptr< SimulationInterceptor > > interceptors
List of interceptors to be used in the simulation.
size_t start_iteration
Start iteration of the simulation.
std::vector< std::shared_ptr< SimpleForceSource > > simple_forces
Simple Forces to be applied to the particles.
std::vector< std::shared_ptr< TargettedForceSource > > targetted_forces
Targetted Forces to be applied to the particles.
double end_time
End time of the simulation.

Friends And Related Function Documentation

◆ DiffusionFunctionInterceptor

friend class DiffusionFunctionInterceptor
friend

Definition at line 66 of file Simulation.h.

◆ FrameWriterInterceptor

friend class FrameWriterInterceptor
friend

Definition at line 61 of file Simulation.h.

◆ ParticleUpdateCounterInterceptor

friend class ParticleUpdateCounterInterceptor
friend

Definition at line 64 of file Simulation.h.

◆ ProgressBarInterceptor

friend class ProgressBarInterceptor
friend

Befriend the interceptors to allow them to access the private members of this class

Definition at line 60 of file Simulation.h.

◆ RadialDistributionFunctionInterceptor

Definition at line 65 of file Simulation.h.

◆ TemperatureSensorInterceptor

friend class TemperatureSensorInterceptor
friend

Definition at line 62 of file Simulation.h.

◆ ThermostatInterceptor

friend class ThermostatInterceptor
friend

Definition at line 63 of file Simulation.h.

◆ VelocityProfileInterceptor

friend class VelocityProfileInterceptor
friend

Definition at line 67 of file Simulation.h.

Member Data Documentation

◆ integration_functor

std::unique_ptr<IntegrationFunctor> Simulation::integration_functor
private

Functor used to integrate the particles.

Definition at line 55 of file Simulation.h.

◆ params

const SimulationParams& Simulation::params
private

Reference to the simulation parameters object.

Definition at line 45 of file Simulation.h.

◆ particle_container

std::unique_ptr<ParticleContainer> Simulation::particle_container
private

Reference to the ParticleContainer on whose content the simulation is performed.

Definition at line 50 of file Simulation.h.


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