ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
core/system/System.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2014-2022 The ESPResSo project
3 *
4 * This file is part of ESPResSo.
5 *
6 * ESPResSo is free software: you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation, either version 3 of the License, or
9 * (at your option) any later version.
10 *
11 * ESPResSo is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program. If not, see <http://www.gnu.org/licenses/>.
18 */
19
20#pragma once
21
22#include "config/config.hpp"
23
24#include "GpuParticleData.hpp"
25#include "ResourceCleanup.hpp"
26
29
30#include "ek/Solver.hpp"
31#include "lb/Solver.hpp"
32
34
35#include <utils/Vector.hpp>
36
37#include <memory>
38
39class BoxGeometry;
40class LocalBox;
41struct CellStructure;
42class Propagation;
44namespace Thermostat {
45class Thermostat;
46}
47class ComFixed;
48class Galilei;
49class Observable_stat;
50namespace BondBreakage {
51class BondBreakage;
52}
53namespace LeesEdwards {
54class LeesEdwards;
55}
56
57namespace System {
58
59/**
60 * @brief Main system class.
61 *
62 * Most components follow the composite pattern and the opaque pointer pattern.
63 * See @ref SystemClassDesign for more details.
64 */
65class System : public std::enable_shared_from_this<System> {
66private:
67 struct Private {};
68 void initialize();
69
70public:
71 System(Private);
72
73 static std::shared_ptr<System> create();
74
75#ifdef CUDA
77#endif
79
80 /** @brief Get @ref time_step. */
81 auto get_time_step() const { return time_step; }
82
83 /** @brief Set @ref time_step. */
84 void set_time_step(double value);
85
86 /** @brief Get @ref sim_time. */
87 auto get_sim_time() const { return sim_time; }
88
89 /** @brief Set @ref sim_time. */
90 void set_sim_time(double value);
91
92 /** @brief Get @ref force_cap. */
93 auto get_force_cap() const { return force_cap; }
94
95 /** @brief Set @ref force_cap. */
96 void set_force_cap(double value);
97
98 /** @brief Get @ref min_global_cut. */
99 auto get_min_global_cut() const { return min_global_cut; }
100
101 /** @brief Set @ref min_global_cut. */
102 void set_min_global_cut(double value);
103
104 /** @brief Change the box dimensions. */
105 void set_box_l(Utils::Vector3d const &box_l);
106
107 /**
108 * @brief Tune the Verlet skin.
109 * Choose the optimal Verlet list skin between @p min_skin and @p max_skin
110 * by bisection to tolerance @p tol.
111 */
112 void tune_verlet_skin(double min_skin, double max_skin, double tol,
113 int int_steps, bool adjust_max_skin);
114
115 /** @brief Change cell structure topology. */
116 void set_cell_structure_topology(CellStructureType topology);
117
118 /** @brief Rebuild cell lists. Use e.g. after a skin change. */
119 void rebuild_cell_structure();
120
121 /** @brief Calculate the maximal cutoff of all interactions. */
122 double maximal_cutoff() const;
123
124 /** @brief Get the interaction range. */
125 double get_interaction_range() const;
126
127 unsigned get_global_ghost_flags() const;
128
129 /** Check electrostatic and magnetostatic methods are properly initialized.
130 * @return true if sanity checks failed.
131 */
132 bool long_range_interactions_sanity_checks() const;
133
134 /** @brief Calculate the total energy. */
135 std::shared_ptr<Observable_stat> calculate_energy();
136
137 /** @brief Calculate the pressure from a virial expansion. */
138 std::shared_ptr<Observable_stat> calculate_pressure();
139
140 /** @brief Calculate all forces. */
141 void calculate_forces();
142
143#ifdef DIPOLE_FIELD_TRACKING
144 /** @brief Calculate dipole fields. */
145 void calculate_long_range_fields();
146#endif
147
148 /**
149 * @brief Compute the short-range energy of a particle.
150 *
151 * Iterate through particles inside cell and neighboring cells and compute
152 * energy contribution for a specific particle.
153 *
154 * @param pid Particle id
155 * @return Non-bonded energy of the particle.
156 */
157 double particle_short_range_energy_contribution(int pid);
158
159 /** Integrate equations of motion
160 * @param n_steps Number of integration steps, can be zero
161 * @param reuse_forces Decide when to re-calculate forces
162 *
163 * @details This function calls two hooks for propagation kernels such as
164 * velocity verlet, velocity verlet + npt box changes, and steepest_descent.
165 * One hook is called before and one after the force calculation.
166 * It is up to the propagation kernels to increment the simulation time.
167 *
168 * This function propagates the system according to the choice of integrator
169 * stored in @ref Propagation::integ_switch. The general structure is:
170 * - if reuse_forces is zero, recalculate the forces based on the current
171 * state of the system
172 * - Loop over the number of simulation steps:
173 * -# initialization (e.g., RATTLE)
174 * -# First hook for propagation kernels
175 * -# Update dependent particles and properties (RATTLE, virtual sites)
176 * -# Calculate forces for the current state of the system. This includes
177 * forces added by the Langevin thermostat and the
178 * Lattice-Boltzmann-particle coupling
179 * -# Second hook for propagation kernels
180 * -# Update dependent properties (Virtual sites, RATTLE)
181 * -# Run single step algorithms (Lattice-Boltzmann propagation, collision
182 * detection, NpT update)
183 * - Final update of dependent properties and statistics/counters
184 *
185 * High-level documentation of the integration and thermostatting schemes
186 * can be found in doc/sphinx/system_setup.rst and /doc/sphinx/running.rst
187 *
188 * @return number of steps that have been integrated, or a negative error
189 * code
190 */
191 int integrate(int n_steps, int reuse_forces);
192
193 int integrate_with_signal_handler(int n_steps, int reuse_forces,
194 bool update_accumulators);
195
196 /** @brief Calculate initial particle forces from active thermostats. */
197 void thermostat_force_init();
198 /** @brief Calculate particle-lattice interactions. */
199 void lb_couple_particles(double time_step);
200
201 /** \name Hook procedures
202 * These procedures are called if several significant changes to
203 * the system happen which may make a reinitialization of subsystems
204 * necessary.
205 */
206 /**@{*/
207 /**
208 * @brief Called when the box length has changed. This routine is relatively
209 * fast, and changing the box length every time step as for example necessary
210 * for NpT is more or less ok.
211 *
212 * @param skip_method_adaption skip the long-range methods adaptions
213 */
214 void on_boxl_change(bool skip_method_adaption = false);
215 void on_node_grid_change();
216 void on_periodicity_change();
217 void on_cell_structure_change();
218 void on_thermostat_param_change();
219 void on_temperature_change();
220 void on_verlet_skin_change();
221 void on_timestep_change();
222 void on_integration_start();
223 void on_short_range_ia_change();
224 void on_non_bonded_ia_change();
225 void on_coulomb_change();
226 void on_dipoles_change();
227 /** @brief Called every time a constraint is changed. */
228 void on_constraint_change();
229 /** @brief Called when the LB boundary conditions change
230 * (geometry, slip velocity, or both).
231 */
232 void on_lb_boundary_conditions_change();
233 /** @brief Called every time a particle local property changes. */
234 void on_particle_local_change();
235 /** @brief Called every time a particle property changes. */
236 void on_particle_change();
237 /** @brief Called every time a particle charge changes. */
238 void on_particle_charge_change();
239 /** called before calculating observables, i.e. energy, pressure or
240 * the integrator (forces). Initialize any methods here which are not
241 * initialized immediately (P3M etc.).
242 */
243 void on_observable_calc();
244 void veto_boxl_change(bool skip_particle_checks = false) const;
245 /**@}*/
246
247 /**
248 * @brief Update particles with properties depending on other particles,
249 * namely virtual sites and ICC charges.
250 */
251 void update_dependent_particles();
252 /**
253 * @brief Update the global propagation bitmask.
254 */
255 void update_used_propagations();
256 /**
257 * @brief Veto temperature change.
258 */
259 void check_kT(double value) const;
260
265 std::shared_ptr<BoxGeometry> box_geo;
266 std::shared_ptr<LocalBox> local_geo;
267 std::shared_ptr<CellStructure> cell_structure;
268 std::shared_ptr<Propagation> propagation;
269 std::shared_ptr<InteractionsNonBonded> nonbonded_ias;
270 std::shared_ptr<Thermostat::Thermostat> thermostat;
271 std::shared_ptr<ComFixed> comfixed;
272 std::shared_ptr<Galilei> galilei;
273 std::shared_ptr<BondBreakage::BondBreakage> bond_breakage;
274 std::shared_ptr<LeesEdwards::LeesEdwards> lees_edwards;
275
276protected:
277 /** @brief Whether the thermostat has to be reinitialized before integration.
278 */
280 /** @brief Molecular dynamics integrator time step. */
281 double time_step;
282 /** @brief Molecular dynamics integrator simulation time. */
283 double sim_time;
284 /** @brief Molecular dynamics integrator force capping. */
285 double force_cap;
286 /**
287 * @brief Minimal global interaction cutoff.
288 * Particles with a distance smaller than this are guaranteed
289 * to be available on the same node (through ghosts).
290 */
292
293 void update_local_geo();
294#ifdef ELECTROSTATICS
295 void update_icc_particles();
296#endif // ELECTROSTATICS
297
298private:
299 /**
300 * @brief Check integrator parameters and incompatibilities between
301 * the integrator and the currently active thermostat(s).
302 */
303 void integrator_sanity_checks() const;
304};
305
307void set_system(std::shared_ptr<System> new_instance);
308void reset_system();
309bool is_system_set();
310
311} // namespace System
CellStructureType
Cell structure topology.
Vector implementation and trait types for boost qvm interoperability.
Particle data communication manager for the GPU.
Observable for the pressure and energy.
Queue to deallocate resources before normal program termination.
auto get_time_step() const
Get time_step.
std::shared_ptr< LeesEdwards::LeesEdwards > lees_edwards
double sim_time
Molecular dynamics integrator simulation time.
std::shared_ptr< LocalBox > local_geo
double min_global_cut
Minimal global interaction cutoff.
auto get_min_global_cut() const
Get min_global_cut.
GpuParticleData gpu
bool reinit_thermo
Whether the thermostat has to be reinitialized before integration.
std::shared_ptr< ComFixed > comfixed
std::shared_ptr< BondBreakage::BondBreakage > bond_breakage
Dipoles::Solver dipoles
double force_cap
Molecular dynamics integrator force capping.
ResourceCleanup cleanup_queue
std::shared_ptr< Propagation > propagation
double time_step
Molecular dynamics integrator time step.
std::shared_ptr< Thermostat::Thermostat > thermostat
std::shared_ptr< CellStructure > cell_structure
auto get_force_cap() const
Get force_cap.
auto get_sim_time() const
Get sim_time.
Coulomb::Solver coulomb
std::shared_ptr< BoxGeometry > box_geo
std::shared_ptr< Galilei > galilei
std::shared_ptr< InteractionsNonBonded > nonbonded_ias
This file contains the defaults for ESPResSo.
void reset_system()
System & get_system()
bool is_system_set()
void set_system(std::shared_ptr< System > new_instance)
Describes a cell structure / cell system.