ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
energy_inline.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2010-2022 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 * Energy calculation.
26 */
27
28#include "config/config.hpp"
29
50
51#include "BoxGeometry.hpp"
52#include "Observable_stat.hpp"
53#include "Particle.hpp"
54#include "bond_error.hpp"
55#include "errorhandling.hpp"
56#include "exclusions.hpp"
57
58#include <utils/Span.hpp>
59#include <utils/Vector.hpp>
60
61#include <boost/optional.hpp>
62#include <boost/range/algorithm/find_if.hpp>
63#include <boost/variant.hpp>
64
65/** Calculate non-bonded energies between a pair of particles.
66 * @param p1 particle 1.
67 * @param p2 particle 2.
68 * @param ia_params the interaction parameters between the two particles
69 * @param d vector between p1 and p2.
70 * @param dist distance between p1 and p2.
71 * @param coulomb_kernel Coulomb energy kernel.
72 * @return the short-range interaction energy between the two particles
73 */
75 Particle const &p1, Particle const &p2, IA_parameters const &ia_params,
76 Utils::Vector3d const &d, double const dist,
78
79 double ret = 0;
80
81#ifdef LENNARD_JONES
82 /* Lennard-Jones */
83 ret += lj_pair_energy(ia_params, dist);
84#endif
85#ifdef WCA
86 /* WCA */
87 ret += wca_pair_energy(ia_params, dist);
88#endif
89
90#ifdef LENNARD_JONES_GENERIC
91 /* Generic Lennard-Jones */
92 ret += ljgen_pair_energy(ia_params, dist);
93#endif
94
95#ifdef SMOOTH_STEP
96 /* smooth step */
97 ret += SmSt_pair_energy(ia_params, dist);
98#endif
99
100#ifdef HERTZIAN
101 /* Hertzian potential */
102 ret += hertzian_pair_energy(ia_params, dist);
103#endif
104
105#ifdef GAUSSIAN
106 /* Gaussian potential */
107 ret += gaussian_pair_energy(ia_params, dist);
108#endif
109
110#ifdef BMHTF_NACL
111 /* BMHTF NaCl */
112 ret += BMHTF_pair_energy(ia_params, dist);
113#endif
114
115#ifdef MORSE
116 /* Morse */
117 ret += morse_pair_energy(ia_params, dist);
118#endif
119
120#ifdef BUCKINGHAM
121 /* Buckingham */
122 ret += buck_pair_energy(ia_params, dist);
123#endif
124
125#ifdef SOFT_SPHERE
126 /* soft-sphere */
127 ret += soft_pair_energy(ia_params, dist);
128#endif
129
130#ifdef HAT
131 /* hat */
132 ret += hat_pair_energy(ia_params, dist);
133#endif
134
135#ifdef LJCOS2
136 /* Lennard-Jones */
137 ret += ljcos2_pair_energy(ia_params, dist);
138#endif
139
140#ifdef THOLE
141 /* Thole damping */
142 ret += thole_pair_energy(p1, p2, ia_params, d, dist, coulomb_kernel);
143#endif
144
145#ifdef TABULATED
146 /* tabulated */
147 ret += tabulated_pair_energy(ia_params, dist);
148#endif
149
150#ifdef LJCOS
151 /* Lennard-Jones cosine */
152 ret += ljcos_pair_energy(ia_params, dist);
153#endif
154
155#ifdef GAY_BERNE
156 /* Gay-Berne */
157 ret += gb_pair_energy(p1.quat(), p2.quat(), ia_params, d, dist);
158#endif
159
160 return ret;
161}
162
163/** Add non-bonded and short-range Coulomb energies between a pair of particles
164 * to the energy observable.
165 * @param p1 particle 1.
166 * @param p2 particle 2.
167 * @param d vector between p1 and p2.
168 * @param dist distance between p1 and p2.
169 * @param dist2 distance squared between p1 and p2.
170 * @param[in] ia_params non-bonded interaction kernels.
171 * @param[in] coulomb_kernel Coulomb energy kernel.
172 * @param[in] dipoles_kernel Dipolar energy kernel.
173 * @param[in,out] obs_energy energy observable.
174 */
176 Particle const &p1, Particle const &p2, Utils::Vector3d const &d,
177 double const dist, double const dist2, IA_parameters const &ia_params,
180 Observable_stat &obs_energy) {
181
182#ifdef EXCLUSIONS
183 if (do_nonbonded(p1, p2))
184#endif
186 p1.type(), p2.type(), p1.mol_id(), p2.mol_id(),
187 calc_non_bonded_pair_energy(p1, p2, ia_params, d, dist,
188 coulomb_kernel));
189
190#ifdef ELECTROSTATICS
191 if (!obs_energy.coulomb.empty() and coulomb_kernel != nullptr) {
192 auto const q1q2 = p1.q() * p2.q();
193 obs_energy.coulomb[0] += (*coulomb_kernel)(p1, p2, q1q2, d, dist);
194 }
195#endif
196
197#ifdef DIPOLES
198 if (!obs_energy.dipolar.empty() and dipoles_kernel != nullptr)
199 obs_energy.dipolar[0] += (*dipoles_kernel)(p1, p2, d, dist, dist2);
200#endif
201}
202
203inline boost::optional<double>
205 Utils::Span<Particle *> partners, BoxGeometry const &box_geo,
207 auto const n_partners = static_cast<int>(partners.size());
208
209 auto p2 = (n_partners > 0) ? partners[0] : nullptr;
210 auto p3 = (n_partners > 1) ? partners[1] : nullptr;
211 auto p4 = (n_partners > 2) ? partners[2] : nullptr;
212
213 if (n_partners == 1) {
214 auto const dx = box_geo.get_mi_vector(p1.pos(), p2->pos());
215 if (auto const *iap = boost::get<FeneBond>(&iaparams)) {
216 return iap->energy(dx);
217 }
218 if (auto const *iap = boost::get<HarmonicBond>(&iaparams)) {
219 return iap->energy(dx);
220 }
221 if (auto const *iap = boost::get<QuarticBond>(&iaparams)) {
222 return iap->energy(dx);
223 }
224#ifdef ELECTROSTATICS
225 if (auto const *iap = boost::get<BondedCoulomb>(&iaparams)) {
226 return iap->energy(p1.q() * p2->q(), dx);
227 }
228 if (auto const *iap = boost::get<BondedCoulombSR>(&iaparams)) {
229 return iap->energy(p1, *p2, dx, *kernel);
230 }
231#endif
232#ifdef BOND_CONSTRAINT
233 if (boost::get<RigidBond>(&iaparams)) {
234 return {0.};
235 }
236#endif
237#ifdef TABULATED
238 if (auto const *iap = boost::get<TabulatedDistanceBond>(&iaparams)) {
239 return iap->energy(dx);
240 }
241#endif
242 if (boost::get<VirtualBond>(&iaparams)) {
243 return {0.};
244 }
245 throw BondUnknownTypeError();
246 } // 1 partner
247 if (n_partners == 2) {
248 auto const vec1 = box_geo.get_mi_vector(p2->pos(), p1.pos());
249 auto const vec2 = box_geo.get_mi_vector(p3->pos(), p1.pos());
250 if (auto const *iap = boost::get<AngleHarmonicBond>(&iaparams)) {
251 return iap->energy(vec1, vec2);
252 }
253 if (auto const *iap = boost::get<AngleCosineBond>(&iaparams)) {
254 return iap->energy(vec1, vec2);
255 }
256 if (auto const *iap = boost::get<AngleCossquareBond>(&iaparams)) {
257 return iap->energy(vec1, vec2);
258 }
259 if (auto const *iap = boost::get<TabulatedAngleBond>(&iaparams)) {
260 return iap->energy(vec1, vec2);
261 }
262 if (boost::get<IBMTriel>(&iaparams)) {
263 runtimeWarningMsg() << "Unsupported bond type " +
264 std::to_string(iaparams.which()) +
265 " in energy calculation.";
266 return 0.;
267 }
268 throw BondUnknownTypeError();
269 } // 2 partners
270 if (n_partners == 3) {
271 // note: particles in a dihedral bond are ordered as p2-p1-p3-p4
272 auto const v12 = box_geo.get_mi_vector(p1.pos(), p2->pos());
273 auto const v23 = box_geo.get_mi_vector(p3->pos(), p1.pos());
274 auto const v34 = box_geo.get_mi_vector(p4->pos(), p3->pos());
275 if (auto const *iap = boost::get<DihedralBond>(&iaparams)) {
276 return iap->energy(v12, v23, v34);
277 }
278 if (auto const *iap = boost::get<TabulatedDihedralBond>(&iaparams)) {
279 return iap->energy(v12, v23, v34);
280 }
281 if (boost::get<IBMTribend>(&iaparams)) {
282 runtimeWarningMsg() << "Unsupported bond type " +
283 std::to_string(iaparams.which()) +
284 " in energy calculation.";
285 return 0.;
286 }
287 throw BondUnknownTypeError();
288 } // 3 partners
289 if (n_partners == 0) {
290 return 0.;
291 }
292
293 throw BondInvalidSizeError(n_partners);
294}
295
296/** Calculate kinetic energies from translation for one particle.
297 * @param p particle for which to calculate energies
298 */
299inline double translational_kinetic_energy(Particle const &p) {
300 return p.is_virtual() ? 0. : 0.5 * p.mass() * p.v().norm2();
301}
302
303/** Calculate kinetic energies from rotation for one particle.
304 * @param p particle for which to calculate energies
305 */
306inline double rotational_kinetic_energy(Particle const &p) {
307#ifdef ROTATION
308 return (p.can_rotate() and not p.is_virtual())
309 ? 0.5 * (hadamard_product(p.omega(), p.omega()) * p.rinertia())
310 : 0.0;
311#else
312 return 0.0;
313#endif
314}
315
316/** Calculate kinetic energies for one particle.
317 * @param p particle for which to calculate energies
318 */
Vector implementation and trait types for boost qvm interoperability.
Routines to calculate the Born-Meyer-Huggins-Tosi-Fumi potential between particle pairs.
double BMHTF_pair_energy(IA_parameters const &ia_params, double dist)
Calculate BMHTF potential energy.
Data structures for bonded interactions.
boost::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_energy(IA_parameters const &ia_params, double dist)
Calculate Buckingham energy.
Utils::Vector< T, 3 > get_mi_vector(const Utils::Vector< T, 3 > &a, const Utils::Vector< T, 3 > &b) const
Get the minimum-image vector between two coordinates.
Observable for the pressure and energy.
Utils::Span< double > dipolar
Contribution(s) from dipolar interactions.
void add_non_bonded_contribution(int type1, int type2, int molid1, int molid2, Utils::Span< const double > data)
Utils::Span< double > coulomb
Contribution(s) from Coulomb interactions.
A stripped-down version of std::span from C++17.
Definition Span.hpp:38
DEVICE_QUALIFIER constexpr size_type size() const
Definition Span.hpp:86
DEVICE_QUALIFIER constexpr bool empty() const
Definition Span.hpp:87
This file contains the defaults for ESPResSo.
Routines to calculate the Gaussian potential between particle pairs.
double gaussian_pair_energy(IA_parameters const &ia_params, double dist)
Calculate Gaussian energy.
double calc_non_bonded_pair_energy(Particle const &p1, Particle const &p2, IA_parameters const &ia_params, Utils::Vector3d const &d, double const dist, Coulomb::ShortRangeEnergyKernel::kernel_type const *coulomb_kernel)
Calculate non-bonded energies between a pair of particles.
double translational_kinetic_energy(Particle const &p)
Calculate kinetic energies from translation for one particle.
double calc_kinetic_energy(Particle const &p)
Calculate kinetic energies for one particle.
double rotational_kinetic_energy(Particle const &p)
Calculate kinetic energies from rotation for one particle.
void add_non_bonded_pair_energy(Particle const &p1, Particle const &p2, Utils::Vector3d const &d, double const dist, double const dist2, IA_parameters const &ia_params, Coulomb::ShortRangeEnergyKernel::kernel_type const *coulomb_kernel, Dipoles::ShortRangeEnergyKernel::kernel_type const *dipoles_kernel, Observable_stat &obs_energy)
Add non-bonded and short-range Coulomb energies between a pair of particles to the energy observable.
boost::optional< double > calc_bonded_energy(Bonded_IA_Parameters const &iaparams, Particle const &p1, Utils::Span< Particle * > partners, BoxGeometry const &box_geo, Coulomb::ShortRangeEnergyKernel::kernel_type const *kernel)
This file contains the errorhandling code for severe errors, like a broken bond or illegal parameter ...
#define runtimeWarningMsg()
bool do_nonbonded(Particle const &p1, Particle const &p2)
Determine if the non-bonded interactions between p1 and p2 should be calculated.
Routines to calculate the Gay-Berne potential between particle pairs.
double gb_pair_energy(Utils::Quaternion< double > const &qi, Utils::Quaternion< double > const &qj, IA_parameters const &ia_params, Utils::Vector3d const &d, double dist)
Calculate Gay-Berne energy.
Routines to calculate the hat potential between particle pairs.
double hat_pair_energy(IA_parameters const &ia_params, double dist)
Calculate hat energy.
Definition hat.hpp:46
Routines to calculate the Hertzian potential between particle pairs.
double hertzian_pair_energy(IA_parameters const &ia_params, double dist)
Calculate Hertzian energy.
Definition hertzian.hpp:49
Routines to calculate the Lennard-Jones potential between particle pairs.
double lj_pair_energy(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones energy.
Definition lj.hpp:52
Routines to calculate the Lennard-Jones with cosine tail potential between particle pairs.
double ljcos2_pair_energy(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones cosine squared energy.
Definition ljcos2.hpp:67
Routines to calculate the Lennard-Jones+cosine potential between particle pairs.
double ljcos_pair_energy(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones cosine energy.
Definition ljcos.hpp:63
Routines to calculate the generalized Lennard-Jones potential between particle pairs.
double ljgen_pair_energy(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones energy.
Definition ljgen.hpp:81
Routines to calculate the Morse potential between particle pairs.
double morse_pair_energy(IA_parameters const &ia_params, double dist)
Calculate Morse energy.
Definition morse.hpp:53
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_energy(IA_parameters const &ia_params, double dist)
Calculate a non-bonded pair energy by linear interpolation from a table.
Routines to calculate the smooth step potential between particle pairs.
double SmSt_pair_energy(IA_parameters const &ia_params, double dist)
Calculate smooth step energy.
Routines to calculate the soft-sphere potential between particle pairs.
double soft_pair_energy(IA_parameters const &ia_params, double dist)
Calculate soft-sphere energy.
Exception indicating that a bond with an unexpected number of partners was encountered.
Exception indicating that a bond type was unknown.
Solver::ShortRangeEnergyKernel kernel_type
Solver::ShortRangeEnergyKernel kernel_type
Parameters for non-bonded interactions.
Struct holding all information for one particle.
Definition Particle.hpp:393
bool can_rotate() const
Definition Particle.hpp:458
auto const & rinertia() const
Definition Particle.hpp:500
auto is_virtual() const
Definition Particle.hpp:516
auto const & mass() const
Definition Particle.hpp:450
auto const & quat() const
Definition Particle.hpp:475
auto const & q() const
Definition Particle.hpp:506
auto const & v() const
Definition Particle.hpp:431
auto const & omega() const
Definition Particle.hpp:479
auto const & type() const
Definition Particle.hpp:416
auto const & mol_id() const
Definition Particle.hpp:414
auto const & pos() const
Definition Particle.hpp:429
Routines to calculate the Thole damping potential between particle pairs.
double thole_pair_energy(Particle const &p1, Particle const &p2, IA_parameters const &ia_params, Utils::Vector3d const &d, double dist, Coulomb::ShortRangeEnergyKernel::kernel_type const *kernel)
Calculate Thole energy.
Definition thole.hpp:67
Routines to calculate the Weeks-Chandler-Andersen potential between particle pairs.
double wca_pair_energy(IA_parameters const &ia_params, double dist)
Calculate WCA energy.
Definition wca.hpp:49