ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
ljgen.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#ifndef CORE_NB_IA_LJGEN_HPP
22#define CORE_NB_IA_LJGEN_HPP
23
24/** \file
25 * Routines to calculate the generalized Lennard-Jones potential between
26 * particle pairs. "Generalized" here means that the LJ energy is of the form
27 * @f[
28 * \varepsilon \cdot
29 * \left[
30 * b_1 \left(\frac{\sigma}{r-r_{\text{offset}}}\right)^{a_1}
31 * - b_2 \left(\frac{\sigma}{r-r_{\text{offset}}}\right)^{a_2}
32 * + \text{shift}
33 * \right]
34 * @f]
35 *
36 * Implementation in \ref ljgen.cpp.
37 */
38
39#include "config/config.hpp"
40
41#ifdef LENNARD_JONES_GENERIC
42
44
45#include <utils/Vector.hpp>
46#include <utils/math/sqr.hpp>
47
48#include <cmath>
49
50/** Calculate Lennard-Jones force factor */
51inline double ljgen_pair_force_factor(IA_parameters const &ia_params,
52 double dist) {
53 if (dist < ia_params.ljgen.max_cutoff()) {
54 auto r_off = dist - ia_params.ljgen.offset;
55
56#ifdef LJGEN_SOFTCORE
57 r_off *= r_off;
58 r_off += Utils::sqr(ia_params.ljgen.sig) * (1.0 - ia_params.ljgen.lambda) *
59 ia_params.ljgen.softrad;
60 r_off = sqrt(r_off);
61#else
62 r_off = std::abs(r_off);
63#endif
64 auto const frac = ia_params.ljgen.sig / r_off;
65 auto const fac = ia_params.ljgen.eps
66#ifdef LJGEN_SOFTCORE
67 * ia_params.ljgen.lambda *
68 (dist - ia_params.ljgen.offset) / r_off
69#endif
70 * (ia_params.ljgen.b1 * ia_params.ljgen.a1 *
71 pow(frac, ia_params.ljgen.a1) -
72 ia_params.ljgen.b2 * ia_params.ljgen.a2 *
73 pow(frac, ia_params.ljgen.a2)) /
74 (r_off * dist);
75 return fac;
76 }
77 return 0.0;
78}
79
80/** Calculate Lennard-Jones energy */
81inline double ljgen_pair_energy(IA_parameters const &ia_params, double dist) {
82 if (dist < ia_params.ljgen.max_cutoff()) {
83 auto r_off = dist - ia_params.ljgen.offset;
84#ifdef LJGEN_SOFTCORE
85 r_off *= r_off;
86 r_off += pow(ia_params.ljgen.sig, 2) * (1.0 - ia_params.ljgen.lambda) *
87 ia_params.ljgen.softrad;
88 r_off = sqrt(r_off);
89#else
90 r_off = std::abs(r_off);
91#endif
92 auto const frac = ia_params.ljgen.sig / r_off;
93 auto const fac = ia_params.ljgen.eps
94#ifdef LJGEN_SOFTCORE
95 * ia_params.ljgen.lambda
96#endif
97 * (ia_params.ljgen.b1 * pow(frac, ia_params.ljgen.a1) -
98 ia_params.ljgen.b2 * pow(frac, ia_params.ljgen.a2) +
99 ia_params.ljgen.shift);
100 return fac;
101 }
102 return 0.0;
103}
104
105#endif // LENNARD_JONES_GENERIC
106#endif
Vector implementation and trait types for boost qvm interoperability.
This file contains the defaults for ESPResSo.
double ljgen_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones force factor.
Definition ljgen.hpp:51
double ljgen_pair_energy(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones energy.
Definition ljgen.hpp:81
DEVICE_QUALIFIER constexpr T sqr(T x)
Calculates the SQuaRe of x.
Definition sqr.hpp:28
Various procedures concerning interactions between particles.
Parameters for non-bonded interactions.