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-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 * 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
76#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
78#endif
80 Utils::Vector3d const &d,
81 double const dist) {
82
83 auto force_factor = 0.;
84/* Lennard-Jones */
85#ifdef ESPRESSO_LENNARD_JONES
86 force_factor += lj_pair_force_factor(ia_params, dist);
87#endif
88/* WCA */
89#ifdef ESPRESSO_WCA
90 force_factor += wca_pair_force_factor(ia_params, dist);
91#endif
92/* Lennard-Jones generic */
93#ifdef ESPRESSO_LENNARD_JONES_GENERIC
94 force_factor += ljgen_pair_force_factor(ia_params, dist);
95#endif
96/* smooth step */
97#ifdef ESPRESSO_SMOOTH_STEP
98 force_factor += SmSt_pair_force_factor(ia_params, dist);
99#endif
100/* Hertzian force */
101#ifdef ESPRESSO_HERTZIAN
102 force_factor += hertzian_pair_force_factor(ia_params, dist);
103#endif
104/* Gaussian force */
105#ifdef ESPRESSO_GAUSSIAN
106 force_factor += gaussian_pair_force_factor(ia_params, dist);
107#endif
108/* BMHTF NaCl */
109#ifdef ESPRESSO_BMHTF_NACL
110 force_factor += BMHTF_pair_force_factor(ia_params, dist);
111#endif
112/* Buckingham*/
113#ifdef ESPRESSO_BUCKINGHAM
114 force_factor += buck_pair_force_factor(ia_params, dist);
115#endif
116/* Morse*/
117#ifdef ESPRESSO_MORSE
118 force_factor += morse_pair_force_factor(ia_params, dist);
119#endif
120/*soft-sphere potential*/
121#ifdef ESPRESSO_SOFT_SPHERE
122 force_factor += soft_pair_force_factor(ia_params, dist);
123#endif
124/*hat potential*/
125#ifdef ESPRESSO_HAT
126 force_factor += hat_pair_force_factor(ia_params, dist);
127#endif
128/* Lennard-Jones cosine */
129#ifdef ESPRESSO_LJCOS
130 force_factor += ljcos_pair_force_factor(ia_params, dist);
131#endif
132/* Lennard-Jones cosine */
133#ifdef ESPRESSO_LJCOS2
134 force_factor += ljcos2_pair_force_factor(ia_params, dist);
135#endif
136/* tabulated */
137#ifdef ESPRESSO_TABULATED
138 force_factor += tabulated_pair_force_factor(ia_params, dist);
139#endif
140 return force_factor * d;
141}
142
144 Particle const &p2,
145 IA_parameters const &ia_params,
146 Utils::Vector3d const &d,
147 double const dist) {
148
149 ParticleForce pf{};
150/* Gay-Berne */
151#ifdef ESPRESSO_GAY_BERNE
152 pf += gb_pair_force(p1.quat(), p2.quat(), ia_params, d, dist);
153#endif
154 return pf;
155}
156
158 Utils::Vector3d const &d) {
159 ParticleForce out{-pf.f};
160#ifdef ESPRESSO_ROTATION
161 // if torque is a null vector, the opposing torque is a null vector too
162 // (this check guards from returning a small yet non-null opposing
163 // torque due to numerical imprecision)
164 if (pf.torque[0] != 0. || pf.torque[1] != 0. || pf.torque[2] != 0.) {
165 out.torque = -(pf.torque + vector_product(d, pf.f));
166 }
167#endif
168 return out;
169}
170
171/**
172 * @brief For interactions which need particle information.
173 */
175 Particle &p1, Particle &p2, ParticleForce &pf, ParticleForce &p1f_asym,
176 ParticleForce &p2f_asym, Utils::Vector3d const &d, double dist,
177 double dist2, double q1q2, IA_parameters const &ia_params,
178 [[maybe_unused]] bool do_nonbonded_flag,
179 Thermostat::Thermostat const &thermostat, BoxGeometry const &box_geo,
180 [[maybe_unused]] BondedInteractionsMap const &bonded_ias,
181 [[maybe_unused]] Utils::Vector3d *const virial,
185 Coulomb::ShortRangeEnergyKernel::kernel_type const *coulomb_u_kernel) {
186
187 /***********************************************/
188 /* non-bonded pair potentials */
189 /***********************************************/
190
191 if (dist < ia_params.max_cut) {
192#ifdef ESPRESSO_EXCLUSIONS
193 if (do_nonbonded_flag) {
194#endif
195#ifdef ESPRESSO_THOLE
196 pf.f += thole_pair_force(p1, p2, ia_params, d, dist, bonded_ias,
197 coulomb_kernel);
198#endif
199 pf += calc_non_central_force(p1, p2, ia_params, d, dist);
200#ifdef ESPRESSO_EXCLUSIONS
201 }
202#endif
203 }
204
205 /*********************************************************************/
206 /* everything before this contributes to the virial pressure in NpT, */
207 /* but nothing afterwards, since the contribution to pressure from */
208 /* electrostatic is calculated by energy */
209 /*********************************************************************/
210#ifdef ESPRESSO_NPT
211 if (virial) {
212 *virial += hadamard_product(pf.f, d);
213 }
214#endif // ESPRESSO_NPT
215
216 /***********************************************/
217 /* short-range electrostatics */
218 /***********************************************/
219
220#ifdef ESPRESSO_ELECTROSTATICS
221 // real-space electrostatic charge-charge interaction
222 if (q1q2 != 0. and coulomb_kernel != nullptr) {
223 pf.f += (*coulomb_kernel)(q1q2, d, dist);
224#ifdef ESPRESSO_NPT
225 if (virial) {
226 (*virial)[0] += (*coulomb_u_kernel)(p1.pos(), p2.pos(), q1q2, d, dist);
227 }
228#endif // ESPRESSO_NPT
229 if (elc_kernel) {
230 (*elc_kernel)(p1.pos(), p2.pos(), p1f_asym.f, p2f_asym.f, q1q2);
231 }
232 }
233#endif // ESPRESSO_ELECTROSTATICS
234
235 /***********************************************/
236 /* thermostat */
237 /***********************************************/
238
239 /* The inter dpd force should not be part of the virial */
240#ifdef ESPRESSO_DPD
241 if (thermostat.thermo_switch & THERMO_DPD) {
242 auto const force =
243 dpd_pair_force(p1.pos(), p1.v(), p1.id(), p2.pos(), p2.v(), p2.id(),
244 *thermostat.dpd, box_geo, ia_params, d, dist, dist2);
245 pf += force;
246 }
247#endif
248
249 /***********************************************/
250 /* short-range magnetostatics */
251 /***********************************************/
252
253#ifdef ESPRESSO_DIPOLES
254 // real-space magnetic dipole-dipole
255 if (dipoles_kernel) {
256 auto const d1d2 = p1.dipm() * p2.dipm();
257 if (d1d2 != 0.) {
258 pf +=
259 (*dipoles_kernel)(d1d2, p1.calc_dip(), p2.calc_dip(), d, dist, dist2);
260 }
261 }
262#endif
263}
264
265/** Calculate non-bonded forces between a pair of particles and update their
266 * forces and torques.
267 * @param[in,out] p1 particle 1.
268 * @param[in,out] p2 particle 2.
269 * @param[in] d vector between @p p1 and @p p2.
270 * @param[in] dist distance between @p p1 and @p p2.
271 * @param[in] dist2 distance squared between @p p1 and @p p2.
272 * @param[in] q1q2 charge x charge between @p p1 and @p p2.
273 * @param[in] ia_params non-bonded interaction kernels.
274 * @param[in] thermostat thermostat.
275 * @param[in] box_geo box geometry.
276 * @param[in] bonded_ias bonded interaction kernels.
277 * @param[out] virial NpT virial.
278 * @param[in] coulomb_kernel Coulomb force kernel.
279 * @param[in] dipoles_kernel Dipolar force kernel.
280 * @param[in] elc_kernel ELC force correction kernel.
281 * @param[in] coulomb_u_kernel Coulomb energy kernel.
282 */
284 Particle &p1, Particle &p2, Utils::Vector3d const &d, double dist,
285 double dist2, double q1q2, IA_parameters const &ia_params,
286 Thermostat::Thermostat const &thermostat, BoxGeometry const &box_geo,
287 [[maybe_unused]] BondedInteractionsMap const &bonded_ias,
288 [[maybe_unused]] Utils::Vector3d *const virial,
292 Coulomb::ShortRangeEnergyKernel::kernel_type const *coulomb_u_kernel) {
293
294 ParticleForce pf{};
295 ParticleForce p1f_asym{};
296 ParticleForce p2f_asym{};
297
298#ifdef ESPRESSO_EXCLUSIONS
299 auto const do_nonbonded_flag = do_nonbonded(p1, p2);
300#else
301 auto constexpr do_nonbonded_flag = true;
302#endif
303
304 if (dist < ia_params.max_cut) {
305#ifdef ESPRESSO_EXCLUSIONS
306 if (do_nonbonded_flag) {
307#endif
308 pf.f += calc_central_radial_force(ia_params, d, dist);
309#ifdef ESPRESSO_EXCLUSIONS
310 }
311#endif
312 }
313
315 p1, p2, pf, p1f_asym, p2f_asym, d, dist, dist2, q1q2, ia_params,
316 do_nonbonded_flag, thermostat, box_geo, bonded_ias, virial,
317 coulomb_kernel, dipoles_kernel, elc_kernel, coulomb_u_kernel);
318
319 /***********************************************/
320 /* add total non-bonded forces to particles */
321 /***********************************************/
322
323 p1.force_and_torque() += pf + p1f_asym;
324 p2.force_and_torque() += calc_opposing_force(pf, d) + p2f_asym;
325}
326
327/** Compute the bonded interaction force between particle pairs.
328 *
329 * @param[in] iaparams Bonded parameters for the interaction.
330 * @param[in] p1 First particle.
331 * @param[in] p2 Second particle.
332 * @param[in] dx Vector between @p p1 and @p p2.
333 * @param[in] kernel Coulomb force kernel.
334 */
335inline std::optional<Utils::Vector3d> calc_bond_pair_force(
336 Bonded_IA_Parameters const &iaparams, Particle const &p1,
337 Particle const &p2, Utils::Vector3d const &dx,
339 if (auto const *iap = std::get_if<FeneBond>(&iaparams)) {
340 return iap->force(dx);
341 }
342 if (auto const *iap = std::get_if<HarmonicBond>(&iaparams)) {
343 return iap->force(dx);
344 }
345 if (auto const *iap = std::get_if<QuarticBond>(&iaparams)) {
346 return iap->force(dx);
347 }
348#ifdef ESPRESSO_ELECTROSTATICS
349 if (auto const *iap = std::get_if<BondedCoulomb>(&iaparams)) {
350 return iap->force(p1.q() * p2.q(), dx);
351 }
352 if (auto const *iap = std::get_if<BondedCoulombSR>(&iaparams)) {
353 return iap->force(dx, *kernel);
354 }
355#endif
356#ifdef ESPRESSO_BOND_CONSTRAINT
357 if (std::get_if<RigidBond>(&iaparams)) {
358 return Utils::Vector3d{};
359 }
360#endif
361#ifdef ESPRESSO_TABULATED
362 if (auto const *iap = std::get_if<TabulatedDistanceBond>(&iaparams)) {
363 return iap->force(dx);
364 }
365#endif
366 if (std::get_if<VirtualBond>(&iaparams)) {
367 return Utils::Vector3d{};
368 }
369 throw BondUnknownTypeError();
370}
371
373 Bonded_IA_Parameters const &iaparams, BoxGeometry const &box_geo,
374 Particle &p1, Particle &p2, [[maybe_unused]] Utils::Vector3d *const virial,
376 auto const dx = box_geo.get_mi_vector(p1.pos(), p2.pos());
377
378 if (auto const *iap = std::get_if<ThermalizedBond>(&iaparams)) {
379 auto result = iap->forces(p1, p2, dx);
380 if (result) {
381 auto const &forces = result.value();
382
383 p1.force() += std::get<0>(forces);
384 p2.force() += std::get<1>(forces);
385
386 return false;
387 }
388 } else {
389 auto result = calc_bond_pair_force(iaparams, p1, p2, dx, kernel);
390 if (result) {
391 p1.force() += result.value();
392 p2.force() -= result.value();
393
394#ifdef ESPRESSO_NPT
395 if (virial) {
396 *virial += hadamard_product(result.value(), dx);
397 }
398#endif
399 return false;
400 }
401 }
402 return true;
403}
404
405inline std::optional<
406 std::tuple<Utils::Vector3d, Utils::Vector3d, Utils::Vector3d>>
408 BoxGeometry const &box_geo, Particle const &p1,
409 Particle const &p2, Particle const &p3) {
410 auto const vec1 = box_geo.get_mi_vector(p2.pos(), p1.pos());
411 auto const vec2 = box_geo.get_mi_vector(p3.pos(), p1.pos());
412 if (auto const *iap = std::get_if<AngleHarmonicBond>(&iaparams)) {
413 return iap->forces(vec1, vec2);
414 }
415 if (auto const *iap = std::get_if<AngleCosineBond>(&iaparams)) {
416 return iap->forces(vec1, vec2);
417 }
418 if (auto const *iap = std::get_if<AngleCossquareBond>(&iaparams)) {
419 return iap->forces(vec1, vec2);
420 }
421#ifdef ESPRESSO_TABULATED
422 if (auto const *iap = std::get_if<TabulatedAngleBond>(&iaparams)) {
423 return iap->forces(vec1, vec2);
424 }
425#endif
426 if (auto const *iap = std::get_if<IBMTriel>(&iaparams)) {
427 return iap->calc_forces(vec1, vec2);
428 }
429 throw BondUnknownTypeError();
430}
431
433 BoxGeometry const &box_geo,
434 Particle &p1, Particle &p2,
435 Particle &p3) {
436 if (std::get_if<OifGlobalForcesBond>(&iaparams)) {
437 return false;
438 }
439 auto const result =
440 calc_bonded_three_body_force(iaparams, box_geo, p1, p2, p3);
441 if (result) {
442 auto const &forces = result.value();
443
444 p1.force() += std::get<0>(forces);
445 p2.force() += std::get<1>(forces);
446 p3.force() += std::get<2>(forces);
447
448 return false;
449 }
450 return true;
451}
452
453inline std::optional<std::tuple<Utils::Vector3d, Utils::Vector3d,
456 BoxGeometry const &box_geo, Particle const &p1,
457 Particle const &p2, Particle const &p3,
458 Particle const &p4) {
459 if (auto const *iap = std::get_if<OifLocalForcesBond>(&iaparams)) {
460 return iap->calc_forces(box_geo, p1, p2, p3, p4);
461 }
462 if (auto const *iap = std::get_if<IBMTribend>(&iaparams)) {
463 return iap->calc_forces(box_geo, p1, p2, p3, p4);
464 }
465 // note: particles in a dihedral bond are ordered as p2-p1-p3-p4
466 auto const v12 = box_geo.get_mi_vector(p1.pos(), p2.pos());
467 auto const v23 = box_geo.get_mi_vector(p3.pos(), p1.pos());
468 auto const v34 = box_geo.get_mi_vector(p4.pos(), p3.pos());
469 if (auto const *iap = std::get_if<DihedralBond>(&iaparams)) {
470 return iap->forces(v12, v23, v34);
471 }
472#ifdef ESPRESSO_TABULATED
473 if (auto const *iap = std::get_if<TabulatedDihedralBond>(&iaparams)) {
474 return iap->forces(v12, v23, v34);
475 }
476#endif
477 throw BondUnknownTypeError();
478}
479
481 BoxGeometry const &box_geo, Particle &p1,
482 Particle &p2, Particle &p3,
483 Particle &p4) {
484 auto const result =
485 calc_bonded_four_body_force(iaparams, box_geo, p1, p2, p3, p4);
486 if (result) {
487 auto const &forces = result.value();
488
489 p1.force() += std::get<0>(forces);
490 p2.force() += std::get<1>(forces);
491 p3.force() += std::get<2>(forces);
492 p4.force() += std::get<3>(forces);
493
494 return false;
495 }
496
497 return true;
498}
499
500inline bool
501add_bonded_force(Particle &p1, int bond_id, std::span<Particle *> partners,
502 BondedInteractionsMap const &bonded_ia_params,
503 BondBreakage::BondBreakage &bond_breakage,
504 BoxGeometry const &box_geo,
505 [[maybe_unused]] Utils::Vector3d *const virial,
507
508 // Consider for bond breakage
509 if (partners.size() == 1u) { // pair bonds
510 auto d = box_geo.get_mi_vector(p1.pos(), partners[0]->pos()).norm();
511 if (bond_breakage.check_and_handle_breakage(
512 p1.id(), {{partners[0]->id(), std::nullopt}}, bond_id, d)) {
513 return false;
514 }
515 }
516 if (partners.size() == 2u) { // angle bond
517 auto d =
518 box_geo.get_mi_vector(partners[0]->pos(), partners[1]->pos()).norm();
519 if (bond_breakage.check_and_handle_breakage(
520 p1.id(), {{partners[0]->id(), partners[1]->id()}}, bond_id, d)) {
521 return false;
522 }
523 }
524
525 auto const &iaparams = *bonded_ia_params.at(bond_id);
526
527 switch (number_of_partners(iaparams)) {
528 case 0:
529 return false;
530 case 1:
531 return add_bonded_two_body_force(iaparams, box_geo, p1, *partners[0],
532 virial, kernel);
533 case 2:
534 return add_bonded_three_body_force(iaparams, box_geo, p1, *partners[0],
535 *partners[1]);
536 case 3:
537 return add_bonded_four_body_force(iaparams, box_geo, p1, *partners[0],
538 *partners[1], *partners[2]);
539 default:
541 }
542}
#define ESPRESSO_ATTR_ALWAYS_INLINE
@ THERMO_DPD
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_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.
int number_of_partners(Bonded_IA_Parameters const &iaparams)
Get the number of bonded partners for the specified bond.
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.
bool check_and_handle_breakage(int particle_id, BondPartners const &bond_partners, int bond_type, double distance)
Check if the bond between the particles should break, if yes, queue it.
container for bonded interactions.
ESPRESSO_ATTR_ALWAYS_INLINE 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.
std::shared_ptr< DPDThermostat > dpd
int thermo_switch
Bitmask of currently active thermostats.
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)
Utils::Vector3d dpd_pair_force(Utils::Vector3d const &p1_position, Utils::Vector3d const &p1_velocity, int const &p1_id, Utils::Vector3d const &p2_position, Utils::Vector3d const &p2_velocity, int const &p2_id, DPDThermostat const &dpd, BoxGeometry const &box_geo, IA_parameters const &ia_params, Utils::Vector3d const &d, double dist, double dist2)
Definition dpd.cpp:79
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 ...
bool do_nonbonded(Particle const &p1, Particle const &p2)
Determine if the non-bonded interactions between p1 and p2 should be calculated.
void add_non_bonded_pair_force_with_p(Particle &p1, Particle &p2, ParticleForce &pf, ParticleForce &p1f_asym, ParticleForce &p2f_asym, Utils::Vector3d const &d, double dist, double dist2, double q1q2, IA_parameters const &ia_params, bool do_nonbonded_flag, Thermostat::Thermostat const &thermostat, BoxGeometry const &box_geo, BondedInteractionsMap const &bonded_ias, Utils::Vector3d *const virial, Coulomb::ShortRangeForceKernel::kernel_type const *coulomb_kernel, Dipoles::ShortRangeForceKernel::kernel_type const *dipoles_kernel, Coulomb::ShortRangeForceCorrectionsKernel::kernel_type const *elc_kernel, Coulomb::ShortRangeEnergyKernel::kernel_type const *coulomb_u_kernel)
For interactions which need particle information.
void add_non_bonded_pair_force(Particle &p1, Particle &p2, Utils::Vector3d const &d, double dist, double dist2, double q1q2, IA_parameters const &ia_params, Thermostat::Thermostat const &thermostat, BoxGeometry const &box_geo, BondedInteractionsMap const &bonded_ias, Utils::Vector3d *const virial, Coulomb::ShortRangeForceKernel::kernel_type const *coulomb_kernel, Dipoles::ShortRangeForceKernel::kernel_type const *dipoles_kernel, Coulomb::ShortRangeForceCorrectionsKernel::kernel_type const *elc_kernel, Coulomb::ShortRangeEnergyKernel::kernel_type const *coulomb_u_kernel)
Calculate non-bonded forces between a pair of particles and update their forces and torques.
std::optional< Utils::Vector3d > calc_bond_pair_force(Bonded_IA_Parameters const &iaparams, Particle const &p1, Particle const &p2, Utils::Vector3d const &dx, Coulomb::ShortRangeForceKernel::kernel_type const *kernel)
Compute the bonded interaction force between particle pairs.
ParticleForce calc_opposing_force(ParticleForce const &pf, Utils::Vector3d const &d)
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, Particle const &p1, Particle const &p2, Particle const &p3, Particle const &p4)
bool add_bonded_three_body_force(Bonded_IA_Parameters const &iaparams, BoxGeometry const &box_geo, Particle &p1, Particle &p2, Particle &p3)
bool add_bonded_four_body_force(Bonded_IA_Parameters const &iaparams, BoxGeometry const &box_geo, Particle &p1, Particle &p2, Particle &p3, Particle &p4)
ParticleForce calc_non_central_force(Particle const &p1, Particle const &p2, IA_parameters const &ia_params, Utils::Vector3d const &d, double const dist)
bool add_bonded_two_body_force(Bonded_IA_Parameters const &iaparams, BoxGeometry const &box_geo, Particle &p1, Particle &p2, Utils::Vector3d *const virial, Coulomb::ShortRangeForceKernel::kernel_type const *kernel)
std::optional< std::tuple< Utils::Vector3d, Utils::Vector3d, Utils::Vector3d > > calc_bonded_three_body_force(Bonded_IA_Parameters const &iaparams, BoxGeometry const &box_geo, Particle const &p1, Particle const &p2, Particle const &p3)
ESPRESSO_ATTR_ALWAYS_INLINE Utils::Vector3d calc_central_radial_force(IA_parameters const &ia_params, Utils::Vector3d const &d, double const dist)
bool add_bonded_force(Particle &p1, int bond_id, std::span< Particle * > partners, BondedInteractionsMap const &bonded_ia_params, BondBreakage::BondBreakage &bond_breakage, BoxGeometry const &box_geo, Utils::Vector3d *const virial, Coulomb::ShortRangeForceKernel::kernel_type const *kernel)
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:185
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 with an unexpected number of partners was encountered.
Exception indicating that a bond type was unknown.
Solver::ShortRangeEnergyKernel kernel_type
Solver::ShortRangeForceCorrectionsKernel kernel_type
Solver::ShortRangeForceKernel kernel_type
Solver::ShortRangeForceKernel kernel_type
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
Utils::Vector3d f
force.
Definition Particle.hpp:371
Struct holding all information for one particle.
Definition Particle.hpp:450
auto const & quat() const
Definition Particle.hpp:532
auto const & q() const
Definition Particle.hpp:593
auto const & force_and_torque() const
Definition Particle.hpp:492
auto calc_dip() const
Definition Particle.hpp:550
auto const & v() const
Definition Particle.hpp:488
auto const & dipm() const
Definition Particle.hpp:548
auto const & pos() const
Definition Particle.hpp:486
auto const & force() const
Definition Particle.hpp:490
auto const & id() const
Definition Particle.hpp:469
Routines to calculate the Thole damping potential between particle pairs.
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:45
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