ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
ljcos2.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_LJCOS2_HPP
22#define CORE_NB_IA_LJCOS2_HPP
23
24/** \file
25 * Routines to calculate the Lennard-Jones with cosine tail potential
26 * between particle pairs. Cosine tail is different from that in
27 * \ref ljcos.hpp. Used for attractive tail/tail interactions in lipid
28 * bilayer calculations.
29 *
30 * Implementation in \ref ljcos2.cpp.
31 */
32
33#include "config/config.hpp"
34
35#ifdef LJCOS2
36
38
39#include <utils/constants.hpp>
41#include <utils/math/sqr.hpp>
42
43#include <cmath>
44
45/** Calculate Lennard-Jones cosine squared force factor */
46inline double ljcos2_pair_force_factor(IA_parameters const &ia_params,
47 double dist) {
48 if (dist < (ia_params.ljcos2.cut + ia_params.ljcos2.offset)) {
49 auto const r_off = dist - ia_params.ljcos2.offset;
50 auto fac = 0.0;
51 if (r_off < ia_params.ljcos2.rchange) {
52 auto const frac6 = Utils::int_pow<6>(ia_params.ljcos2.sig / r_off);
53 fac =
54 48.0 * ia_params.ljcos2.eps * frac6 * (frac6 - 0.5) / (r_off * dist);
55 } else if (r_off < ia_params.ljcos2.rchange + ia_params.ljcos2.w) {
56 fac = -ia_params.ljcos2.eps * Utils::pi() / 2 / ia_params.ljcos2.w /
57 dist *
58 sin(Utils::pi() * (r_off - ia_params.ljcos2.rchange) /
59 ia_params.ljcos2.w);
60 }
61 return fac;
62 }
63 return 0.0;
64}
65
66/** Calculate Lennard-Jones cosine squared energy */
67inline double ljcos2_pair_energy(IA_parameters const &ia_params, double dist) {
68 if (dist < (ia_params.ljcos2.cut + ia_params.ljcos2.offset)) {
69 auto const r_off = dist - ia_params.ljcos2.offset;
70 if (r_off < ia_params.ljcos2.rchange) {
71 auto const frac6 = Utils::int_pow<6>(ia_params.ljcos2.sig / r_off);
72 return 4.0 * ia_params.ljcos2.eps * (Utils::sqr(frac6) - frac6);
73 }
74 if (r_off < (ia_params.ljcos2.rchange + ia_params.ljcos2.w)) {
75 return -ia_params.ljcos2.eps / 2 *
76 (cos(Utils::pi() * (r_off - ia_params.ljcos2.rchange) /
77 ia_params.ljcos2.w) +
78 1);
79 }
80 }
81 return 0.0;
82}
83
84#endif /* ifdef LJCOS2 */
85#endif
This file contains the defaults for ESPResSo.
double ljcos2_pair_energy(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones cosine squared energy.
Definition ljcos2.hpp:67
double ljcos2_pair_force_factor(IA_parameters const &ia_params, double dist)
Calculate Lennard-Jones cosine squared force factor.
Definition ljcos2.hpp:46
DEVICE_QUALIFIER constexpr T pi()
Ratio of diameter and circumference of a circle.
Definition constants.hpp:36
DEVICE_QUALIFIER constexpr T sqr(T x)
Calculates the SQuaRe of x.
Definition sqr.hpp:26
Various procedures concerning interactions between particles.
Parameters for non-bonded interactions.