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