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 "forces.hpp"
31
32#include "BoxGeometry.hpp"
33#include "actor/visitors.hpp"
60
61#ifdef DPD
62#include "dpd.hpp"
63#endif
64
65#include "Particle.hpp"
66#include "bond_error.hpp"
67#include "errorhandling.hpp"
68#include "exclusions.hpp"
69#include "thermostat.hpp"
70
71#include <utils/Span.hpp>
72#include <utils/Vector.hpp>
73
74#include <boost/optional.hpp>
75#include <boost/variant.hpp>
76
77#include <optional>
78#include <tuple>
79
81 Utils::Vector3d const &d,
82 double const dist) {
83
84 ParticleForce pf{};
85 auto force_factor = 0.;
86/* Lennard-Jones */
87#ifdef LENNARD_JONES
88 force_factor += lj_pair_force_factor(ia_params, dist);
89#endif
90/* WCA */
91#ifdef WCA
92 force_factor += wca_pair_force_factor(ia_params, dist);
93#endif
94/* Lennard-Jones generic */
95#ifdef LENNARD_JONES_GENERIC
96 force_factor += ljgen_pair_force_factor(ia_params, dist);
97#endif
98/* smooth step */
99#ifdef SMOOTH_STEP
100 force_factor += SmSt_pair_force_factor(ia_params, dist);
101#endif
102/* Hertzian force */
103#ifdef HERTZIAN
104 force_factor += hertzian_pair_force_factor(ia_params, dist);
105#endif
106/* Gaussian force */
107#ifdef GAUSSIAN
108 force_factor += gaussian_pair_force_factor(ia_params, dist);
109#endif
110/* BMHTF NaCl */
111#ifdef BMHTF_NACL
112 force_factor += BMHTF_pair_force_factor(ia_params, dist);
113#endif
114/* Buckingham*/
115#ifdef BUCKINGHAM
116 force_factor += buck_pair_force_factor(ia_params, dist);
117#endif
118/* Morse*/
119#ifdef MORSE
120 force_factor += morse_pair_force_factor(ia_params, dist);
121#endif
122/*soft-sphere potential*/
123#ifdef SOFT_SPHERE
124 force_factor += soft_pair_force_factor(ia_params, dist);
125#endif
126/*hat potential*/
127#ifdef HAT
128 force_factor += hat_pair_force_factor(ia_params, dist);
129#endif
130/* Lennard-Jones cosine */
131#ifdef LJCOS
132 force_factor += ljcos_pair_force_factor(ia_params, dist);
133#endif
134/* Lennard-Jones cosine */
135#ifdef LJCOS2
136 force_factor += ljcos2_pair_force_factor(ia_params, dist);
137#endif
138/* tabulated */
139#ifdef TABULATED
140 force_factor += tabulated_pair_force_factor(ia_params, dist);
141#endif
142 pf.f += force_factor * d;
143 return pf;
144}
145
147 Particle const &p1, Particle const &p2, IA_parameters const &ia_params,
148 Utils::Vector3d const &d, double const dist,
149 Coulomb::ShortRangeForceKernel::kernel_type const *coulomb_kernel) {
150
151 ParticleForce pf{};
152/* Thole damping */
153#ifdef THOLE
154 pf.f += thole_pair_force(p1, p2, ia_params, d, dist, coulomb_kernel);
155#endif
156
157 return pf;
158}
159
161 Particle const &p2,
162 IA_parameters const &ia_params,
163 Utils::Vector3d const &d,
164 double const dist) {
165
166 ParticleForce pf{};
167/* Gay-Berne */
168#ifdef GAY_BERNE
169 pf += gb_pair_force(p1.quat(), p2.quat(), ia_params, d, dist);
170#endif
171 return pf;
172}
173
175 Utils::Vector3d const &d) {
176 ParticleForce out{-pf.f};
177#ifdef ROTATION
178 // if torque is a null vector, the opposing torque is a null vector too
179 // (this check guards from returning a small yet non-null opposing
180 // torque due to numerical imprecision)
181 if (pf.torque[0] != 0. || pf.torque[1] != 0. || pf.torque[2] != 0.) {
182 out.torque = -(pf.torque + vector_product(d, pf.f));
183 }
184#endif
185 return out;
186}
187
188/** Calculate non-bonded forces between a pair of particles and update their
189 * forces and torques.
190 * @param[in,out] p1 particle 1.
191 * @param[in,out] p2 particle 2.
192 * @param[in] d vector between @p p1 and @p p2.
193 * @param[in] dist distance between @p p1 and @p p2.
194 * @param[in] dist2 distance squared between @p p1 and @p p2.
195 * @param[in] ia_params non-bonded interaction kernels.
196 * @param[in] thermostat thermostat.
197 * @param[in] box_geo box geometry.
198 * @param[in] coulomb_kernel Coulomb force kernel.
199 * @param[in] dipoles_kernel Dipolar force kernel.
200 * @param[in] elc_kernel ELC force correction kernel.
201 */
203 Particle &p1, Particle &p2, Utils::Vector3d const &d, double dist,
204 double dist2, IA_parameters const &ia_params,
205 Thermostat::Thermostat const &thermostat, BoxGeometry const &box_geo,
209
210 ParticleForce pf{};
211
212 /***********************************************/
213 /* non-bonded pair potentials */
214 /***********************************************/
215
216 if (dist < ia_params.max_cut) {
217#ifdef EXCLUSIONS
218 if (do_nonbonded(p1, p2)) {
219#endif
220 pf += calc_central_radial_force(ia_params, d, dist);
221 pf += calc_central_radial_charge_force(p1, p2, ia_params, d, dist,
222 coulomb_kernel);
223 pf += calc_non_central_force(p1, p2, ia_params, d, dist);
224#ifdef EXCLUSIONS
225 }
226#endif
227 }
228
229 /***********************************************/
230 /* short-range electrostatics */
231 /***********************************************/
232
233#ifdef ELECTROSTATICS
234 // real-space electrostatic charge-charge interaction
235 auto const q1q2 = p1.q() * p2.q();
236 if (q1q2 != 0. and coulomb_kernel != nullptr) {
237 pf.f += (*coulomb_kernel)(q1q2, d, dist);
238#ifdef P3M
239 if (elc_kernel)
240 (*elc_kernel)(p1, p2, q1q2);
241#endif // P3M
242 }
243#endif // ELECTROSTATICS
244
245 /*********************************************************************/
246 /* everything before this contributes to the virial pressure in NpT, */
247 /* but nothing afterwards */
248 /*********************************************************************/
249#ifdef NPT
251#endif
252
253 /***********************************************/
254 /* thermostat */
255 /***********************************************/
256
257 /* The inter dpd force should not be part of the virial */
258#ifdef DPD
259 if (thermostat.thermo_switch & THERMO_DPD) {
260 auto const force = dpd_pair_force(p1, p2, *thermostat.dpd, box_geo,
261 ia_params, d, dist, dist2);
262 p1.force() += force;
263 p2.force() -= force;
264 }
265#endif
266
267 /***********************************************/
268 /* short-range magnetostatics */
269 /***********************************************/
270
271#ifdef DIPOLES
272 // real-space magnetic dipole-dipole
273 if (dipoles_kernel) {
274 pf += (*dipoles_kernel)(p1, p2, d, dist, dist2);
275 }
276#endif
277
278 /***********************************************/
279 /* add total non-bonded forces to particles */
280 /***********************************************/
281
282 p1.force_and_torque() += pf;
284}
285
286/** Compute the bonded interaction force between particle pairs.
287 *
288 * @param[in] p1 First particle.
289 * @param[in] p2 Second particle.
290 * @param[in] iaparams Bonded parameters for the interaction.
291 * @param[in] dx Vector between @p p1 and @p p2.
292 * @param[in] kernel Coulomb force kernel.
293 */
294inline boost::optional<Utils::Vector3d> calc_bond_pair_force(
295 Particle const &p1, Particle const &p2,
296 Bonded_IA_Parameters const &iaparams, Utils::Vector3d const &dx,
298 if (auto const *iap = boost::get<FeneBond>(&iaparams)) {
299 return iap->force(dx);
300 }
301 if (auto const *iap = boost::get<HarmonicBond>(&iaparams)) {
302 return iap->force(dx);
303 }
304 if (auto const *iap = boost::get<QuarticBond>(&iaparams)) {
305 return iap->force(dx);
306 }
307#ifdef ELECTROSTATICS
308 if (auto const *iap = boost::get<BondedCoulomb>(&iaparams)) {
309 return iap->force(p1.q() * p2.q(), dx);
310 }
311 if (auto const *iap = boost::get<BondedCoulombSR>(&iaparams)) {
312 return iap->force(dx, *kernel);
313 }
314#endif
315#ifdef BOND_CONSTRAINT
316 if (boost::get<RigidBond>(&iaparams)) {
317 return Utils::Vector3d{};
318 }
319#endif
320#ifdef TABULATED
321 if (auto const *iap = boost::get<TabulatedDistanceBond>(&iaparams)) {
322 return iap->force(dx);
323 }
324#endif
325 if (boost::get<VirtualBond>(&iaparams)) {
326 return Utils::Vector3d{};
327 }
328 throw BondUnknownTypeError();
329}
330
332 Bonded_IA_Parameters const &iaparams, Particle &p1, Particle &p2,
333 BoxGeometry const &box_geo,
335 auto const dx = box_geo.get_mi_vector(p1.pos(), p2.pos());
336
337 if (auto const *iap = boost::get<ThermalizedBond>(&iaparams)) {
338 auto result = iap->forces(p1, p2, dx);
339 if (result) {
340 auto const &forces = result.get();
341
342 p1.force() += std::get<0>(forces);
343 p2.force() += std::get<1>(forces);
344
345 return false;
346 }
347 } else {
348 auto result = calc_bond_pair_force(p1, p2, iaparams, dx, kernel);
349 if (result) {
350 p1.force() += result.get();
351 p2.force() -= result.get();
352
353#ifdef NPT
354 npt_add_virial_force_contribution(result.get(), dx);
355#endif
356 return false;
357 }
358 }
359 return true;
360}
361
362inline boost::optional<
363 std::tuple<Utils::Vector3d, Utils::Vector3d, Utils::Vector3d>>
365 BoxGeometry const &box_geo, Particle const &p1,
366 Particle const &p2, Particle const &p3) {
367 auto const vec1 = box_geo.get_mi_vector(p2.pos(), p1.pos());
368 auto const vec2 = box_geo.get_mi_vector(p3.pos(), p1.pos());
369 if (auto const *iap = boost::get<AngleHarmonicBond>(&iaparams)) {
370 return iap->forces(vec1, vec2);
371 }
372 if (auto const *iap = boost::get<AngleCosineBond>(&iaparams)) {
373 return iap->forces(vec1, vec2);
374 }
375 if (auto const *iap = boost::get<AngleCossquareBond>(&iaparams)) {
376 return iap->forces(vec1, vec2);
377 }
378#ifdef TABULATED
379 if (auto const *iap = boost::get<TabulatedAngleBond>(&iaparams)) {
380 return iap->forces(vec1, vec2);
381 }
382#endif
383 if (auto const *iap = boost::get<IBMTriel>(&iaparams)) {
384 return iap->calc_forces(vec1, vec2);
385 }
386 throw BondUnknownTypeError();
387}
388
390 BoxGeometry const &box_geo,
391 Particle &p1, Particle &p2,
392 Particle &p3) {
393 if (boost::get<OifGlobalForcesBond>(&iaparams)) {
394 return false;
395 }
396 auto const result =
397 calc_bonded_three_body_force(iaparams, box_geo, p1, p2, p3);
398 if (result) {
399 auto const &forces = result.get();
400
401 p1.force() += std::get<0>(forces);
402 p2.force() += std::get<1>(forces);
403 p3.force() += std::get<2>(forces);
404
405 return false;
406 }
407 return true;
408}
409
410inline boost::optional<std::tuple<Utils::Vector3d, Utils::Vector3d,
413 BoxGeometry const &box_geo, Particle const &p1,
414 Particle const &p2, Particle const &p3,
415 Particle const &p4) {
416 if (auto const *iap = boost::get<OifLocalForcesBond>(&iaparams)) {
417 return iap->calc_forces(box_geo, p1, p2, p3, p4);
418 }
419 if (auto const *iap = boost::get<IBMTribend>(&iaparams)) {
420 return iap->calc_forces(box_geo, p1, p2, p3, p4);
421 }
422 // note: particles in a dihedral bond are ordered as p2-p1-p3-p4
423 auto const v12 = box_geo.get_mi_vector(p1.pos(), p2.pos());
424 auto const v23 = box_geo.get_mi_vector(p3.pos(), p1.pos());
425 auto const v34 = box_geo.get_mi_vector(p4.pos(), p3.pos());
426 if (auto const *iap = boost::get<DihedralBond>(&iaparams)) {
427 return iap->forces(v12, v23, v34);
428 }
429#ifdef TABULATED
430 if (auto const *iap = boost::get<TabulatedDihedralBond>(&iaparams)) {
431 return iap->forces(v12, v23, v34);
432 }
433#endif
434 throw BondUnknownTypeError();
435}
436
438 BoxGeometry const &box_geo, Particle &p1,
439 Particle &p2, Particle &p3,
440 Particle &p4) {
441 auto const result =
442 calc_bonded_four_body_force(iaparams, box_geo, p1, p2, p3, p4);
443 if (result) {
444 auto const &forces = result.get();
445
446 p1.force() += std::get<0>(forces);
447 p2.force() += std::get<1>(forces);
448 p3.force() += std::get<2>(forces);
449 p4.force() += std::get<3>(forces);
450
451 return false;
452 }
453
454 return true;
455}
456
457inline bool
459 BondBreakage::BondBreakage &bond_breakage,
460 BoxGeometry const &box_geo,
462
463 // Consider for bond breakage
464 if (partners.size() == 1) { // pair bonds
465 auto d = box_geo.get_mi_vector(p1.pos(), partners[0]->pos()).norm();
466 if (bond_breakage.check_and_handle_breakage(
467 p1.id(), {{partners[0]->id(), std::nullopt}}, bond_id, d)) {
468 return false;
469 }
470 }
471 if (partners.size() == 2) { // angle bond
472 auto d =
473 box_geo.get_mi_vector(partners[0]->pos(), partners[1]->pos()).norm();
474 if (bond_breakage.check_and_handle_breakage(
475 p1.id(), {{partners[0]->id(), partners[1]->id()}}, bond_id, d)) {
476 return false;
477 }
478 }
479
480 auto const &iaparams = *bonded_ia_params.at(bond_id);
481
482 switch (number_of_partners(iaparams)) {
483 case 0:
484 return false;
485 case 1:
486 return add_bonded_two_body_force(iaparams, p1, *partners[0], box_geo,
487 kernel);
488 case 2:
489 return add_bonded_three_body_force(iaparams, box_geo, p1, *partners[0],
490 *partners[1]);
491 case 3:
492 return add_bonded_four_body_force(iaparams, box_geo, p1, *partners[0],
493 *partners[1], *partners[2]);
494 default:
496 }
497}
@ THERMO_DPD
Vector implementation and trait types for boost qvm interoperability.
__global__ float * force
__shared__ int pos[MAXDEPTH *THREADS5/WARPSIZE]
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.
BondedInteractionsMap bonded_ia_params
Field containing the parameters of the bonded ia types.
Data structures for bonded interactions.
int number_of_partners(Bonded_IA_Parameters const &iaparams)
Return the number of bonded partners for the specified bond.
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_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.
mapped_type at(key_type const &key) const
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.
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
This file contains the defaults for ESPResSo.
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(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 .
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 npt_add_virial_force_contribution(const Utils::Vector3d &force, const Utils::Vector3d &d)
Update the NpT virial.
Definition forces.cpp:269
Force calculation.
boost::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)
bool add_bonded_force(Particle &p1, int bond_id, Utils::Span< Particle * > partners, BondBreakage::BondBreakage &bond_breakage, BoxGeometry const &box_geo, Coulomb::ShortRangeForceKernel::kernel_type const *kernel)
ParticleForce calc_central_radial_force(IA_parameters const &ia_params, Utils::Vector3d const &d, double const dist)
void add_non_bonded_pair_force(Particle &p1, Particle &p2, Utils::Vector3d const &d, double dist, double dist2, IA_parameters const &ia_params, Thermostat::Thermostat const &thermostat, BoxGeometry const &box_geo, Coulomb::ShortRangeForceKernel::kernel_type const *coulomb_kernel, Dipoles::ShortRangeForceKernel::kernel_type const *dipoles_kernel, Coulomb::ShortRangeForceCorrectionsKernel::kernel_type const *elc_kernel)
Calculate non-bonded forces between a pair of particles and update their forces and torques.
ParticleForce calc_opposing_force(ParticleForce const &pf, Utils::Vector3d const &d)
bool add_bonded_two_body_force(Bonded_IA_Parameters const &iaparams, Particle &p1, Particle &p2, BoxGeometry const &box_geo, Coulomb::ShortRangeForceKernel::kernel_type const *kernel)
boost::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)
boost::optional< Utils::Vector3d > calc_bond_pair_force(Particle const &p1, Particle const &p2, Bonded_IA_Parameters const &iaparams, Utils::Vector3d const &dx, Coulomb::ShortRangeForceKernel::kernel_type const *kernel)
Compute the bonded interaction force between particle pairs.
ParticleForce calc_central_radial_charge_force(Particle const &p1, Particle const &p2, IA_parameters const &ia_params, Utils::Vector3d const &d, double const dist, Coulomb::ShortRangeForceKernel::kernel_type const *coulomb_kernel)
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)
Routines to calculate the Gay-Berne potential between particle pairs.
ParticleForce gb_pair_force(Utils::Quaternion< double > const &qi, Utils::Quaternion< double > const &qj, 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:42
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:157
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 energy or/and and force for a particle triple (triangle f...
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::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:290
Utils::Vector3d torque
torque.
Definition Particle.hpp:318
Utils::Vector3d f
force.
Definition Particle.hpp:314
Struct holding all information for one particle.
Definition Particle.hpp:393
auto const & quat() const
Definition Particle.hpp:475
auto const & q() const
Definition Particle.hpp:506
auto const & force_and_torque() const
Definition Particle.hpp:435
auto const & pos() const
Definition Particle.hpp:429
auto const & force() const
Definition Particle.hpp:433
auto const & id() const
Definition Particle.hpp:412
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, Coulomb::ShortRangeForceKernel::kernel_type const *kernel)
Calculate Thole force.
Definition thole.hpp:44
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:39