ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
ShapeBasedConstraint.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2010-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
21
22#include "BoxGeometry.hpp"
23#include "Observable_stat.hpp"
25#include "communication.hpp"
26#include "config/config.hpp"
27#include "dpd.hpp"
28#include "energy_inline.hpp"
29#include "errorhandling.hpp"
30#include "forces_inline.hpp"
32#include "system/System.hpp"
33#include "thermostat.hpp"
34
35#include <utils/Vector.hpp>
36
37#include <boost/mpi/collectives.hpp>
38
39#include <algorithm>
40#include <functional>
41#include <limits>
42#include <numeric>
43
44namespace Constraints {
45/** Check if a non-bonded interaction is defined */
46static bool is_active(IA_parameters const &data) {
47 return data.max_cut != inactive_cutoff;
48}
49
51 part_rep.type() = type;
52 m_system.lock()->nonbonded_ias->make_particle_type_exist(type);
53}
54
55IA_parameters const &ShapeBasedConstraint::get_ia_param(int type) const {
56 return m_system.lock()->nonbonded_ias->get_ia_param(type, part_rep.type());
57}
58
60 return all_reduce(comm_cart, m_local_force, std::plus<>());
61}
62
64 return all_reduce(comm_cart, m_outer_normal_force, std::plus<double>());
65}
66
68 ParticleRange const &particles) const {
69 auto global_mindist = std::numeric_limits<double>::infinity();
70
71 auto const local_mindist = std::accumulate(
72 particles.begin(), particles.end(),
73 std::numeric_limits<double>::infinity(),
74 [this, &box_geo](double acc, Particle const &p) {
75 auto const &ia_params = get_ia_param(p.type());
76 if (is_active(ia_params)) {
77 double dist;
78 Utils::Vector3d vec;
79 m_shape->calculate_dist(box_geo.folded_position(p.pos()), dist, vec);
80 acc = std::min(acc, dist);
81 }
82 return acc;
83 });
84 boost::mpi::reduce(comm_cart, local_mindist, global_mindist,
85 boost::mpi::minimum<double>(), 0);
86 return global_mindist;
87}
88
89ParticleForce ShapeBasedConstraint::force(Particle const &p,
90 Utils::Vector3d const &folded_pos,
91 double) {
92 ParticleForce pf{};
93 auto const &ia_params = get_ia_param(p.type());
94
95 if (is_active(ia_params)) {
96 double dist = 0.;
97 Utils::Vector3d dist_vec;
98 m_shape->calculate_dist(folded_pos, dist, dist_vec);
99 auto &system = *m_system.lock();
100 auto const coulomb_kernel = system.coulomb.pair_force_kernel();
101
102#ifdef ESPRESSO_DPD
103 Utils::Vector3d dpd_force{};
104#endif
105 Utils::Vector3d outer_normal_vec{};
106
107 if (dist > 0) {
108 outer_normal_vec = -dist_vec / dist;
109 pf.f = calc_central_radial_force(ia_params, dist_vec, dist);
110#ifdef ESPRESSO_THOLE
111 pf.f += thole_pair_force(p, part_rep, ia_params, dist_vec, dist,
112 *system.bonded_ias, get_ptr(coulomb_kernel));
113#endif
114 pf += calc_non_central_force(p, part_rep, ia_params, dist_vec, dist);
115
116#ifdef ESPRESSO_DPD
117 if (system.thermostat->thermo_switch & THERMO_DPD) {
118 dpd_force = dpd_pair_force(p.pos(), p.v(), p.id(), part_rep.pos(),
119 part_rep.v(), part_rep.id(),
120 *system.thermostat->dpd, *system.box_geo,
121 ia_params, dist_vec, dist, dist * dist);
122 // Additional use of DPD here requires counter increase
123 system.thermostat->dpd->rng_increment();
124 }
125#endif
126 } else if (m_penetrable && (dist <= 0)) {
127 if ((!m_only_positive) && (dist < 0)) {
128 pf.f = calc_central_radial_force(ia_params, dist_vec, -dist);
129#ifdef ESPRESSO_THOLE
130 pf.f += thole_pair_force(p, part_rep, ia_params, dist_vec, -dist,
131 *system.bonded_ias, get_ptr(coulomb_kernel));
132#endif
133 pf += calc_non_central_force(p, part_rep, ia_params, dist_vec, -dist);
134
135#ifdef ESPRESSO_DPD
136 if (system.thermostat->thermo_switch & THERMO_DPD) {
137 dpd_force = dpd_pair_force(p.pos(), p.v(), p.id(), part_rep.pos(),
138 part_rep.v(), part_rep.id(),
139 *system.thermostat->dpd, *system.box_geo,
140 ia_params, dist_vec, dist, dist * dist);
141 // Additional use of DPD here requires counter increase
142 system.thermostat->dpd->rng_increment();
143 }
144#endif
145 }
146 } else {
147 runtimeErrorMsg() << "Constraint violated by particle " << p.id()
148 << " dist " << dist;
149 }
150
151#ifdef ESPRESSO_ROTATION
152 part_rep.torque() += calc_opposing_force(pf, dist_vec).torque;
153#endif
154#ifdef ESPRESSO_DPD
155 pf.f += dpd_force;
156#endif
157 m_local_force -= pf.f;
158 m_outer_normal_force -= outer_normal_vec * pf.f;
159 }
160 return pf;
161}
162
163void ShapeBasedConstraint::add_energy(const Particle &p,
164 const Utils::Vector3d &folded_pos, double,
165 Observable_stat &obs_energy) const {
166 double energy = 0.0;
167 auto const &ia_params = get_ia_param(p.type());
168
169 if (is_active(ia_params)) {
170 auto &system = *m_system.lock();
171 auto const coulomb_kernel = system.coulomb.pair_energy_kernel();
172 double dist = 0.0;
173 Utils::Vector3d vec;
174 m_shape->calculate_dist(folded_pos, dist, vec);
175 if (dist > 0.) {
176 energy = calc_non_bonded_pair_energy(p, part_rep, ia_params, vec, dist,
177 *system.bonded_ias,
178 get_ptr(coulomb_kernel));
179 } else if (dist <= 0. and m_penetrable) {
180 if (!m_only_positive and dist < 0.) {
181 energy = calc_non_bonded_pair_energy(p, part_rep, ia_params, vec, -dist,
182 *system.bonded_ias,
183 get_ptr(coulomb_kernel));
184 }
185 } else {
186 runtimeErrorMsg() << "Constraint violated by particle " << p.id();
187 }
188 }
189 // NOLINTNEXTLINE(clang-analyzer-cplusplus.NewDeleteLeaks)
190 if (part_rep.type() >= 0) {
192 p.type(), part_rep.type(), p.mol_id(), part_rep.mol_id(), energy);
193 }
194}
195} // namespace Constraints
@ THERMO_DPD
Vector implementation and trait types for boost qvm interoperability.
Data structures for bonded interactions.
double min_dist(BoxGeometry const &box_geo, ParticleRange const &particles) const
Calculate the minimum distance between all particle pairs.
Observable for the pressure and energy.
void add_non_bonded_contribution(int type1, int type2, int molid1, int molid2, std::span< const double > data)
A range of particles.
boost::mpi::communicator comm_cart
The communicator.
constexpr double inactive_cutoff
Special cutoff value for an inactive interaction.
Definition config.hpp:44
const T * get_ptr(std::optional< T > const &opt)
Utils::Vector3d dpd_pair_force(DPDParameters const &params, Utils::Vector3d const &v, double dist, Utils::Vector3d const &noise)
Definition dpd.cpp:86
Routines to use DPD as thermostat or pair force .
Energy calculation.
double calc_non_bonded_pair_energy(Particle const &p1, Particle const &p2, IA_parameters const &ia_params, Utils::Vector3d const &d, double const dist, BondedInteractionsMap const &bonded_ias, Coulomb::ShortRangeEnergyKernel::kernel_type const *coulomb_kernel)
Calculate non-bonded energies between a pair of particles.
This file contains the errorhandling code for severe errors, like a broken bond or illegal parameter ...
#define runtimeErrorMsg()
Force calculation.
ParticleForce calc_opposing_force(ParticleForce const &pf, Utils::Vector3d const &d)
ParticleForce calc_non_central_force(Particle const &p1, Particle const &p2, IA_parameters const &ia_params, Utils::Vector3d const &d, double const dist)
ESPRESSO_ATTR_ALWAYS_INLINE Utils::Vector3d calc_central_radial_force(IA_parameters const &ia_params, Utils::Vector3d const &d, double const dist)
static bool is_active(IA_parameters const &data)
Check if a non-bonded interaction is defined.
Various procedures concerning interactions between particles.
Parameters for non-bonded interactions.
double max_cut
maximal cutoff for this pair of particle types.
Force information on a particle.
Definition Particle.hpp:345
Utils::Vector3d torque
torque.
Definition Particle.hpp:375
Struct holding all information for one particle.
Definition Particle.hpp:450
auto const & v() const
Definition Particle.hpp:488
auto const & type() const
Definition Particle.hpp:473
auto const & mol_id() const
Definition Particle.hpp:471
auto const & pos() const
Definition Particle.hpp:486
auto const & id() const
Definition Particle.hpp:469
Utils::Vector3d thole_pair_force(Particle const &p1, Particle const &p2, IA_parameters const &ia_params, Utils::Vector3d const &d, double dist, BondedInteractionsMap const &bonded_ias, Coulomb::ShortRangeForceKernel::kernel_type const *kernel)
Calculate Thole force.
Definition thole.hpp:44