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-2026 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 "ResourceCleanup.hpp"
25
28
29#include "ek/Solver.hpp"
30#include "lb/Solver.hpp"
31
33
34#include <utils/Vector.hpp>
35
36#include <memory>
37#include <optional>
38#include <vector>
39
40class BoxGeometry;
41class LocalBox;
42class CellStructure;
43class Propagation;
46namespace Thermostat {
47class Thermostat;
48}
49class ComFixed;
50class Galilei;
51class Observable_stat;
52class OifGlobal;
54namespace CollisionDetection {
56}
57namespace BondBreakage {
58class BondBreakage;
59}
60namespace LeesEdwards {
61class LeesEdwards;
62}
63namespace Accumulators {
64class AutoUpdateAccumulators;
65}
66namespace Constraints {
67class Constraints;
68}
69struct SteepestDescent;
71struct NptIsoParameters;
73#ifdef ESPRESSO_CUDA
74class GpuParticleData;
75#endif
76
77namespace System {
78
79/**
80 * @brief Main system class.
81 *
82 * Most components follow the composite pattern and the opaque pointer pattern.
83 * See @ref SystemClassDesign for more details.
84 */
85class System : public std::enable_shared_from_this<System> {
86private:
87 struct Private {};
88 void initialize();
89
90public:
91 System(Private);
92
93 static std::shared_ptr<System> create();
94
95#ifdef ESPRESSO_CUDA
96 std::shared_ptr<GpuParticleData> gpu;
97#endif
99
100 /** @brief Get @ref time_step. */
101 auto get_time_step() const { return time_step; }
102
103 /** @brief Set @ref time_step. */
104 void set_time_step(double value);
105
106 /** @brief Get @ref sim_time. */
107 auto get_sim_time() const { return sim_time; }
108
109 /** @brief Set @ref sim_time. */
110 void set_sim_time(double value);
111
112 /** @brief Get @ref force_cap. */
113 auto get_force_cap() const { return force_cap; }
114
115 /** @brief Set @ref force_cap. */
116 void set_force_cap(double value);
117
118 /** @brief Get @ref min_global_cut. */
119 auto get_min_global_cut() const { return min_global_cut; }
120
121 /** @brief Set @ref min_global_cut. */
122 void set_min_global_cut(double value);
123
124 /** @brief Change the box dimensions. */
125 void set_box_l(Utils::Vector3d const &box_l);
126
127 /**
128 * @brief Tune the Verlet skin.
129 * Choose the optimal Verlet list skin between @p min_skin and @p max_skin
130 * by bisection to tolerance @p tol.
131 */
132 void tune_verlet_skin(double min_skin, double max_skin, double tol,
133 int int_steps, bool adjust_max_skin);
134
135 /** @brief Change cell structure topology. */
136 void set_cell_structure_topology(CellStructureType topology);
137
138 /** @brief Rebuild cell lists. Use e.g. after a skin change. */
139 void rebuild_cell_structure();
140 void rebuild_aosoa();
141
142 /** @brief Calculate the maximal cutoff of all interactions. */
143 double maximal_cutoff() const;
144
145 /** @brief Get the interaction range. */
146 double get_interaction_range() const;
147
148 unsigned get_global_ghost_flags() const;
149
150 /** Check electrostatic and magnetostatic methods are properly initialized.
151 * @return true if sanity checks failed.
152 */
153 bool long_range_interactions_sanity_checks() const;
154
155 /** @brief Calculate the total energy. */
156 std::shared_ptr<Observable_stat> calculate_energy();
157
158 /** @brief Calculate the pressure from a virial expansion. */
159 std::shared_ptr<Observable_stat> calculate_pressure();
160
161#ifdef ESPRESSO_NPT
162 /** @brief get the instantaneous pressure with (q(t+dt), p(t+dt/2))*/
164
165 /** @brief get the instantaneous virial pressure with q(t+dt)*/
167
168 /** @brief Synchronize NpT state such as instantaneous and average pressure */
169 void synchronize_npt_state();
170 /** @brief Reinitialize the NpT state. */
171 void npt_ensemble_init(bool recalc_forces);
172 void npt_add_virial_contribution(double energy);
173 bool has_npt_enabled() const;
174#endif // ESPRESSO_NPT
175 Utils::Vector3d *get_npt_virial() const;
176
177 /** @brief Calculate all forces. */
178 void calculate_forces();
179
180#ifdef ESPRESSO_DIPOLE_FIELD_TRACKING
181 /** @brief Calculate dipole fields. */
183#endif
184
185#ifdef ESPRESSO_COLLISION_DETECTION
186 bool has_collision_detection_enabled() const;
187#endif
188
189 /**
190 * @brief Compute the short-range energy of a particle.
191 *
192 * Iterate through particles inside cell and neighboring cells and compute
193 * energy contribution for a specific particle.
194 *
195 * @param pid Particle id
196 * @return Non-bonded energy of the particle.
197 */
198 double particle_short_range_energy_contribution(int pid);
199 /**
200 * @brief Compute the energy of a given bond which has to exist on the given
201 * particle.
202 *
203 * Requires that bond partners are visible on the same MPI rank as the
204 * primary particle.
205 * Returns nothing if the primary particle is not owned by this MPI rank.
206 *
207 * @param pid Particle id
208 * @param bond_id Bond id
209 * @param partners Particle ids of the bond partners
210 *
211 * @return energy of the bond given the primary particle and bond partners
212 */
213 std::optional<double> particle_bond_energy(int pid, int bond_id,
214 std::vector<int> partners);
215
216 /** Integrate equations of motion
217 * @param n_steps Number of integration steps, can be zero
218 * @param reuse_forces Decide when to re-calculate forces
219 *
220 * @details This function calls two hooks for propagation kernels such as
221 * velocity Verlet, velocity Verlet + NpT, or steepest descent.
222 * One hook is called before and one after the force calculation.
223 * It is up to the propagation kernels to increment the simulation time.
224 *
225 * This function propagates the system according to the choice of integrator
226 * stored in @ref Propagation::integ_switch. The general structure is:
227 * - if reuse_forces is zero, recalculate the forces based on the current
228 * state of the system
229 * - Loop over the number of simulation steps:
230 * -# initialization (e.g., RATTLE)
231 * -# First hook for propagation kernels
232 * -# Update dependent particles and properties (RATTLE, virtual sites)
233 * -# Calculate forces for the current state of the system. This includes
234 * forces added by the Langevin thermostat and the
235 * Lattice-Boltzmann-particle coupling
236 * -# Second hook for propagation kernels
237 * -# Update dependent properties (Virtual sites, RATTLE)
238 * -# Run single step algorithms (Lattice-Boltzmann propagation, collision
239 * detection, NpT update)
240 * - Final update of dependent properties and statistics/counters
241 *
242 * High-level documentation of the integration and thermostatting schemes
243 * can be found in doc/sphinx/system_setup.rst and /doc/sphinx/running.rst
244 *
245 * @return number of steps that have been integrated, or a negative error
246 * code
247 */
248 int integrate(int n_steps, int reuse_forces);
249
250 int integrate_with_signal_handler(int n_steps, int reuse_forces,
251 bool update_accumulators);
252
253 /** @brief Calculate particle-lattice interactions. */
254 void lb_couple_particles();
255
256 /** \name Hook procedures
257 * These procedures are called if several significant changes to
258 * the system happen which may make a reinitialization of subsystems
259 * necessary.
260 */
261 /**@{*/
262 /**
263 * @brief Called when the box length has changed. This routine is relatively
264 * fast, and changing the box length every time step as for example necessary
265 * for NpT is more or less ok.
266 *
267 * @param skip_method_adaption skip the long-range methods adaptions
268 */
269 void on_boxl_change(bool skip_method_adaption = false);
270 void on_node_grid_change();
271 void on_periodicity_change();
272 void on_cell_structure_change();
273 void on_thermostat_param_change();
274 void on_temperature_change();
275 void on_verlet_skin_change();
276 void on_timestep_change();
277 void on_integration_start();
278 void on_short_range_ia_change();
279 void on_non_bonded_ia_change();
280 void on_coulomb_change();
281 void on_dipoles_change();
282 /** @brief Called every time a constraint is changed. */
283 void on_constraint_change();
284 /** @brief Called when the LB boundary conditions change
285 * (geometry, slip velocity, or both).
286 */
287 void on_lb_boundary_conditions_change();
288 /** @brief Called every time a particle local property changes. */
289 void on_particle_local_change();
290 /** @brief Called every time a particle property changes. */
291 void on_particle_change();
292 /** @brief Called every time a particle charge changes. */
293 void on_particle_charge_change();
294 /** called before calculating observables, i.e. energy, pressure or
295 * the integrator (forces). Initialize any methods here which are not
296 * initialized immediately (P3M etc.).
297 */
298 void on_observable_calc();
299 void on_lees_edwards_change();
300 void veto_boxl_change(bool skip_particle_checks = false) const;
301 /**@}*/
302
303 /**
304 * @brief Update particles with properties depending on other particles,
305 * namely virtual sites and ICC charges.
306 */
307 void update_dependent_particles();
308 /**
309 * @brief Update the global propagation bitmask.
310 */
311 void update_used_propagations();
312 /**
313 * @brief Veto temperature change.
314 */
315 void check_kT(double value) const;
316
321 std::shared_ptr<BoxGeometry> box_geo;
322 std::shared_ptr<LocalBox> local_geo;
323 std::shared_ptr<CellStructure> cell_structure;
324 std::shared_ptr<Propagation> propagation;
325 std::shared_ptr<BondedInteractionsMap> bonded_ias;
326 std::shared_ptr<InteractionsNonBonded> nonbonded_ias;
327 std::shared_ptr<Thermostat::Thermostat> thermostat;
328 std::shared_ptr<ComFixed> comfixed;
329 std::shared_ptr<Galilei> galilei;
330 std::shared_ptr<OifGlobal> oif_global;
331 std::shared_ptr<ImmersedBoundaries> immersed_boundaries;
332#ifdef ESPRESSO_COLLISION_DETECTION
333 std::shared_ptr<CollisionDetection::CollisionDetection> collision_detection;
334#endif
335 std::shared_ptr<BondBreakage::BondBreakage> bond_breakage;
336 std::shared_ptr<LeesEdwards::LeesEdwards> lees_edwards;
337 std::shared_ptr<Accumulators::AutoUpdateAccumulators>
339 std::shared_ptr<Constraints::Constraints> constraints;
340 std::shared_ptr<SteepestDescent> steepest_descent;
341#ifdef ESPRESSO_STOKESIAN_DYNAMICS
342 std::shared_ptr<StokesianDynamics> stokesian_dynamics;
343#endif
344#ifdef ESPRESSO_NPT
345 std::shared_ptr<NptIsoParameters> nptiso;
346 std::shared_ptr<InstantaneousPressure> npt_inst_pressure;
347#endif
348
349protected:
350 /** @brief Whether the thermostat has to be reinitialized before integration.
351 */
353 /** @brief Molecular dynamics integrator time step. */
354 double time_step;
355 /** @brief Molecular dynamics integrator simulation time. */
356 double sim_time;
357 /** @brief Molecular dynamics integrator force capping. */
358 double force_cap;
359 /**
360 * @brief Minimal global interaction cutoff.
361 * Particles with a distance smaller than this are guaranteed
362 * to be available on the same node (through ghosts).
363 */
365
366 void update_local_geo();
367#ifdef ESPRESSO_ELECTROSTATICS
368 void update_icc_particles();
369 bool has_icc_enabled() const;
370#endif // ESPRESSO_ELECTROSTATICS
371#ifdef ESPRESSO_THERMAL_STONER_WOHLFARTH
372 void integrate_magnetodynamics();
373#endif
374
375private:
376 /**
377 * @brief Check integrator parameters and incompatibilities between
378 * the integrator and the currently active thermostat(s).
379 */
380 void integrator_sanity_checks() const;
381};
382
384void set_system(std::shared_ptr<System> new_instance);
385void reset_system();
386
387} // namespace System
CellStructureType
Cell structure topology.
Vector implementation and trait types for boost qvm interoperability.
container for bonded interactions.
Describes a cell structure / cell system.
Particle data communication manager for the GPU.
Observable for the pressure and energy.
Queue to deallocate resources before normal program termination.
double get_instantaneous_pressure_virial()
get the instantaneous virial pressure with q(t+dt)
auto get_time_step() const
Get time_step.
std::shared_ptr< StokesianDynamics > stokesian_dynamics
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.
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
double get_instantaneous_pressure()
get the instantaneous pressure with (q(t+dt), p(t+dt/2))
std::shared_ptr< InstantaneousPressure > npt_inst_pressure
std::shared_ptr< OifGlobal > oif_global
std::shared_ptr< NptIsoParameters > nptiso
double force_cap
Molecular dynamics integrator force capping.
ResourceCleanup cleanup_queue
std::shared_ptr< Propagation > propagation
double time_step
Molecular dynamics integrator time step.
void calculate_long_range_fields()
Calculate dipole fields.
std::shared_ptr< SteepestDescent > steepest_descent
std::shared_ptr< GpuParticleData > gpu
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
System & get_system()
void set_system(std::shared_ptr< System > new_instance)
Instantaneous pressure during force calculation for NPT integration.
Definition npt.hpp:92
Parameters of the isotropic NpT-integration scheme.
Definition npt.hpp:42
Steepest descent algorithm.