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
script_interface/system/System.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2013-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#include "System.hpp"
21
22#include "config/config.hpp"
23
24#include "core/BoxGeometry.hpp"
25#include "core/Particle.hpp"
29#include "core/cells.hpp"
31#include "core/exclusions.hpp"
33#include "core/npt.hpp"
36#include "core/propagation.hpp"
37#include "core/rotation.hpp"
40
59
60#include <utils/Vector.hpp>
61#include <utils/demangle.hpp>
64
65#include <boost/mpi/collectives.hpp>
66
67#include <algorithm>
68#include <array>
69#include <cassert>
70#include <cmath>
71#include <functional>
72#include <initializer_list>
73#include <memory>
74#include <stdexcept>
75#include <string>
76#include <type_traits>
77#include <vector>
78
79namespace ScriptInterface {
80namespace System {
81
82static bool system_created = false;
83
84#ifdef EXCLUSIONS
86 Variant const &exclusions) {
87 p.call_method("set_exclusions", {{"p_ids", exclusions}});
88}
89#endif // EXCLUSIONS
90
91static void set_bonds(Particles::ParticleHandle &p, Variant const &bonds) {
93 for (auto const &bond_flat : bond_list_flat) {
94 auto const bond_id = bond_flat[0];
95 auto const part_id =
96 std::vector<int>{bond_flat.begin() + 1, bond_flat.end()};
97 p.call_method("add_bond",
98 {{"bond_id", bond_id}, {"part_id", std::move(part_id)}});
99 }
100}
101
102/** @brief Container for leaves of the system class. */
104 Leaves() = default;
105 std::shared_ptr<CellSystem::CellSystem> cell_system;
106 std::shared_ptr<Integrators::IntegratorHandle> integrator;
107 std::shared_ptr<Interactions::BondedInteractions> bonded_interactions;
108#ifdef COLLISION_DETECTION
109 std::shared_ptr<CollisionDetection::CollisionDetection> collision_detection;
110#endif
111 std::shared_ptr<Thermostat::Thermostat> thermostat;
112 std::shared_ptr<Analysis::Analysis> analysis;
113 std::shared_ptr<Galilei::ComFixed> comfixed;
114 std::shared_ptr<Galilei::Galilei> galilei;
115 std::shared_ptr<BondBreakage::BreakageSpecs> bond_breakage;
116 std::shared_ptr<LeesEdwards::LeesEdwards> lees_edwards;
117 std::shared_ptr<Accumulators::AutoUpdateAccumulators>
119 std::shared_ptr<Constraints::Constraints> constraints;
120 std::shared_ptr<Interactions::NonBondedInteractions> non_bonded_inter;
121#ifdef ELECTROSTATICS
122 std::shared_ptr<Coulomb::Container> electrostatics;
123#endif
124#ifdef DIPOLES
125 std::shared_ptr<Dipoles::Container> magnetostatics;
126#endif
127 std::shared_ptr<Particles::ParticleList> part;
128
130 // Clear containers whose elements call MPI callbacks upon destruction.
131 // The containers lifetime is extended by the global context shared object
132 // registry on the head node, which can cause MPI deadlocks if they still
133 // contain elements.
135 bonded_interactions->clear();
136 }
137 if (constraints) {
138 constraints->clear();
139 }
140 }
141};
142
143System::System() : m_instance{}, m_leaves{std::make_unique<Leaves>()} {
144 auto const add_parameter =
145 [this, ptr = m_leaves.get()](std::string key, auto(Leaves::*member)) {
147 key.c_str(),
148 [this, ptr, member, key](Variant const &val) {
149 auto &dst = ptr->*member;
150 if (dst != nullptr) {
151 throw WriteError(key);
152 }
153 dst = get_value<std::remove_reference_t<decltype(dst)>>(val);
154 dst->bind_system(m_instance);
155 },
156 [ptr, member]() { return ptr->*member; })});
157 };
158
159 add_parameters({
160 {"box_l",
161 [this](Variant const &v) {
162 context()->parallel_try_catch([&]() {
165 throw std::domain_error("Attribute 'box_l' must be > 0");
166 }
167 m_instance->veto_boxl_change();
168 m_instance->box_geo->set_length(new_value);
169 m_instance->on_boxl_change();
170 });
171 },
172 [this]() { return m_instance->box_geo->length(); }},
173 {"periodicity",
174 [this](Variant const &v) {
176 for (unsigned int i = 0u; i < 3u; ++i) {
177 m_instance->box_geo->set_periodic(i, periodicity[i]);
178 }
179 context()->parallel_try_catch(
180 [&]() { m_instance->on_periodicity_change(); });
181 },
182 [this]() {
183 return Utils::Vector3b{m_instance->box_geo->periodic(0),
184 m_instance->box_geo->periodic(1),
185 m_instance->box_geo->periodic(2)};
186 }},
187 {"min_global_cut",
188 [this](Variant const &v) {
189 context()->parallel_try_catch([&]() {
190 auto const new_value = get_value<double>(v);
192 throw std::domain_error("Attribute 'min_global_cut' must be >= 0");
193 }
194 m_instance->set_min_global_cut(new_value);
195 });
196 },
197 [this]() { return m_instance->get_min_global_cut(); }},
198 {"max_oif_objects",
199 [this](Variant const &v) {
200 m_instance->oif_global->max_oif_objects = get_value<int>(v);
201 },
202 [this]() { return m_instance->oif_global->max_oif_objects; }},
203 });
204 // note: the order of leaves matters! e.g. bonds depend on thermostats,
205 // and thus a thermostat object must be instantiated before the bonds
206 add_parameter("cell_system", &Leaves::cell_system);
207 add_parameter("integrator", &Leaves::integrator);
208 add_parameter("thermostat", &Leaves::thermostat);
209 add_parameter("analysis", &Leaves::analysis);
210 add_parameter("comfixed", &Leaves::comfixed);
211 add_parameter("galilei", &Leaves::galilei);
212 add_parameter("bonded_inter", &Leaves::bonded_interactions);
213#ifdef COLLISION_DETECTION
214 add_parameter("collision_detection", &Leaves::collision_detection);
215#endif
216 add_parameter("bond_breakage", &Leaves::bond_breakage);
217 add_parameter("lees_edwards", &Leaves::lees_edwards);
218 add_parameter("auto_update_accumulators", &Leaves::auto_update_accumulators);
219 add_parameter("constraints", &Leaves::constraints);
220 add_parameter("non_bonded_inter", &Leaves::non_bonded_inter);
221#ifdef ELECTROSTATICS
222 add_parameter("electrostatics", &Leaves::electrostatics);
223#endif
224#ifdef DIPOLES
225 add_parameter("magnetostatics", &Leaves::magnetostatics);
226#endif
227 add_parameter("part", &Leaves::part);
228}
229
230template <typename LeafType>
231void System::do_set_default_parameter(std::string const &name) {
232 assert(context()->is_head_node());
233 auto const so_name = Utils::demangle<LeafType>().substr(17);
234 set_parameter(name, Variant{context()->make_shared(so_name, {})});
235}
236
237void System::do_construct(VariantMap const &params) {
238 /* When reloading the system state from a checkpoint file,
239 * the order of global variables instantiation matters.
240 * The @c box_l must be set before any other global variable.
241 * All these globals re-initialize the cell system, and we
242 * cannot use the default-constructed @c box_geo when e.g.
243 * long-range interactions exist in the system, otherwise
244 * runtime errors about the local geometry being smaller
245 * than the interaction range would be raised.
246 */
247 context()->parallel_try_catch([&]() {
248 if (not params.contains("box_l")) {
249 throw std::domain_error("Required argument 'box_l' not provided.");
250 }
251 if (params.contains("_regular_constructor") and system_created) {
252 throw std::runtime_error(
253 "You can only have one instance of the system class at a time");
254 }
255 });
256 m_instance = ::System::System::create();
257 ::System::set_system(m_instance);
258
259 // domain decomposition can only be set after box_l is set
260 m_instance->set_cell_structure_topology(CellStructureType::NSQUARE);
261 do_set_parameter("box_l", params.at("box_l"));
262 m_instance->set_cell_structure_topology(CellStructureType::REGULAR);
263
264 m_instance->lb.bind_system(m_instance);
265 m_instance->ek.bind_system(m_instance);
266
267 if (params.contains("_regular_constructor")) {
268 std::set<std::string> const setable_properties = {
269 "box_l", "min_global_cut",
270 "periodicity", "time",
271 "time_step", "force_cap",
272 "max_oif_objects", "_regular_constructor"};
273 for (auto const &kv : params) {
274 if (not setable_properties.contains(kv.first)) {
275 context()->parallel_try_catch([&kv]() {
276 throw std::domain_error(
277 "Property '" + kv.first +
278 "' cannot be set via argument to System class");
279 });
280 }
281 }
282 for (std::string attr :
283 {"min_global_cut", "periodicity", "max_oif_objects"}) {
284 if (params.contains(attr)) {
285 do_set_parameter(attr, params.at(attr));
286 }
287 }
288 if (not context()->is_head_node()) {
289 return;
290 }
291 auto integrator = std::dynamic_pointer_cast<Integrators::IntegratorHandle>(
292 context()->make_shared("Integrators::IntegratorHandle", {}));
293 set_parameter("integrator", integrator);
294 for (std::string attr : {"time", "time_step", "force_cap"}) {
295 if (params.contains(attr)) {
296 integrator->set_parameter(attr, params.at(attr));
297 }
298 }
299 // note: the order of leaves matters! e.g. bonds depend on thermostats,
300 // and thus a thermostat object must be instantiated before the bonds
304#ifdef COLLISION_DETECTION
306 "collision_detection");
307#endif
314 "auto_update_accumulators");
317 "non_bonded_inter");
318#ifdef ELECTROSTATICS
320#endif
321#ifdef DIPOLES
323#endif
325 } else {
326 for (auto const &key : get_parameter_insertion_order()) {
327 if (key != "box_l" and params.contains(key)) {
328 do_set_parameter(key, params.at(key));
329 }
330 }
331 }
332 if (not context()->is_head_node()) {
333 return;
334 }
335 call_method("internal_attach_leaves", {});
336}
337
338static void rotate_system(CellStructure &cell_structure, double phi,
339 double theta, double alpha) {
340 auto const particles = cell_structure.local_particles();
341
342 // Calculate center of mass
344 double local_mass = 0.0;
345
346 for (auto const &p : particles) {
347 if (not p.is_virtual()) {
348 local_com += p.mass() * p.pos();
349 local_mass += p.mass();
350 }
351 }
352
353 auto const total_mass =
354 boost::mpi::all_reduce(comm_cart, local_mass, std::plus<>());
355 auto const com =
356 boost::mpi::all_reduce(comm_cart, local_com, std::plus<>()) / total_mass;
357
358 // Rotation axis in Cartesian coordinates
359 Utils::Vector3d axis;
360 axis[0] = std::sin(theta) * std::cos(phi);
361 axis[1] = std::sin(theta) * std::sin(phi);
362 axis[2] = std::cos(theta);
363
364 // Rotate particle coordinates
365 for (auto &p : particles) {
366 // Move the center of mass of the system to the origin
367 p.pos() = com + Utils::vec_rotate(axis, alpha, p.pos() - com);
368#ifdef ROTATION
369 local_rotate_particle(p, axis, alpha);
370#endif
371 }
372
376}
377
378/** Rescale all particle positions in direction @p dir by a factor @p scale.
379 * @param cell_structure cell structure
380 * @param dir direction to scale (0/1/2 = x/y/z, 3 = x+y+z isotropically)
381 * @param scale factor by which to rescale (>1: stretch, <1: contract)
382 */
383static void rescale_particles(CellStructure &cell_structure, int dir,
384 double scale) {
385 for (auto &p : cell_structure.local_particles()) {
386 if (dir < 3)
387 p.pos()[dir] *= scale;
388 else {
389 p.pos() *= scale;
390 }
391 }
392}
393
394Variant System::do_call_method(std::string const &name,
395 VariantMap const &parameters) {
396 if (name == "lock_system_creation") {
397 system_created = true;
398 return {};
399 }
400 if (name == "rescale_boxl") {
401 auto &box_geo = *m_instance->box_geo;
402 auto const coord = get_value<int>(parameters, "coord");
403 auto const length = get_value<double>(parameters, "length");
404 assert(coord >= 0);
405 assert(coord != 3 or ((box_geo.length()[0] == box_geo.length()[1]) and
406 (box_geo.length()[1] == box_geo.length()[2])));
407 auto const scale = (coord == 3) ? length * box_geo.length_inv()[0]
408 : length * box_geo.length_inv()[coord];
409 context()->parallel_try_catch([&]() {
410 if (length <= 0.) {
411 throw std::domain_error("Parameter 'd_new' must be > 0");
412 }
413 m_instance->veto_boxl_change(true);
414 });
415 auto new_value = Utils::Vector3d{};
416 if (coord == 3) {
418 } else {
419 new_value = box_geo.length();
420 new_value[static_cast<unsigned>(coord)] = length;
421 }
422 // when shrinking, rescale the particles first
423 if (scale <= 1.) {
424 rescale_particles(*m_instance->cell_structure, coord, scale);
425 m_instance->on_particle_change();
426 }
427 m_instance->box_geo->set_length(new_value);
428 m_instance->on_boxl_change();
429 if (scale > 1.) {
430 rescale_particles(*m_instance->cell_structure, coord, scale);
431 m_instance->on_particle_change();
432 }
433 return {};
434 }
435 if (name == "setup_type_map") {
436 auto const types = get_value<std::vector<int>>(parameters, "type_list");
437 for (auto const type : types) {
438 ::init_type_map(type);
439 }
440 return {};
441 }
442 if (name == "number_of_particles") {
443 auto const type = get_value<int>(parameters, "type");
444 return ::number_of_particles_with_type(type);
445 }
446 if (name == "velocity_difference") {
447 auto const pos1 = get_value<Utils::Vector3d>(parameters, "pos1");
448 auto const pos2 = get_value<Utils::Vector3d>(parameters, "pos2");
449 auto const v1 = get_value<Utils::Vector3d>(parameters, "v1");
450 auto const v2 = get_value<Utils::Vector3d>(parameters, "v2");
451 return m_instance->box_geo->velocity_difference(pos2, pos1, v2, v1);
452 }
453 if (name == "distance_vec") {
454 auto const pos1 = get_value<Utils::Vector3d>(parameters, "pos1");
455 auto const pos2 = get_value<Utils::Vector3d>(parameters, "pos2");
456 return m_instance->box_geo->get_mi_vector(pos2, pos1);
457 }
458 if (name == "rotate_system") {
459 rotate_system(*m_instance->cell_structure,
462 get_value<double>(parameters, "alpha"));
463 m_instance->on_particle_change();
464 m_instance->update_dependent_particles();
465 return {};
466 }
467 if (name == "get_propagation_modes_enum") {
469 }
470 if (name == "session_shutdown") {
471 if (m_instance) {
472 if (&::System::get_system() == m_instance.get()) {
474 }
475 assert(m_instance.use_count() == 1l);
476 m_leaves.reset();
477 m_instance.reset();
478 }
479 return {};
480 }
481 if (name == "internal_attach_leaves") {
482 m_leaves->part->attach(m_leaves->cell_system,
483 m_leaves->bonded_interactions);
484#ifdef COLLISION_DETECTION
485 m_leaves->collision_detection->attach(m_leaves->bonded_interactions);
486#endif
487 return {};
488 }
489 return {};
490}
491
492/**
493 * @brief Serialize particles.
494 * Particles need to be serialized here to reduce overhead,
495 * and also to guarantee particles get instantiated after the cell structure
496 * was instantiated (since they store a weak pointer to it).
497 */
498std::string System::get_internal_state() const {
499 auto const p_ids = get_particle_ids();
500 std::vector<std::string> object_states(p_ids.size());
501
502 std::ranges::transform(p_ids, object_states.begin(), [this](auto const p_id) {
503 auto p_obj = context()->make_shared(
504 "Particles::ParticleHandle",
505 {{"id", p_id}, {"__cell_structure", m_leaves->cell_system}});
506 auto &p_handle = dynamic_cast<Particles::ParticleHandle &>(*p_obj);
507 auto const packed_state = p_handle.serialize();
508 // custom particle serialization
509 auto state = Utils::unpack<ObjectState>(packed_state);
510 state.name = "Particles::ParticleHandle";
511 auto const bonds_view = p_handle.call_method("get_bonds_view", {});
512 state.params.emplace_back(std::string{"bonds"}, pack(bonds_view));
513#ifdef EXCLUSIONS
514 auto const exclusions = p_handle.call_method("get_exclusions", {});
515 state.params.emplace_back(std::string{"exclusions"}, pack(exclusions));
516#endif // EXCLUSIONS
517 state.params.emplace_back(std::string{"__cpt_sentinel"}, pack(None{}));
518 return Utils::pack(state);
519 });
520
522}
523
524void System::set_internal_state(std::string const &state) {
525 auto const object_states = Utils::unpack<std::vector<std::string>>(state);
526#ifdef EXCLUSIONS
527 std::unordered_map<int, Variant> exclusions = {};
528#endif // EXCLUSIONS
529 std::unordered_map<int, Variant> bonds = {};
530
531 for (auto const &packed_object : object_states) {
532 auto state = Utils::unpack<ObjectState>(packed_object);
533 VariantMap params = {};
534 for (auto const &kv : state.params) {
535 params[kv.first] = unpack(kv.second, {});
536 }
537 auto const p_id = get_value<int>(params.at("id"));
538 bonds[p_id] = params.extract("bonds").mapped();
539#ifdef EXCLUSIONS
540 exclusions[p_id] = params.extract("exclusions").mapped();
541#endif // EXCLUSIONS
542 params["__cell_structure"] = get_parameter("cell_system");
543 context()->make_shared("Particles::ParticleHandle", params);
544 }
545
546 for (auto const p_id : get_particle_ids()) {
547 auto p_obj = context()->make_shared(
548 "Particles::ParticleHandle",
549 {{"id", p_id}, {"__cell_structure", m_leaves->cell_system}});
550 auto &p_handle = dynamic_cast<Particles::ParticleHandle &>(*p_obj);
551 set_bonds(p_handle, bonds[p_id]);
552#ifdef EXCLUSIONS
553 set_exclusions(p_handle, exclusions[p_id]);
554#endif // EXCLUSIONS
555 }
556}
557
558} // namespace System
559} // namespace ScriptInterface
@ NSQUARE
Atom decomposition (N-square).
@ REGULAR
Regular decomposition.
static int coord(std::string const &s)
Vector implementation and trait types for boost qvm interoperability.
This file contains everything related to the global cell structure / cell system.
void add_parameters(std::vector< AutoParameter > &&params)
Type to indicate no value in Variant.
std::string serialize() const
Variant call_method(const std::string &name, const VariantMap &params)
Call a method on the object.
static std::shared_ptr< System > create()
static DEVICE_QUALIFIER constexpr Vector< T, N > broadcast(typename Base::value_type const &value) noexcept
Create a vector that has all entries set to the same value.
Definition Vector.hpp:111
boost::mpi::communicator comm_cart
The communicator.
This file contains the defaults for ESPResSo.
This file contains the asynchronous MPI communication.
@ DATA_PART_PROPERTIES
Particle::p.
@ DATA_PART_POSITION
Particle::r.
static void set_bonds(Particles::ParticleHandle &p, Variant const &bonds)
static void set_exclusions(Particles::ParticleHandle &p, Variant const &exclusions)
static void rotate_system(CellStructure &cell_structure, double phi, double theta, double alpha)
static void rescale_particles(CellStructure &cell_structure, int dir, double scale)
Rescale all particle positions in direction dir by a factor scale.
PackedVariant pack(const Variant &v)
Transform a Variant to a PackedVariant.
T get_value(Variant const &v)
Extract value of specific type T from a Variant.
std::unordered_map< std::string, Variant > VariantMap
Definition Variant.hpp:69
auto make_unordered_map_of_variants(std::unordered_map< K, V > const &v)
Definition Variant.hpp:80
boost::make_recursive_variant< None, bool, int, std::size_t, double, std::string, ObjectRef, Utils::Vector3b, Utils::Vector3i, Utils::Vector2d, Utils::Vector3d, Utils::Vector4d, std::vector< int >, std::vector< double >, std::vector< boost::recursive_variant_ >, std::unordered_map< int, boost::recursive_variant_ >, std::unordered_map< std::string, boost::recursive_variant_ > >::type Variant
Possible types for parameters.
Definition Variant.hpp:67
Variant unpack(const PackedVariant &v, std::unordered_map< ObjectId, ObjectRef > const &objects)
Unpack a PackedVariant.
System & get_system()
void set_system(std::shared_ptr< System > new_instance)
Vector3d vec_rotate(const Vector3d &axis, double angle, const Vector3d &vector)
Rotate a vector around an axis.
std::string pack(T const &v)
Pack a serialize type into a string.
Definition pack.hpp:38
Various procedures concerning interactions between particles.
constexpr double INACTIVE_CUTOFF
Cutoff for deactivated interactions.
Exports for the NpT code.
Routines to calculate the OIF global forces for a particle triple (triangle from mesh).
void init_type_map(int type)
std::vector< int > get_particle_ids()
Get all particle ids.
Particles creation and deletion.
std::unordered_map< std::string, int > propagation_flags_map()
Convert PropagationMode::PropagationMode to name/value pairs.
This file contains all subroutines required to process rotational motion.
void local_rotate_particle(Particle &p, const Utils::Vector3d &axis_space_frame, const double phi)
Rotate the particle p around the NORMALIZED axis aSpaceFrame by amount phi.
Definition rotation.hpp:133
static SteepestDescentParameters params
Currently active steepest descent instance.
Describes a cell structure / cell system.
void update_ghosts_and_resort_particle(unsigned data_parts)
Update ghost particles, with particle resort if needed.
void set_resort_particles(Cells::Resort level)
Increase the local resort level at least to level.
ParticleRange local_particles() const
Description and getter/setter for a parameter.
Container for leaves of the system class.
std::shared_ptr< Thermostat::Thermostat > thermostat
std::shared_ptr< Integrators::IntegratorHandle > integrator
std::shared_ptr< Dipoles::Container > magnetostatics
std::shared_ptr< CellSystem::CellSystem > cell_system
std::shared_ptr< Particles::ParticleList > part
std::shared_ptr< BondBreakage::BreakageSpecs > bond_breakage
std::shared_ptr< Constraints::Constraints > constraints
std::shared_ptr< Interactions::BondedInteractions > bonded_interactions
std::shared_ptr< CollisionDetection::CollisionDetection > collision_detection
std::shared_ptr< Interactions::NonBondedInteractions > non_bonded_inter
std::shared_ptr< Coulomb::Container > electrostatics
std::shared_ptr< Accumulators::AutoUpdateAccumulators > auto_update_accumulators
std::shared_ptr< LeesEdwards::LeesEdwards > lees_edwards