ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
forces_inline.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2010-2026 The ESPResSo project
3 * Copyright (C) 2002,2003,2004,2005,2006,2007,2008,2009,2010
4 * Max-Planck-Institute for Polymer Research, Theory Group
5 *
6 * This file is part of ESPResSo.
7 *
8 * ESPResSo is free software: you can redistribute it and/or modify
9 * it under the terms of the GNU General Public License as published by
10 * the Free Software Foundation, either version 3 of the License, or
11 * (at your option) any later version.
12 *
13 * ESPResSo is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public License
19 * along with this program. If not, see <http://www.gnu.org/licenses/>.
20 */
21
22#pragma once
23
24/** \file
25 * Force calculation.
26 */
27
28#include <config/config.hpp>
29
30#include "BoxGeometry.hpp"
31#include "actor/visitors.hpp"
58
59#ifdef ESPRESSO_DPD
60#include "dpd.hpp"
61#endif
62
63#include "Particle.hpp"
64#include "bond_error.hpp"
65#include "errorhandling.hpp"
66#include "exclusions.hpp"
67#include "thermostat.hpp"
68
69#include <utils/Vector.hpp>
70
71#include <optional>
72#include <span>
73#include <tuple>
74#include <variant>
75
78 Utils::Vector3d const &d,
79 double const dist) {
80
81 auto force_factor = 0.;
82/* Lennard-Jones */
83#ifdef ESPRESSO_LENNARD_JONES
84 force_factor += lj_pair_force_factor(ia_params, dist);
85#endif
86/* WCA */
87#ifdef ESPRESSO_WCA
88 force_factor += wca_pair_force_factor(ia_params, dist);
89#endif
90/* Lennard-Jones generic */
91#ifdef ESPRESSO_LENNARD_JONES_GENERIC
92 force_factor += ljgen_pair_force_factor(ia_params, dist);
93#endif
94/* smooth step */
95#ifdef ESPRESSO_SMOOTH_STEP
96 force_factor += SmSt_pair_force_factor(ia_params, dist);
97#endif
98/* Hertzian force */
99#ifdef ESPRESSO_HERTZIAN
100 force_factor += hertzian_pair_force_factor(ia_params, dist);
101#endif
102/* Gaussian force */
103#ifdef ESPRESSO_GAUSSIAN
104 force_factor += gaussian_pair_force_factor(ia_params, dist);
105#endif
106/* BMHTF NaCl */
107#ifdef ESPRESSO_BMHTF_NACL
108 force_factor += BMHTF_pair_force_factor(ia_params, dist);
109#endif
110/* Buckingham*/
111#ifdef ESPRESSO_BUCKINGHAM
112 force_factor += buck_pair_force_factor(ia_params, dist);
113#endif
114/* Morse*/
115#ifdef ESPRESSO_MORSE
116 force_factor += morse_pair_force_factor(ia_params, dist);
117#endif
118/*soft-sphere potential*/
119#ifdef ESPRESSO_SOFT_SPHERE
120 force_factor += soft_pair_force_factor(ia_params, dist);
121#endif
122/*hat potential*/
123#ifdef ESPRESSO_HAT
124 force_factor += hat_pair_force_factor(ia_params, dist);
125#endif
126/* Lennard-Jones cosine */
127#ifdef ESPRESSO_LJCOS
128 force_factor += ljcos_pair_force_factor(ia_params, dist);
129#endif
130/* Lennard-Jones cosine */
131#ifdef ESPRESSO_LJCOS2
132 force_factor += ljcos2_pair_force_factor(ia_params, dist);
133#endif
134/* tabulated */
135#ifdef ESPRESSO_TABULATED
136 force_factor += tabulated_pair_force_factor(ia_params, dist);
137#endif
138 return force_factor * d;
139}
140
142 Particle const &p2,
143 IA_parameters const &ia_params,
144 Utils::Vector3d const &d,
145 double const dist) {
146
147 ParticleForce pf{};
148/* Gay-Berne */
149#ifdef ESPRESSO_GAY_BERNE
150 pf += gb_pair_force(p1.quat(), p2.quat(), ia_params, d, dist);
151#endif
152 return pf;
153}
154
156 Utils::Vector3d const &d) {
157 ParticleForce out{-pf.f};
158#ifdef ESPRESSO_ROTATION
159 // if torque is a null vector, the opposing torque is a null vector too
160 // (this check guards from returning a small yet non-null opposing
161 // torque due to numerical imprecision)
162 if (pf.torque[0] != 0. || pf.torque[1] != 0. || pf.torque[2] != 0.) {
163 out.torque = -(pf.torque + vector_product(d, pf.f));
164 }
165#endif
166 return out;
167}
168
169/** Compute the bonded interaction force between particle pairs.
170 *
171 * @param[in] iaparams Bonded parameters for the interaction.
172 * @param[in] q1q2 Product of the particle charges.
173 * @param[in] dx Vector between @p p1 and @p p2.
174 * @param[in] kernel Coulomb force kernel.
175 */
177inline std::optional<Utils::Vector3d> calc_bond_pair_force(
178 Bonded_IA_Parameters const &iaparams, Utils::Vector3d const &dx,
179 double const q1q2,
181 if (auto const *iap = std::get_if<FeneBond>(&iaparams)) {
182 return iap->force(dx);
183 }
184 if (auto const *iap = std::get_if<HarmonicBond>(&iaparams)) {
185 return iap->force(dx);
186 }
187 if (auto const *iap = std::get_if<QuarticBond>(&iaparams)) {
188 return iap->force(dx);
189 }
190#ifdef ESPRESSO_ELECTROSTATICS
191 if (auto const *iap = std::get_if<BondedCoulomb>(&iaparams)) {
192 return iap->force(q1q2, dx);
193 }
194 if (auto const *iap = std::get_if<BondedCoulombSR>(&iaparams)) {
195 return iap->force(dx, *kernel);
196 }
197#endif
198#ifdef ESPRESSO_BOND_CONSTRAINT
199 if (std::get_if<RigidBond>(&iaparams)) {
200 return Utils::Vector3d{};
201 }
202#endif
203#ifdef ESPRESSO_TABULATED
204 if (auto const *iap = std::get_if<TabulatedDistanceBond>(&iaparams)) {
205 return iap->force(dx);
206 }
207#endif
208 if (std::get_if<VirtualBond>(&iaparams)) {
209 return Utils::Vector3d{};
210 }
211 throw BondUnknownTypeError();
212}
213
215inline std::optional<
216 std::tuple<Utils::Vector3d, Utils::Vector3d, Utils::Vector3d>>
218 Utils::Vector3d const &vec1,
219 Utils::Vector3d const &vec2) {
220 if (auto const *iap = std::get_if<AngleHarmonicBond>(&iaparams)) {
221 return iap->forces(vec1, vec2);
222 }
223 if (auto const *iap = std::get_if<AngleCosineBond>(&iaparams)) {
224 return iap->forces(vec1, vec2);
225 }
226 if (auto const *iap = std::get_if<AngleCossquareBond>(&iaparams)) {
227 return iap->forces(vec1, vec2);
228 }
229#ifdef ESPRESSO_TABULATED
230 if (auto const *iap = std::get_if<TabulatedAngleBond>(&iaparams)) {
231 return iap->forces(vec1, vec2);
232 }
233#endif
234 if (auto const *iap = std::get_if<IBMTriel>(&iaparams)) {
235 return iap->calc_forces(vec1, vec2);
236 }
237 throw BondUnknownTypeError();
238}
239
241inline std::optional<std::tuple<Utils::Vector3d, Utils::Vector3d,
244 Bonded_IA_Parameters const &iaparams, BoxGeometry const &box_geo,
245 Utils::Vector3d const &pos1, Utils::Vector3d const &pos2,
246 Utils::Vector3d const &pos3, Utils::Vector3d const &pos4,
247 Utils::Vector3d const &vel1, Utils::Vector3d const &vel3,
248 Utils::Vector3i const &image1) {
249 if (auto const *iap = std::get_if<OifLocalForcesBond>(&iaparams)) {
250 // note: particles in a dihedral bond are ordered as p2-p1-p3-p4
251 auto const fp2 = box_geo.unfolded_position(pos1, image1);
252 auto const fp1 = fp2 + box_geo.get_mi_vector(pos2, fp2);
253 auto const fp3 = fp2 + box_geo.get_mi_vector(pos3, fp2);
254 auto const fp4 = fp2 + box_geo.get_mi_vector(pos4, fp2);
255 return iap->calc_forces(fp2, fp1, fp3, fp4, vel1, vel3);
256 }
257 if (auto const *iap = std::get_if<IBMTribend>(&iaparams)) {
258 return iap->calc_forces(box_geo, pos1, pos2, pos3, pos4);
259 }
260 // note: particles in a dihedral bond are ordered as p2-p1-p3-p4
261 auto const v12 = box_geo.get_mi_vector(pos1, pos2);
262 auto const v23 = box_geo.get_mi_vector(pos3, pos1);
263 auto const v34 = box_geo.get_mi_vector(pos4, pos3);
264 if (auto const *iap = std::get_if<DihedralBond>(&iaparams)) {
265 return iap->forces(v12, v23, v34);
266 }
267#ifdef ESPRESSO_TABULATED
268 if (auto const *iap = std::get_if<TabulatedDihedralBond>(&iaparams)) {
269 return iap->forces(v12, v23, v34);
270 }
271#endif
272 throw BondUnknownTypeError();
273}
Vector implementation and trait types for boost qvm interoperability.
#define ESPRESSO_ATTR_ALWAYS_INLINE
Routines to calculate the Born-Meyer-Huggins-Tosi-Fumi potential between particle pairs.
double BMHTF_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate BMHTF force factor.
Data structures for bonded interactions.
std::variant< NoneBond, FeneBond, HarmonicBond, QuarticBond, BondedCoulomb, BondedCoulombSR, AngleHarmonicBond, AngleCosineBond, AngleCossquareBond, DihedralBond, TabulatedDistanceBond, TabulatedAngleBond, TabulatedDihedralBond, ThermalizedBond, RigidBond, IBMTriel, IBMVolCons, IBMTribend, OifGlobalForcesBond, OifLocalForcesBond, VirtualBond > Bonded_IA_Parameters
Variant in which to store the parameters of an individual bonded interaction.
Routines to calculate the Buckingham potential between particle pairs.
double buck_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate Buckingham force factor.
auto unfolded_position(Utils::Vector3d const &pos, Utils::Vector3i const &image_box) const
Unfold particle coordinates to image box.
ESPRESSO_ATTR_ALWAYS_INLINE Utils::Vector3< T > get_mi_vector(Utils::Vector3< T > const &a, Utils::Vector3< T > const &b) const
Get the minimum-image vector between two coordinates.
Routines to calculate the Gaussian potential between particle pairs.
double gaussian_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate Gaussian force factor.
__device__ void vector_product(float const *a, float const *b, float *out)
Routines to use DPD as thermostat or pair force .
This file contains the errorhandling code for severe errors, like a broken bond or illegal parameter ...
ESPRESSO_ATTR_ALWAYS_INLINE std::optional< Utils::Vector3d > calc_bond_pair_force(Bonded_IA_Parameters const &iaparams, Utils::Vector3d const &dx, double const q1q2, Coulomb::ShortRangeForceKernel::kernel_type const *kernel)
Compute the bonded interaction force between particle pairs.
ESPRESSO_ATTR_ALWAYS_INLINE std::optional< std::tuple< Utils::Vector3d, Utils::Vector3d, Utils::Vector3d, Utils::Vector3d > > calc_bonded_four_body_force(Bonded_IA_Parameters const &iaparams, BoxGeometry const &box_geo, Utils::Vector3d const &pos1, Utils::Vector3d const &pos2, Utils::Vector3d const &pos3, Utils::Vector3d const &pos4, Utils::Vector3d const &vel1, Utils::Vector3d const &vel3, Utils::Vector3i const &image1)
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)
ESPRESSO_ATTR_ALWAYS_INLINE std::optional< std::tuple< Utils::Vector3d, Utils::Vector3d, Utils::Vector3d > > calc_bonded_three_body_force(Bonded_IA_Parameters const &iaparams, Utils::Vector3d const &vec1, Utils::Vector3d const &vec2)
Routines to calculate the Gay-Berne potential between particle pairs.
ParticleForce gb_pair_force(Utils::Vector3d const &ui, Utils::Vector3d const &uj, IA_parameters const &ia_params, Utils::Vector3d const &d, double dist)
Calculate Gay-Berne force and torques.
Definition gay_berne.hpp:49
Routines to calculate the hat potential between particle pairs.
double hat_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate hat force factor.
Definition hat.hpp:37
Routines to calculate the Hertzian potential between particle pairs.
double hertzian_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate Hertzian force factor.
Definition hertzian.hpp:39
Routines to calculate the Lennard-Jones potential between particle pairs.
double lj_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones force factor.
Definition lj.hpp:41
Routines to calculate the Lennard-Jones with cosine tail potential between particle pairs.
double ljcos2_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones cosine squared force factor.
Definition ljcos2.hpp:46
Routines to calculate the Lennard-Jones+cosine potential between particle pairs.
double ljcos_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones cosine force factor.
Definition ljcos.hpp:43
Routines to calculate the generalized Lennard-Jones potential between particle pairs.
double ljgen_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones force factor.
Definition ljgen.hpp:51
Routines to calculate the Morse potential between particle pairs.
double morse_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate Morse force factor.
Definition morse.hpp:41
VectorXd< 3 > Vector3d
Definition Vector.hpp:184
Various procedures concerning interactions between particles.
Routines to calculate the energy and/or force for particle pairs via interpolation of lookup tables.
double tabulated_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate a non-bonded pair force factor by linear interpolation from a table.
Routines to calculate the OIF global forces for a particle triple (triangle from mesh).
Routines to calculate the OIF local forces for a particle quadruple (two neighboring triangles with c...
Routines to calculate the smooth step potential between particle pairs.
double SmSt_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate smooth step force factor.
Routines to calculate the soft-sphere potential between particle pairs.
double soft_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate soft-sphere force factor.
Exception indicating that a bond type was unknown.
Solver::ShortRangeForceKernel kernel_type
Parameters for non-bonded interactions.
Force information on a particle.
Definition Particle.hpp:330
Utils::Vector3d torque
torque.
Definition Particle.hpp:360
Utils::Vector3d f
force.
Definition Particle.hpp:356
Struct holding all information for one particle.
Definition Particle.hpp:435
auto const & quat() const
Definition Particle.hpp:517
Routines to calculate the Thole damping potential between particle pairs.
Routines to calculate the Weeks-Chandler-Andersen potential between particle pairs.
double wca_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate WCA force factor.
Definition wca.hpp:40