ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
ParticleHandle.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 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 "config/config.hpp"
21
22#include "ParticleHandle.hpp"
23
29
30#include "core/BoxGeometry.hpp"
33#include "core/bonds.hpp"
35#include "core/exclusions.hpp"
38#include "core/propagation.hpp"
39#include "core/rotation.hpp"
42
43#include <utils/Vector.hpp>
45
46#include <boost/format.hpp>
47#include <boost/mpi/collectives/all_reduce.hpp>
48#include <boost/mpi/collectives/broadcast.hpp>
49#include <boost/mpi/communicator.hpp>
50
51#include <algorithm>
52#include <cassert>
53#include <cmath>
54#include <cstddef>
55#include <memory>
56#include <optional>
57#include <set>
58#include <sstream>
59#include <stdexcept>
60#include <string>
61#include <tuple>
62#include <type_traits>
63#include <vector>
64
65namespace ScriptInterface {
66namespace Particles {
67
68static void particle_checks(int p_id, Utils::Vector3d const &pos) {
69 if (p_id < 0) {
70 throw std::domain_error("Invalid particle id: " + std::to_string(p_id));
71 }
72#ifndef __FAST_MATH__
73 if (std::isnan(pos[0]) or std::isnan(pos[1]) or std::isnan(pos[2]) or
74 std::isinf(pos[0]) or std::isinf(pos[1]) or std::isinf(pos[2])) {
75 throw std::domain_error("Particle position must be finite");
76 }
77#endif // __FAST_MATH__
78}
79
81 auto bitfield = static_cast<uint8_t>(0u);
82 if (flag[0])
83 bitfield |= static_cast<uint8_t>(1u);
84 if (flag[1])
85 bitfield |= static_cast<uint8_t>(2u);
86 if (flag[2])
87 bitfield |= static_cast<uint8_t>(4u);
88 return bitfield;
89}
90
91static auto error_msg(std::string const &name, std::string const &reason) {
92 std::stringstream msg;
93 msg << "attribute '" << name << "' of 'ParticleHandle' " << reason;
94 return msg.str();
95}
96
98 return Utils::Vector4d{{q[0], q[1], q[2], q[3]}};
99}
100
101static auto get_quaternion_safe(std::string const &name, Variant const &value) {
102 auto const q = get_value<Utils::Vector4d>(value);
103 if (q.norm2() == 0.) {
104 throw std::domain_error(error_msg(name, "must be non-zero"));
105 }
106 return Utils::Quaternion<double>{{q[0], q[1], q[2], q[3]}};
107}
108
109#ifdef THERMOSTAT_PER_PARTICLE
110static auto get_gamma_safe(Variant const &value) {
111#ifdef PARTICLE_ANISOTROPY
112 try {
114 } catch (...) {
115 return get_value<Utils::Vector3d>(value);
116 }
117#else // PARTICLE_ANISOTROPY
118 return get_value<double>(value);
119#endif // PARTICLE_ANISOTROPY
120}
121#endif // THERMOSTAT_PER_PARTICLE
122
123static auto get_real_particle(boost::mpi::communicator const &comm, int p_id,
124 ::CellStructure &cell_structure) {
125 if (p_id < 0) {
126 throw std::domain_error("Invalid particle id: " + std::to_string(p_id));
127 }
128 auto ptr = cell_structure.get_local_particle(p_id);
129 if (ptr != nullptr and ptr->is_ghost()) {
130 ptr = nullptr;
131 }
132 auto const n_found = boost::mpi::all_reduce(
133 comm, static_cast<int>(ptr != nullptr), std::plus<>());
134 if (n_found == 0) {
135 throw std::runtime_error("Particle with id " + std::to_string(p_id) +
136 " not found");
137 }
138 return ptr;
139}
140
141template <typename T, class F>
142T ParticleHandle::get_particle_property(F const &fun) const {
143 auto cell_structure_si = get_cell_structure();
144 auto &cell_structure = cell_structure_si->get_cell_structure();
145 auto const &comm = context()->get_comm();
146 auto const ptr = const_cast<Particle const *>(
147 get_real_particle(comm, m_pid, cell_structure));
148 std::optional<T> ret;
149 if (ptr == nullptr) {
150 ret = {};
151 } else {
152 ret = {fun(*ptr)};
153 }
154 return Utils::Mpi::reduce_optional(comm, ret);
155}
156
157template <typename T>
158T ParticleHandle::get_particle_property(T const &(Particle::*getter)()
159 const) const {
160 return get_particle_property<T>(
161 [getter](Particle const &p) { return (p.*getter)(); });
162}
163
164template <class F>
165void ParticleHandle::set_particle_property(F const &fun) const {
166 auto cell_structure_si = get_cell_structure();
167 auto &cell_structure = cell_structure_si->get_cell_structure();
168 auto const &comm = context()->get_comm();
169 auto const ptr = get_real_particle(comm, m_pid, cell_structure);
170 if (ptr != nullptr) {
171 fun(*ptr);
172 }
174}
175
176template <typename T>
177void ParticleHandle::set_particle_property(T &(Particle::*setter)(),
178 Variant const &value) const {
179 set_particle_property(
180 [&value, setter](Particle &p) { (p.*setter)() = get_value<T>(value); });
181}
182
183ParticleHandle::ParticleHandle() {
184 /* Warning: the order of particle property setters matters! Some properties
185 * override each other, e.g. quat/director/dip or dip/dipm.
186 * This is revelant during checkpointing: all particle properties are set at
187 * once, in the order specified in the following call to `add_parameters()`.
188 */
189 add_parameters({
190 {"id", AutoParameter::read_only, [this]() { return m_pid; }},
191 {"type",
192 [this](Variant const &value) {
193 auto const old_type = get_particle_property(&Particle::type);
194 auto const new_type = get_value<int>(value);
195 if (new_type < 0) {
196 throw std::domain_error(
197 error_msg("type", "must be an integer >= 0"));
198 }
199 get_system()->nonbonded_ias->make_particle_type_exist(new_type);
200 on_particle_type_change(m_pid, old_type, new_type);
201 set_particle_property(&Particle::type, value);
202 },
203 [this]() { return get_particle_data(m_pid).type(); }},
204 {"pos",
205 [this](Variant const &value) {
206 auto const pos = get_value<Utils::Vector3d>(value);
207 particle_checks(m_pid, pos);
208 set_particle_pos(m_pid, pos);
209 },
210 [this]() {
211 auto const p = get_particle_data(m_pid);
212 auto const pos = p.pos();
213 auto const image_box = p.image_box();
214 return get_system()->box_geo->unfolded_position(pos, image_box);
215 }},
216 {"v",
217 [this](Variant const &value) {
218 set_particle_property(&Particle::v, value);
219 },
220 [this]() { return get_particle_data(m_pid).v(); }},
221 {"f",
222 [this](Variant const &value) {
223 set_particle_property(&Particle::force, value);
224 },
225 [this]() { return get_particle_data(m_pid).force(); }},
226 {"mass",
227#ifdef MASS
228 [this](Variant const &value) {
229 if (get_value<double>(value) <= 0.) {
230 throw std::domain_error(error_msg("mass", "must be a float > 0"));
231 }
232 set_particle_property(&Particle::mass, value);
233 },
234#else // MASS
235 [](Variant const &value) {
236 auto const default_mass = Particle().mass();
237 if (std::abs(get_value<double>(value) - default_mass) > 1e-10) {
238 throw std::runtime_error("Feature MASS not compiled in");
239 }
240 },
241#endif // MASS
242 [this]() { return get_particle_data(m_pid).mass(); }},
243 {"q",
244#ifdef ELECTROSTATICS
245 [this](Variant const &value) {
246 set_particle_property(&Particle::q, value);
247 },
248#else // ELECTROSTATICS
249 [](Variant const &value) {
250 if (get_value<double>(value) != 0.) {
251 throw std::runtime_error("Feature ELECTROSTATICS not compiled in");
252 }
253 },
254#endif // ELECTROSTATICS
255 [this]() { return get_particle_data(m_pid).q(); }},
256#ifdef DIPOLES
257 {"dip",
258 [this](Variant const &value) {
259 set_particle_property([&value](Particle &p) {
260 auto const dip = get_value<Utils::Vector3d>(value);
261 std::tie(p.quat(), p.dipm()) = convert_dip_to_quat(dip);
262 });
263 },
264 [this]() { return get_particle_data(m_pid).calc_dip(); }},
265 {"dipm",
266 [this](Variant const &value) {
267 set_particle_property(&Particle::dipm, value);
268 },
269 [this]() { return get_particle_data(m_pid).dipm(); }},
270#endif // DIPOLES
271#ifdef DIPOLE_FIELD_TRACKING
272 {"dip_fld",
273 [this](Variant const &value) {
274 set_particle_property(&Particle::dip_fld, value);
275 },
276 [this]() { return get_particle_data(m_pid).dip_fld(); }},
277#endif
278#ifdef ROTATION
279 {"director",
280 [this](Variant const &value) {
281 set_particle_property([&value](Particle &p) {
282 auto const director = get_value<Utils::Vector3d>(value).normalized();
284 });
285 },
286 [this]() {
287 auto const quat = get_particle_data(m_pid).quat();
289 }},
290 {"quat",
291 [this](Variant const &value) {
292 auto const quat = get_quaternion_safe("quat", value);
293 set_particle_property([&quat](Particle &p) { p.quat() = quat; });
294 },
295 [this]() { return quat2vector(get_particle_data(m_pid).quat()); }},
296 {"omega_body",
297 [this](Variant const &value) {
298 set_particle_property(&Particle::omega, value);
299 },
300 [this]() { return get_particle_data(m_pid).omega(); }},
301 {"rotation",
302 [this](Variant const &value) {
303 set_particle_property([&value](Particle &p) {
304 auto const rotation_flag =
305 Utils::Vector3i{(get_value<Utils::Vector3b>(value))};
306 p.rotation() = bitfield_from_flag(rotation_flag);
307 });
308 },
309 [this]() {
310 auto const rotation_bits = get_particle_data(m_pid).rotation();
311 return Utils::Vector3b{{::detail::get_nth_bit(rotation_bits, 0),
312 ::detail::get_nth_bit(rotation_bits, 1),
313 ::detail::get_nth_bit(rotation_bits, 2)}};
314 }},
315 {"omega_lab",
316 [this](Variant const &value) {
317 set_particle_property([&value](Particle &p) {
318 auto const omega = get_value<Utils::Vector3d>(value);
319 p.omega() = convert_vector_space_to_body(p, omega);
320 });
321 },
322 [this]() {
323 auto &p = get_particle_data(m_pid);
324 return convert_vector_body_to_space(p, p.omega());
325 }},
326 {"torque_lab",
327 [this](Variant const &value) {
328 set_particle_property([&value](Particle &p) {
329 auto const torque = get_value<Utils::Vector3d>(value);
330 p.torque() = convert_vector_space_to_body(p, torque);
331 });
332 },
333 [this]() {
334 auto &p = get_particle_data(m_pid);
335 return convert_vector_body_to_space(p, p.torque());
336 }},
337#endif // ROTATION
338#ifdef ROTATIONAL_INERTIA
339 {"rinertia",
340 [this](Variant const &value) {
341 set_particle_property(&Particle::rinertia, value);
342 },
343 [this]() { return get_particle_data(m_pid).rinertia(); }},
344#endif // ROTATIONAL_INERTIA
345#ifdef LB_ELECTROHYDRODYNAMICS
346 {"mu_E",
347 [this](Variant const &value) {
348 set_particle_property(&Particle::mu_E, value);
349 },
350 [this]() { return get_particle_data(m_pid).mu_E(); }},
351#endif // LB_ELECTROHYDRODYNAMICS
352#ifdef EXTERNAL_FORCES
353 {"fix",
354 [this](Variant const &value) {
355 set_particle_property([&value](Particle &p) {
356 auto const fix_flag =
357 Utils::Vector3i{(get_value<Utils::Vector3b>(value))};
358 p.fixed() = bitfield_from_flag(fix_flag);
359 });
360 },
361 [this]() {
362 auto const fixed = get_particle_data(m_pid).fixed();
363 return Utils::Vector3b{{::detail::get_nth_bit(fixed, 0),
364 ::detail::get_nth_bit(fixed, 1),
365 ::detail::get_nth_bit(fixed, 2)}};
366 }},
367 {"ext_force",
368 [this](Variant const &value) {
369 set_particle_property(&Particle::ext_force, value);
370 },
371 [this]() { return get_particle_data(m_pid).ext_force(); }},
372#ifdef ROTATION
373 {"ext_torque",
374 [this](Variant const &value) {
375 set_particle_property(&Particle::ext_torque, value);
376 },
377 [this]() { return get_particle_data(m_pid).ext_torque(); }},
378#endif // ROTATION
379#endif // EXTERNAL_FORCES
380#ifdef THERMOSTAT_PER_PARTICLE
381 {"gamma",
382 [this](Variant const &value) {
383 set_particle_property(&Particle::gamma,
384 Variant{get_gamma_safe(value)});
385 },
386 [this]() { return get_particle_data(m_pid).gamma(); }},
387#ifdef ROTATION
388 {"gamma_rot",
389 [this](Variant const &value) {
390 set_particle_property(&Particle::gamma_rot,
391 Variant{get_gamma_safe(value)});
392 },
393 [this]() { return get_particle_data(m_pid).gamma_rot(); }},
394#endif // ROTATION
395#endif // THERMOSTAT_PER_PARTICLE
396 {"pos_folded", AutoParameter::read_only,
397 [this]() {
398 auto const &box_geo = *get_system()->box_geo;
399 return box_geo.folded_position(get_particle_data(m_pid).pos());
400 }},
401 {"lees_edwards_offset",
402 [this](Variant const &value) {
403 set_particle_property(&Particle::lees_edwards_offset, value);
404 },
405 [this]() { return get_particle_data(m_pid).lees_edwards_offset(); }},
406 {"lees_edwards_flag", AutoParameter::read_only,
407 [this]() { return get_particle_data(m_pid).lees_edwards_flag(); }},
408 {"image_box", AutoParameter::read_only,
409 [this]() {
410 auto const &box_geo = *get_system()->box_geo;
411 auto const p = get_particle_data(m_pid);
412 return box_geo.folded_image_box(p.pos(), p.image_box());
413 }},
414 {"node", AutoParameter::read_only,
415 [this]() {
416 return (context()->is_head_node()) ? get_particle_node(m_pid) : -1;
417 }},
418 {"mol_id",
419 [this](Variant const &value) {
420 auto const mol_id = get_value<int>(value);
421 if (mol_id < 0) {
422 throw std::domain_error(
423 error_msg("mol_id", "must be an integer >= 0"));
424 }
425 set_particle_property(&Particle::mol_id, Variant{mol_id});
426 },
427 [this]() { return get_particle_data(m_pid).mol_id(); }},
428#ifdef VIRTUAL_SITES_RELATIVE
429 {"vs_quat",
430 [this](Variant const &value) {
431 auto const quat = get_quaternion_safe("vs_quat", value);
432 set_particle_property(
433 [&quat](Particle &p) { p.vs_relative().quat = quat; });
434 },
435 [this]() {
436 return quat2vector(get_particle_data(m_pid).vs_relative().quat);
437 }},
438 {"vs_relative",
439 [this](Variant const &value) {
441 try {
442 auto const array = get_value<std::vector<Variant>>(value);
443 if (array.size() != 3) {
444 throw 0;
445 }
446 vs_relative.distance = get_value<double>(array[1]);
447 vs_relative.to_particle_id = get_value<int>(array[0]);
448 vs_relative.rel_orientation =
449 get_quaternion_safe("vs_relative", array[2]);
450 } catch (...) {
451 throw std::invalid_argument(error_msg(
452 "vs_relative", "must take the form [id, distance, quaternion]"));
453 }
454 set_particle_property(
455 [&vs_relative](Particle &p) { p.vs_relative() = vs_relative; });
456 },
457 [this]() {
458 auto const vs_rel = get_particle_data(m_pid).vs_relative();
459 return std::vector<Variant>{{vs_rel.to_particle_id, vs_rel.distance,
460 quat2vector(vs_rel.rel_orientation)}};
461 }},
462#endif // VIRTUAL_SITES_RELATIVE
463 {"propagation",
464 [this](Variant const &value) {
465 auto const propagation = get_value<int>(value);
466 if (!is_valid_propagation_combination(propagation)) {
467 throw std::domain_error(error_msg(
468 "propagation", "propagation combination not accepted: " +
469 propagation_bitmask_to_string(propagation)));
470 }
471 set_particle_property(&Particle::propagation, value);
472 },
473 [this]() { return get_particle_data(m_pid).propagation(); }},
474#ifdef ENGINE
475 {"swimming",
476 [this](Variant const &value) {
477 set_particle_property([&value](Particle &p) {
479 swim.swimming = true;
480 auto const dict = get_value<VariantMap>(value);
481 if (dict.count("f_swim") != 0) {
482 swim.f_swim = get_value<double>(dict.at("f_swim"));
483 }
484 if (dict.count("is_engine_force_on_fluid") != 0) {
485 auto const is_engine_force_on_fluid =
486 get_value<bool>(dict.at("is_engine_force_on_fluid"));
487 swim.is_engine_force_on_fluid = is_engine_force_on_fluid;
488 }
489 p.swimming() = swim;
490 });
491 },
492 [this]() {
493 auto const swim = get_particle_data(m_pid).swimming();
494 return VariantMap{
495 {"f_swim", swim.f_swim},
496 {"is_engine_force_on_fluid", swim.is_engine_force_on_fluid},
497 };
498 }},
499#endif // ENGINE
500 });
501}
502
503#ifdef EXCLUSIONS
504/**
505 * @brief Locally add an exclusion to a particle.
506 * @param pid1 the identity of the first exclusion partner
507 * @param pid2 the identity of the second exclusion partner
508 * @param cell_structure the cell structure
509 */
510static void local_add_exclusion(int pid1, int pid2,
511 ::CellStructure &cell_structure) {
512 if (auto p1 = cell_structure.get_local_particle(pid1)) {
513 add_exclusion(*p1, pid2);
514 }
515 if (auto p2 = cell_structure.get_local_particle(pid2)) {
516 add_exclusion(*p2, pid1);
517 }
518}
519
520/**
521 * @brief Locally remove an exclusion to a particle.
522 * @param pid1 the identity of the first exclusion partner
523 * @param pid2 the identity of the second exclusion partner
524 * @param cell_structure the cell structure
525 */
526static void local_remove_exclusion(int pid1, int pid2,
527 ::CellStructure &cell_structure) {
528 if (auto p1 = cell_structure.get_local_particle(pid1)) {
529 delete_exclusion(*p1, pid2);
530 }
531 if (auto p2 = cell_structure.get_local_particle(pid2)) {
532 delete_exclusion(*p2, pid1);
533 }
534}
535
536void ParticleHandle::particle_exclusion_sanity_checks(int pid1,
537 int pid2) const {
538 if (pid1 == pid2) {
539 throw std::runtime_error("Particles cannot exclude themselves (id " +
540 std::to_string(pid1) + ")");
541 }
542 auto cell_structure_si = get_cell_structure();
543 auto &cell_structure = cell_structure_si->get_cell_structure();
544 std::ignore = get_real_particle(context()->get_comm(), pid1, cell_structure);
545 std::ignore = get_real_particle(context()->get_comm(), pid2, cell_structure);
546}
547#endif // EXCLUSIONS
548
549Variant ParticleHandle::do_call_method(std::string const &name,
550 VariantMap const &params) {
551 if (name == "set_param_parallel") {
552 auto const param_name = get_value<std::string>(params, "name");
553 if (params.count("value") == 0) {
554 throw Exception("Parameter '" + param_name + "' is missing.");
555 }
556 auto const &value = params.at("value");
557 context()->parallel_try_catch(
558 [&]() { do_set_parameter(param_name, value); });
559 return {};
560 }
561 if (name == "get_bond_by_id") {
562 if (not context()->is_head_node()) {
563 return {};
564 }
565 return get_bonded_ias()->call_method("get_bond", params);
566 }
567 if (name == "get_bonds_view") {
568 if (not context()->is_head_node()) {
569 return {};
570 }
571 auto const bond_list = get_particle_data(m_pid).bonds();
572 std::vector<std::vector<Variant>> bonds_flat;
573 for (auto const &&bond_view : bond_list) {
574 std::vector<Variant> bond_flat;
575 bond_flat.emplace_back(bond_view.bond_id());
576 for (auto const pid : bond_view.partner_ids()) {
577 bond_flat.emplace_back(pid);
578 }
579 bonds_flat.emplace_back(std::move(bond_flat));
580 }
581 return make_vector_of_variants(bonds_flat);
582 }
583 if (name == "add_bond") {
584 auto const bond_id = get_value<int>(params, "bond_id");
585 auto const partner_ids = get_value<std::vector<int>>(params, "part_id");
586 std::vector<int> particle_ids = {m_pid};
587 std::ranges::copy(partner_ids, std::back_inserter(particle_ids));
588 ::add_bond(*get_system(), bond_id, particle_ids);
589 get_system()->on_particle_change();
590 } else if (name == "del_bond") {
591 set_particle_property([&params](Particle &p) {
592 auto const bond_id = get_value<int>(params, "bond_id");
593 auto const part_id = get_value<std::vector<int>>(params, "part_id");
594 auto const bond_view =
595 BondView(bond_id, {part_id.data(), part_id.size()});
596 auto &bond_list = p.bonds();
597 auto it = std::find(bond_list.begin(), bond_list.end(), bond_view);
598 if (it != bond_list.end()) {
599 bond_list.erase(it);
600 }
601 });
602 } else if (name == "delete_all_bonds") {
603 set_particle_property([&](Particle &p) { p.bonds().clear(); });
604 } else if (name == "is_valid_bond_id") {
605 auto const bond_id = get_value<int>(params, "bond_id");
606 return get_system()->bonded_ias->get_zero_based_type(bond_id) != 0;
607 }
608 if (name == "remove_particle") {
609 context()->parallel_try_catch([&]() {
610 auto cell_structure_si = get_cell_structure();
611 auto &cell_structure = cell_structure_si->get_cell_structure();
612 std::ignore =
613 get_real_particle(context()->get_comm(), m_pid, cell_structure);
614 remove_particle(m_pid);
615 });
616 } else if (name == "is_virtual") {
617 if (not context()->is_head_node()) {
618 return {};
619 }
620 return get_particle_data(m_pid).is_virtual();
621#ifdef VIRTUAL_SITES_RELATIVE
622 } else if (name == "vs_auto_relate_to") {
623 if (not context()->is_head_node()) {
624 return {};
625 }
626 auto const other_pid = get_value<int>(params, "pid");
627 auto const override_cutoff_check =
628 get_value<bool>(params, "override_cutoff_check");
629 if (m_pid == other_pid) {
630 throw std::invalid_argument("A virtual site cannot relate to itself");
631 }
632 if (other_pid < 0) {
633 throw std::domain_error("Invalid particle id: " +
634 std::to_string(other_pid));
635 }
636 auto const system = get_system();
637 /* note: this code can be rewritten as parallel code, but only with a call
638 * to `cells_update_ghosts(DATA_PART_POSITION | DATA_PART_PROPERTIES)`, as
639 * there is no guarantee the virtual site has visibility of the relative
640 * particle through the ghost layer during particle creation. However,
641 * ghost updates can scramble the particle ordering in the local cells,
642 * which is an issue for checkpointing: the H5MD writer will use the
643 * scrambled ordering before writing to a checkpoint file and the
644 * non-scrambled ordering after reloading from a checkpoint file.
645 */
646 auto const &p_current = get_particle_data(m_pid);
647 auto const &p_relate_to = get_particle_data(other_pid);
648 auto const [quat, dist] = calculate_vs_relate_to_params(
649 p_current, p_relate_to, *system->box_geo, system->get_min_global_cut(),
650 override_cutoff_check);
651 set_parameter("vs_relative", Variant{std::vector<Variant>{
652 {other_pid, dist, quat2vector(quat)}}});
653 set_parameter("propagation",
656#endif // VIRTUAL_SITES_RELATIVE
657#ifdef EXCLUSIONS
658 } else if (name == "has_exclusion") {
659 auto const other_pid = get_value<int>(params, "pid");
660 auto cell_structure_si = get_cell_structure();
661 auto &cell_structure = cell_structure_si->get_cell_structure();
662 auto const p =
663 get_real_particle(context()->get_comm(), m_pid, cell_structure);
664 if (p != nullptr) {
665 return p->has_exclusion(other_pid);
666 }
667 }
668 if (name == "add_exclusion") {
669 auto const other_pid = get_value<int>(params, "pid");
670 auto cell_structure_si = get_cell_structure();
671 auto &cell_structure = cell_structure_si->get_cell_structure();
672 context()->parallel_try_catch(
673 [&]() { particle_exclusion_sanity_checks(m_pid, other_pid); });
674 local_add_exclusion(m_pid, other_pid, cell_structure);
675 get_system()->on_particle_change();
676 } else if (name == "del_exclusion") {
677 auto const other_pid = get_value<int>(params, "pid");
678 auto cell_structure_si = get_cell_structure();
679 auto &cell_structure = cell_structure_si->get_cell_structure();
680 context()->parallel_try_catch(
681 [&]() { particle_exclusion_sanity_checks(m_pid, other_pid); });
682 local_remove_exclusion(m_pid, other_pid, cell_structure);
683 get_system()->on_particle_change();
684 } else if (name == "set_exclusions") {
685 std::vector<int> exclusion_list;
686 try {
687 auto const pid = get_value<int>(params, "p_ids");
688 exclusion_list.push_back(pid);
689 } catch (...) {
690 exclusion_list = get_value<std::vector<int>>(params, "p_ids");
691 }
692 context()->parallel_try_catch([&]() {
693 for (auto const pid : exclusion_list) {
694 particle_exclusion_sanity_checks(m_pid, pid);
695 }
696 });
697 set_particle_property([this, &exclusion_list](Particle &p) {
698 auto cell_structure_si = get_cell_structure();
699 auto &cell_structure = cell_structure_si->get_cell_structure();
700 for (auto const pid : p.exclusions()) {
701 local_remove_exclusion(m_pid, pid, cell_structure);
702 }
703 for (auto const pid : exclusion_list) {
704 if (!p.has_exclusion(pid)) {
705 local_add_exclusion(m_pid, pid, cell_structure);
706 }
707 }
708 });
709 } else if (name == "get_exclusions") {
710 if (not context()->is_head_node()) {
711 return {};
712 }
713 auto const excl_list = get_particle_data(m_pid).exclusions();
714 return Variant{std::vector<int>{excl_list.begin(), excl_list.end()}};
715#endif // EXCLUSIONS
716#ifdef ROTATION
717 }
718 if (name == "rotate_particle") {
719 set_particle_property([&params](Particle &p) {
720 auto const axis = get_value<Utils::Vector3d>(params, "axis");
721 auto const angle = get_value<double>(params, "angle");
722 local_rotate_particle(p, axis, angle);
723 });
724 }
725 if (name == "convert_vector_body_to_space") {
726 return get_particle_property<std::vector<double>>(
727 [&params](Particle const &p) {
728 auto const vec = get_value<Utils::Vector3d>(params, "vec");
730 });
731 }
732 if (name == "convert_vector_space_to_body") {
733 return get_particle_property<std::vector<double>>(
734 [&params](Particle const &p) {
735 auto const vec = get_value<Utils::Vector3d>(params, "vec");
737 });
738#endif // ROTATION
739 }
740 return {};
741}
742
743#ifdef ROTATION
744static auto const contradicting_arguments_quat = std::vector<
745 std::array<std::string, 3>>{{
746 {{"dip", "dipm",
747 "Setting 'dip' is sufficient as it defines the scalar dipole moment."}},
748 {{"quat", "director",
749 "Setting 'quat' is sufficient as it defines the director."}},
750 {{"dip", "quat",
751 "Setting 'dip' would overwrite 'quat'. Set 'quat' and 'dipm' instead."}},
752 {{"dip", "director",
753 "Setting 'dip' would overwrite 'director'. Set 'director' and "
754 "'dipm' instead."}},
755}};
756#endif // ROTATION
757
758void ParticleHandle::do_construct(VariantMap const &params) {
759 auto const n_extra_args = params.size() - params.count("id") -
760 params.count("__cell_structure") -
761 params.count("__bonded_ias");
762 m_pid = (params.contains("id")) ? get_value<int>(params, "id")
764
765#ifndef NDEBUG
766 if (not params.contains("id")) {
767 auto head_node_reference = m_pid;
768 boost::mpi::broadcast(context()->get_comm(), head_node_reference, 0);
769 assert(m_pid == head_node_reference && "global max_seen_pid has diverged");
770 }
771#endif
772
773 if (params.contains("__cell_structure")) {
774 auto so = get_value<std::shared_ptr<CellSystem::CellSystem>>(
775 params, "__cell_structure");
776 so->configure(*this);
777 m_cell_structure = so;
778 }
779 if (params.contains("__bonded_ias")) {
780 m_bonded_ias = get_value<std::shared_ptr<Interactions::BondedInteractions>>(
781 params, "__bonded_ias");
782 }
783
784 // create a new particle if extra arguments were passed
785 if (n_extra_args == 0) {
786 return;
787 }
788
789 auto const pos = get_value<Utils::Vector3d>(params, "pos");
790 context()->parallel_try_catch([&]() {
791 particle_checks(m_pid, pos);
792 auto cell_structure_si = get_cell_structure();
793 auto &cell_structure = cell_structure_si->get_cell_structure();
794 auto ptr = cell_structure.get_local_particle(m_pid);
795 if (ptr != nullptr) {
796 throw std::invalid_argument("Particle " + std::to_string(m_pid) +
797 " already exists");
798 }
799 });
800
801#ifdef ROTATION
802 context()->parallel_try_catch([&]() {
803 // if we are not constructing a particle from a checkpoint file,
804 // check the quaternion is not accidentally set twice by the user
805 if (not params.contains("__cpt_sentinel")) {
806 auto formatter =
807 boost::format("Contradicting particle attributes: '%s' and '%s'. %s");
808 for (auto const &[prop1, prop2, reason] : contradicting_arguments_quat) {
809 if (params.contains(prop1) and params.contains(prop2)) {
810 auto const err_msg = boost::str(formatter % prop1 % prop2 % reason);
811 throw std::invalid_argument(err_msg);
812 }
813 }
814 }
815 });
816#endif // ROTATION
817
818 // create a default-constructed particle
819 make_new_particle(m_pid, pos);
820
821 try {
822 context()->parallel_try_catch([&]() {
823 /* clang-format off */
824 // set particle properties (filter out read-only and deferred properties)
825 std::set<std::string> const skip = {
826 "pos_folded", "pos", "id", "exclusions", "node", "image_box", "bonds",
827 "lees_edwards_flag", "__cpt_sentinel",
828 };
829 /* clang-format on */
830 for (auto const &name : get_parameter_insertion_order()) {
831 if (params.contains(name) and not skip.contains(name)) {
832 do_set_parameter(name, params.at(name));
833 }
834 }
835 if (not params.contains("type")) {
836 do_set_parameter("type", 0);
837 }
838 });
839 } catch (...) {
840 remove_particle(m_pid);
841 throw;
842 }
843}
844
845} // namespace Particles
846} // namespace ScriptInterface
Vector implementation and trait types for boost qvm interoperability.
Data structures for bonded interactions.
bool add_bond(System::System &system, int bond_id, std::vector< int > const &particle_ids)
Add a bond to a particle.
Definition bonds.cpp:25
Immutable view on a bond.
Definition BondList.hpp:44
virtual boost::mpi::communicator const & get_comm() const =0
Context * context() const
Responsible context.
std::shared_ptr< BondedInteractionsMap > bonded_ias
void on_particle_change()
Called every time a particle property changes.
std::shared_ptr< BoxGeometry > box_geo
std::shared_ptr< InteractionsNonBonded > nonbonded_ias
std::vector< T > as_vector() const
Definition Vector.hpp:119
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
This file contains the defaults for ESPResSo.
void delete_exclusion(Particle &p, int p_id)
Remove exclusion from particle if possible.
void add_exclusion(Particle &p, int p_id)
Insert an exclusion if not already set.
static uint8_t bitfield_from_flag(Utils::Vector3i const &flag)
static void local_remove_exclusion(int pid1, int pid2, ::CellStructure &cell_structure)
Locally remove an exclusion to a particle.
static void particle_checks(int p_id, Utils::Vector3d const &pos)
static void local_add_exclusion(int pid1, int pid2, ::CellStructure &cell_structure)
Locally add an exclusion to a particle.
static auto get_quaternion_safe(std::string const &name, Variant const &value)
static auto const contradicting_arguments_quat
static auto get_gamma_safe(Variant const &value)
static auto quat2vector(Utils::Quaternion< double > const &q)
static auto get_real_particle(boost::mpi::communicator const &comm, int p_id, ::CellStructure &cell_structure)
static auto error_msg(std::string const &name, std::string const &reason)
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_vector_of_variants(std::vector< T > const &v)
Definition Variant.hpp:88
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
System & get_system()
T reduce_optional(boost::mpi::communicator const &comm, std::optional< T > const &result)
Reduce an optional on the head node.
Quaternion< T > convert_director_to_quaternion(Vector< T, 3 > const &d)
Convert director to quaternion.
Vector< T, 3 > convert_quaternion_to_director(Quaternion< T > const &quat)
Convert quaternion to director.
Various procedures concerning interactions between particles.
void make_new_particle(int p_id, Utils::Vector3d const &pos)
Create a new particle and attach it to a cell.
const Particle & get_particle_data(int p_id)
Get particle data.
int get_particle_node(int p_id)
Get the MPI rank which owns the a specific particle.
void set_particle_pos(int p_id, Utils::Vector3d const &pos)
Move particle to a new position.
void remove_particle(int p_id)
Remove particle with a given identity.
int get_maximal_particle_id()
Get maximal particle id.
static auto & get_cell_structure()
void on_particle_type_change(int p_id, int old_type, int new_type)
Particles creation and deletion.
std::string propagation_bitmask_to_string(int propagation)
Convert a propagation modes bitmask to a string.
bool is_valid_propagation_combination(int propagation)
Note for developers: when enabling new propagation mode combinations, make sure every single line of ...
This file contains all subroutines required to process rotational motion.
Utils::Vector3d convert_vector_body_to_space(const Particle &p, const Utils::Vector3d &vec)
Definition rotation.hpp:57
std::pair< Utils::Quaternion< double >, double > convert_dip_to_quat(const Utils::Vector3d &dip)
convert a dipole moment to quaternions and dipolar strength
Definition rotation.hpp:108
Utils::Vector3d convert_vector_space_to_body(const Particle &p, const Utils::Vector3d &v)
Definition rotation.hpp:61
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.
Particle * get_local_particle(int id)
Get a local particle by id.
Properties of a self-propelled particle.
Definition Particle.hpp:50
bool swimming
Is the particle a swimmer.
Definition Particle.hpp:54
The following properties define, with respect to which real particle a virtual site is placed and at ...
Definition Particle.hpp:152
Struct holding all information for one particle.
Definition Particle.hpp:395
auto const & dip_fld() const
Definition Particle.hpp:498
bool has_exclusion(int pid) const
Definition Particle.hpp:576
auto const & swimming() const
Definition Particle.hpp:561
auto const & lees_edwards_offset() const
Definition Particle.hpp:446
auto const & propagation() const
Definition Particle.hpp:421
auto const & rinertia() const
Definition Particle.hpp:502
auto is_virtual() const
Definition Particle.hpp:518
auto const & mass() const
Definition Particle.hpp:452
Utils::compact_vector< int > & exclusions()
Definition Particle.hpp:574
auto const & quat() const
Definition Particle.hpp:477
auto const & rotation() const
Definition Particle.hpp:458
auto const & vs_relative() const
Definition Particle.hpp:527
auto const & q() const
Definition Particle.hpp:508
auto const & gamma() const
Definition Particle.hpp:531
auto const & gamma_rot() const
Definition Particle.hpp:534
auto const & lees_edwards_flag() const
Definition Particle.hpp:448
auto calc_dip() const
Definition Particle.hpp:495
auto const & v() const
Definition Particle.hpp:433
auto const & torque() const
Definition Particle.hpp:479
auto const & fixed() const
Definition Particle.hpp:539
auto const & ext_force() const
Definition Particle.hpp:554
auto const & omega() const
Definition Particle.hpp:481
auto const & type() const
Definition Particle.hpp:418
auto const & ext_torque() const
Definition Particle.hpp:484
auto const & bonds() const
Definition Particle.hpp:428
auto const & dipm() const
Definition Particle.hpp:493
auto const & mol_id() const
Definition Particle.hpp:416
auto const & mu_E() const
Definition Particle.hpp:514
auto const & force() const
Definition Particle.hpp:435
Quaternion representation.
std::tuple< Utils::Quaternion< double >, double > calculate_vs_relate_to_params(Particle const &p_vs, Particle const &p_relate_to, BoxGeometry const &box_geo, double min_global_cut, bool override_cutoff_check)
Calculate the rotation quaternion and distance between two particles.