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;
45namespace Thermostat {
46class Thermostat;
47}
48class ComFixed;
49class Galilei;
50class Observable_stat;
51class OifGlobal;
53namespace CollisionDetection {
55}
56namespace BondBreakage {
57class BondBreakage;
58}
59namespace LeesEdwards {
60class LeesEdwards;
61}
62namespace Accumulators {
63class AutoUpdateAccumulators;
64}
65namespace Constraints {
66class Constraints;
67}
68
69namespace System {
70
71/**
72 * @brief Main system class.
73 *
74 * Most components follow the composite pattern and the opaque pointer pattern.
75 * See @ref SystemClassDesign for more details.
76 */
77class System : public std::enable_shared_from_this<System> {
78private:
79 struct Private {};
80 void initialize();
81
82public:
83 System(Private);
84
85 static std::shared_ptr<System> create();
86
87#ifdef CUDA
89#endif
91
92 /** @brief Get @ref time_step. */
93 auto get_time_step() const { return time_step; }
94
95 /** @brief Set @ref time_step. */
96 void set_time_step(double value);
97
98 /** @brief Get @ref sim_time. */
99 auto get_sim_time() const { return sim_time; }
100
101 /** @brief Set @ref sim_time. */
102 void set_sim_time(double value);
103
104 /** @brief Get @ref force_cap. */
105 auto get_force_cap() const { return force_cap; }
106
107 /** @brief Set @ref force_cap. */
108 void set_force_cap(double value);
109
110 /** @brief Get @ref min_global_cut. */
111 auto get_min_global_cut() const { return min_global_cut; }
112
113 /** @brief Set @ref min_global_cut. */
114 void set_min_global_cut(double value);
115
116 /** @brief Change the box dimensions. */
117 void set_box_l(Utils::Vector3d const &box_l);
118
119 /**
120 * @brief Tune the Verlet skin.
121 * Choose the optimal Verlet list skin between @p min_skin and @p max_skin
122 * by bisection to tolerance @p tol.
123 */
124 void tune_verlet_skin(double min_skin, double max_skin, double tol,
125 int int_steps, bool adjust_max_skin);
126
127 /** @brief Change cell structure topology. */
128 void set_cell_structure_topology(CellStructureType topology);
129
130 /** @brief Rebuild cell lists. Use e.g. after a skin change. */
131 void rebuild_cell_structure();
132
133 /** @brief Calculate the maximal cutoff of all interactions. */
134 double maximal_cutoff() const;
135
136 /** @brief Get the interaction range. */
137 double get_interaction_range() const;
138
139 unsigned get_global_ghost_flags() const;
140
141 /** Check electrostatic and magnetostatic methods are properly initialized.
142 * @return true if sanity checks failed.
143 */
144 bool long_range_interactions_sanity_checks() const;
145
146 /** @brief Calculate the total energy. */
147 std::shared_ptr<Observable_stat> calculate_energy();
148
149 /** @brief Calculate the pressure from a virial expansion. */
150 std::shared_ptr<Observable_stat> calculate_pressure();
151
152 /** @brief Calculate all forces. */
153 void calculate_forces();
154
155#ifdef DIPOLE_FIELD_TRACKING
156 /** @brief Calculate dipole fields. */
157 void calculate_long_range_fields();
158#endif
159
160 /**
161 * @brief Compute the short-range energy of a particle.
162 *
163 * Iterate through particles inside cell and neighboring cells and compute
164 * energy contribution for a specific particle.
165 *
166 * @param pid Particle id
167 * @return Non-bonded energy of the particle.
168 */
169 double particle_short_range_energy_contribution(int pid);
170
171 /** Integrate equations of motion
172 * @param n_steps Number of integration steps, can be zero
173 * @param reuse_forces Decide when to re-calculate forces
174 *
175 * @details This function calls two hooks for propagation kernels such as
176 * velocity verlet, velocity verlet + npt box changes, and steepest_descent.
177 * One hook is called before and one after the force calculation.
178 * It is up to the propagation kernels to increment the simulation time.
179 *
180 * This function propagates the system according to the choice of integrator
181 * stored in @ref Propagation::integ_switch. The general structure is:
182 * - if reuse_forces is zero, recalculate the forces based on the current
183 * state of the system
184 * - Loop over the number of simulation steps:
185 * -# initialization (e.g., RATTLE)
186 * -# First hook for propagation kernels
187 * -# Update dependent particles and properties (RATTLE, virtual sites)
188 * -# Calculate forces for the current state of the system. This includes
189 * forces added by the Langevin thermostat and the
190 * Lattice-Boltzmann-particle coupling
191 * -# Second hook for propagation kernels
192 * -# Update dependent properties (Virtual sites, RATTLE)
193 * -# Run single step algorithms (Lattice-Boltzmann propagation, collision
194 * detection, NpT update)
195 * - Final update of dependent properties and statistics/counters
196 *
197 * High-level documentation of the integration and thermostatting schemes
198 * can be found in doc/sphinx/system_setup.rst and /doc/sphinx/running.rst
199 *
200 * @return number of steps that have been integrated, or a negative error
201 * code
202 */
203 int integrate(int n_steps, int reuse_forces);
204
205 int integrate_with_signal_handler(int n_steps, int reuse_forces,
206 bool update_accumulators);
207
208 /** @brief Calculate initial particle forces from active thermostats. */
209 void thermostat_force_init();
210 /** @brief Calculate particle-lattice interactions. */
211 void lb_couple_particles();
212
213 /** \name Hook procedures
214 * These procedures are called if several significant changes to
215 * the system happen which may make a reinitialization of subsystems
216 * necessary.
217 */
218 /**@{*/
219 /**
220 * @brief Called when the box length has changed. This routine is relatively
221 * fast, and changing the box length every time step as for example necessary
222 * for NpT is more or less ok.
223 *
224 * @param skip_method_adaption skip the long-range methods adaptions
225 */
226 void on_boxl_change(bool skip_method_adaption = false);
227 void on_node_grid_change();
228 void on_periodicity_change();
229 void on_cell_structure_change();
230 void on_thermostat_param_change();
231 void on_temperature_change();
232 void on_verlet_skin_change();
233 void on_timestep_change();
234 void on_integration_start();
235 void on_short_range_ia_change();
236 void on_non_bonded_ia_change();
237 void on_coulomb_change();
238 void on_dipoles_change();
239 /** @brief Called every time a constraint is changed. */
240 void on_constraint_change();
241 /** @brief Called when the LB boundary conditions change
242 * (geometry, slip velocity, or both).
243 */
244 void on_lb_boundary_conditions_change();
245 /** @brief Called every time a particle local property changes. */
246 void on_particle_local_change();
247 /** @brief Called every time a particle property changes. */
248 void on_particle_change();
249 /** @brief Called every time a particle charge changes. */
250 void on_particle_charge_change();
251 /** called before calculating observables, i.e. energy, pressure or
252 * the integrator (forces). Initialize any methods here which are not
253 * initialized immediately (P3M etc.).
254 */
255 void on_observable_calc();
256 void on_lees_edwards_change();
257 void veto_boxl_change(bool skip_particle_checks = false) const;
258 /**@}*/
259
260 /**
261 * @brief Update particles with properties depending on other particles,
262 * namely virtual sites and ICC charges.
263 */
264 void update_dependent_particles();
265 /**
266 * @brief Update the global propagation bitmask.
267 */
268 void update_used_propagations();
269 /**
270 * @brief Veto temperature change.
271 */
272 void check_kT(double value) const;
273
278 std::shared_ptr<BoxGeometry> box_geo;
279 std::shared_ptr<LocalBox> local_geo;
280 std::shared_ptr<CellStructure> cell_structure;
281 std::shared_ptr<Propagation> propagation;
282 std::shared_ptr<BondedInteractionsMap> bonded_ias;
283 std::shared_ptr<InteractionsNonBonded> nonbonded_ias;
284 std::shared_ptr<Thermostat::Thermostat> thermostat;
285 std::shared_ptr<ComFixed> comfixed;
286 std::shared_ptr<Galilei> galilei;
287 std::shared_ptr<OifGlobal> oif_global;
288 std::shared_ptr<ImmersedBoundaries> immersed_boundaries;
289#ifdef COLLISION_DETECTION
290 std::shared_ptr<CollisionDetection::CollisionDetection> collision_detection;
291#endif
292 std::shared_ptr<BondBreakage::BondBreakage> bond_breakage;
293 std::shared_ptr<LeesEdwards::LeesEdwards> lees_edwards;
294 std::shared_ptr<Accumulators::AutoUpdateAccumulators>
296 std::shared_ptr<Constraints::Constraints> constraints;
297
298protected:
299 /** @brief Whether the thermostat has to be reinitialized before integration.
300 */
302 /** @brief Molecular dynamics integrator time step. */
303 double time_step;
304 /** @brief Molecular dynamics integrator simulation time. */
305 double sim_time;
306 /** @brief Molecular dynamics integrator force capping. */
307 double force_cap;
308 /**
309 * @brief Minimal global interaction cutoff.
310 * Particles with a distance smaller than this are guaranteed
311 * to be available on the same node (through ghosts).
312 */
314
315 void update_local_geo();
316#ifdef ELECTROSTATICS
317 void update_icc_particles();
318#endif // ELECTROSTATICS
319
320private:
321 /**
322 * @brief Check integrator parameters and incompatibilities between
323 * the integrator and the currently active thermostat(s).
324 */
325 void integrator_sanity_checks() const;
326};
327
329void set_system(std::shared_ptr<System> new_instance);
330void reset_system();
331
332} // namespace System
CellStructureType
Cell structure topology.
Vector implementation and trait types for boost qvm interoperability.
container for bonded interactions.
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
std::shared_ptr< BondedInteractionsMap > bonded_ias
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< ImmersedBoundaries > immersed_boundaries
std::shared_ptr< Accumulators::AutoUpdateAccumulators > auto_update_accumulators
std::shared_ptr< ComFixed > comfixed
std::shared_ptr< BondBreakage::BondBreakage > bond_breakage
Dipoles::Solver dipoles
std::shared_ptr< CollisionDetection::CollisionDetection > collision_detection
std::shared_ptr< OifGlobal > oif_global
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
std::shared_ptr< Constraints::Constraints > constraints
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.
System & get_system()
void set_system(std::shared_ptr< System > new_instance)
Describes a cell structure / cell system.