Rename some paramters in nblib listed forces
[alexxy/gromacs.git] / api / nblib / listed_forces / bondtypes.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2020,2021, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \inpublicapi \file
36  * \brief
37  * Implements nblib supported bondtypes
38  *
39  * We choose to forward comparison operations to the
40  * corresponding std::tuple comparison operations.
41  * In order to do that without temporary copies,
42  * we employ std::tie, which requires lvalues as input.
43  * For this reason, bond type parameter getters are implemented
44  * with a const lvalue reference return.
45  *
46  * \author Victor Holanda <victor.holanda@cscs.ch>
47  * \author Joe Jordan <ejjordan@kth.se>
48  * \author Prashanth Kanduri <kanduri@cscs.ch>
49  * \author Sebastian Keller <keller@cscs.ch>
50  * \author Artem Zhmurov <zhmurov@gmail.com>
51  */
52 #ifndef NBLIB_LISTEDFORCES_BONDTYPES_H
53 #define NBLIB_LISTEDFORCES_BONDTYPES_H
54
55 #include <array>
56
57 #include "nblib/particletype.h"
58 #include "nblib/util/util.hpp"
59
60 namespace nblib
61 {
62 using Name          = std::string;
63 using ForceConstant = real;
64 using EquilConstant = real;
65 using Exponent      = real;
66
67 using Degrees = StrongType<real, struct DegreeParameter>;
68 using Radians = StrongType<real, struct RadianParameter>;
69
70 /*! \brief Basic template for interactions with 2 parameters named forceConstant and equilConstant
71  *
72  * \tparam Phantom unused template parameter for type distinction
73  *
74  * Distinct bond types can be generated from this template with using declarations
75  * and declared, but undefined structs. For example:
76  * using HarmonicBondType = TwoParameterInteraction<struct HarmonicBondTypeParameter>;
77  * Note that HarmonicBondTypeParameter does not have to be defined.
78  */
79 template<class Phantom>
80 class TwoParameterInteraction
81 {
82 public:
83     TwoParameterInteraction() = default;
84     TwoParameterInteraction(ForceConstant f, EquilConstant d) : forceConstant_(f), equilConstant_(d)
85     {
86     }
87
88     [[nodiscard]] const ForceConstant& forceConstant() const { return forceConstant_; }
89     [[nodiscard]] const EquilConstant& equilConstant() const { return equilConstant_; }
90
91 private:
92     ForceConstant forceConstant_;
93     EquilConstant equilConstant_;
94 };
95
96 template<class Phantom>
97 inline bool operator<(const TwoParameterInteraction<Phantom>& a, const TwoParameterInteraction<Phantom>& b)
98 {
99     return std::tie(a.forceConstant(), a.equilConstant())
100            < std::tie(b.forceConstant(), b.equilConstant());
101 }
102
103 template<class Phantom>
104 inline bool operator==(const TwoParameterInteraction<Phantom>& a, const TwoParameterInteraction<Phantom>& b)
105 {
106     return std::tie(a.forceConstant(), a.equilConstant())
107            == std::tie(b.forceConstant(), b.equilConstant());
108 }
109
110 /*! \brief harmonic bond type
111  *
112  *  It represents the interaction of the form
113  *  V(r, forceConstant, equilDistance) = 0.5 * forceConstant * (r - equilConstant)^2
114  */
115 using HarmonicBondType = TwoParameterInteraction<struct HarmonicBondTypeParameter>;
116
117
118 /*! \brief GROMOS bond type
119  *
120  * It represents the interaction of the form
121  * V(r, forceConstant, equilDistance) = 0.25 * forceConstant * (r^2 - equilConstant^2)^2
122  */
123 using G96BondType = TwoParameterInteraction<struct G96BondTypeParameter>;
124
125
126 /*! \brief FENE bond type
127  *
128  * It represents the interaction of the form
129  * V(r, forceConstant, equilConstant) = - 0.5 * forceConstant * equilDistance^2 * log( 1 - (r / equilConstant)^2)
130  */
131 using FENEBondType = TwoParameterInteraction<struct FENEBondTypeParameter>;
132
133
134 /*! \brief Half-attractive quartic bond type
135  *
136  * It represents the interaction of the form
137  * V(r, forceConstant, equilConstant) = 0.5 * forceConstant * (r - equilConstant)^4
138  */
139 using HalfAttractiveQuarticBondType =
140         TwoParameterInteraction<struct HalfAttractiveQuarticBondTypeParameter>;
141
142
143 /*! \brief Cubic bond type
144  *
145  * It represents the interaction of the form
146  * V(r, quadraticForceConstant, cubicForceConstant, equilConstant) = quadraticForceConstant * (r -
147  * equilConstant)^2 + quadraticForceConstant * cubicForceConstant * (r - equilDistance)
148  */
149 struct CubicBondType
150 {
151     CubicBondType() = default;
152     CubicBondType(ForceConstant fq, ForceConstant fc, EquilConstant d) :
153         quadraticForceConstant_(fq), cubicForceConstant_(fc), equilDistance_(d)
154     {
155     }
156
157     [[nodiscard]] const ForceConstant& quadraticForceConstant() const
158     {
159         return quadraticForceConstant_;
160     }
161     [[nodiscard]] const ForceConstant& cubicForceConstant() const { return cubicForceConstant_; }
162     [[nodiscard]] const EquilConstant& equilDistance() const { return equilDistance_; }
163
164 private:
165     ForceConstant quadraticForceConstant_;
166     ForceConstant cubicForceConstant_;
167     EquilConstant equilDistance_;
168 };
169
170 inline bool operator<(const CubicBondType& a, const CubicBondType& b)
171 {
172     return std::tie(a.quadraticForceConstant(), a.cubicForceConstant(), a.equilDistance())
173            < std::tie(b.quadraticForceConstant(), b.cubicForceConstant(), b.equilDistance());
174 }
175
176 inline bool operator==(const CubicBondType& a, const CubicBondType& b)
177 {
178     return std::tie(a.quadraticForceConstant(), a.cubicForceConstant(), a.equilDistance())
179            == std::tie(b.quadraticForceConstant(), b.cubicForceConstant(), b.equilDistance());
180 }
181
182 /*! \brief Morse bond type
183  *
184  * It represents the interaction of the form
185  * V(r, forceConstant, exponent, equilDistance) = forceConstant * ( 1 - exp( -exponent * (r - equilConstant))
186  */
187 class MorseBondType
188 {
189 public:
190     MorseBondType() = default;
191     MorseBondType(ForceConstant f, Exponent e, EquilConstant d) :
192         forceConstant_(f), exponent_(e), equilDistance_(d)
193     {
194     }
195
196     [[nodiscard]] const ForceConstant& forceConstant() const { return forceConstant_; }
197     [[nodiscard]] const Exponent&      exponent() const { return exponent_; }
198     [[nodiscard]] const EquilConstant& equilDistance() const { return equilDistance_; }
199
200 private:
201     ForceConstant forceConstant_;
202     Exponent      exponent_;
203     EquilConstant equilDistance_;
204 };
205
206 inline bool operator<(const MorseBondType& a, const MorseBondType& b)
207 {
208     return std::tie(a.forceConstant(), a.exponent(), a.equilDistance())
209            < std::tie(b.forceConstant(), b.exponent(), b.equilDistance());
210 }
211
212 inline bool operator==(const MorseBondType& a, const MorseBondType& b)
213 {
214     return std::tie(a.forceConstant(), a.exponent(), a.equilDistance())
215            == std::tie(b.forceConstant(), b.exponent(), b.equilDistance());
216 }
217
218
219 /*! \brief Harmonic angle type
220  *
221  * Note: the angle is always stored as radians internally
222  */
223 struct HarmonicAngle : public TwoParameterInteraction<struct HarmonicAngleTypeParameter>
224 {
225     HarmonicAngle() = default;
226     //! \brief construct from angle given in radians
227     HarmonicAngle(Radians angle, ForceConstant f) :
228         TwoParameterInteraction<struct HarmonicAngleTypeParameter>{ f, angle }
229     {
230     }
231
232     //! \brief construct from angle given in degrees
233     HarmonicAngle(Degrees angle, ForceConstant f) :
234         TwoParameterInteraction<struct HarmonicAngleTypeParameter>{ f, angle * DEG2RAD }
235     {
236     }
237 };
238
239 /*! \brief Proper Dihedral Implementation
240  */
241 class ProperDihedral
242 {
243 public:
244     using Multiplicity = int;
245
246     ProperDihedral() = default;
247     ProperDihedral(Radians phi, ForceConstant f, Multiplicity m) :
248         phi_(phi), forceConstant_(f), multiplicity_(m)
249     {
250     }
251     ProperDihedral(Degrees phi, ForceConstant f, Multiplicity m) :
252         phi_(phi * DEG2RAD), forceConstant_(f), multiplicity_(m)
253     {
254     }
255
256     [[nodiscard]] const EquilConstant& equilDistance() const { return phi_; }
257     [[nodiscard]] const ForceConstant& forceConstant() const { return forceConstant_; }
258     [[nodiscard]] const Multiplicity&  multiplicity() const { return multiplicity_; }
259
260 private:
261     EquilConstant phi_;
262     ForceConstant forceConstant_;
263     Multiplicity  multiplicity_;
264 };
265
266 inline bool operator<(const ProperDihedral& a, const ProperDihedral& b)
267 {
268     return std::tie(a.equilDistance(), a.forceConstant(), a.multiplicity())
269            < std::tie(b.equilDistance(), b.forceConstant(), b.multiplicity());
270 }
271
272 inline bool operator==(const ProperDihedral& a, const ProperDihedral& b)
273 {
274     return std::tie(a.equilDistance(), a.forceConstant(), a.multiplicity())
275            == std::tie(b.equilDistance(), b.forceConstant(), b.multiplicity());
276 }
277
278
279 /*! \brief Improper Dihedral Implementation
280  */
281 struct ImproperDihedral : public TwoParameterInteraction<struct ImproperDihdedralParameter>
282 {
283     ImproperDihedral() = default;
284     ImproperDihedral(Radians phi, ForceConstant f) :
285         TwoParameterInteraction<struct ImproperDihdedralParameter>{ f, phi }
286     {
287     }
288     ImproperDihedral(Degrees phi, ForceConstant f) :
289         TwoParameterInteraction<struct ImproperDihdedralParameter>{ f, phi * DEG2RAD }
290     {
291     }
292 };
293
294 /*! \brief Ryckaert-Belleman Dihedral Implementation
295  */
296 class RyckaertBellemanDihedral
297 {
298 public:
299     RyckaertBellemanDihedral() = default;
300     RyckaertBellemanDihedral(real p1, real p2, real p3, real p4, real p5, real p6) :
301         parameters_{ p1, p2, p3, p4, p5, p6 }
302     {
303     }
304
305     const real& operator[](std::size_t i) const { return parameters_[i]; }
306
307     [[nodiscard]] const std::array<real, 6>& parameters() const { return parameters_; }
308
309     [[nodiscard]] std::size_t size() const { return parameters_.size(); }
310
311 private:
312     std::array<real, 6> parameters_;
313 };
314
315 inline bool operator<(const RyckaertBellemanDihedral& a, const RyckaertBellemanDihedral& b)
316 {
317     return a.parameters() < b.parameters();
318 }
319
320 inline bool operator==(const RyckaertBellemanDihedral& a, const RyckaertBellemanDihedral& b)
321 {
322     return a.parameters() == b.parameters();
323 }
324
325
326 /*! \brief Type for 5-center interaction (C-MAP)
327  *
328  *  Note: no kernels currently implemented
329  */
330 class Default5Center
331 {
332 public:
333     Default5Center() = default;
334     Default5Center(Radians phi, Radians psi, ForceConstant fphi, ForceConstant fpsi) :
335         phi_(phi), psi_(psi), fphi_(fphi), fpsi_(fpsi)
336     {
337     }
338
339     [[nodiscard]] const Radians&       phi() const { return phi_; }
340     [[nodiscard]] const Radians&       psi() const { return psi_; }
341     [[nodiscard]] const ForceConstant& fphi() const { return fphi_; }
342     [[nodiscard]] const ForceConstant& fpsi() const { return fpsi_; }
343
344 private:
345     Radians       phi_, psi_;
346     ForceConstant fphi_, fpsi_;
347 };
348
349 inline bool operator<(const Default5Center& a, const Default5Center& b)
350 {
351     return std::tie(a.phi(), a.psi(), a.fphi(), a.fpsi())
352            < std::tie(b.phi(), b.psi(), b.fphi(), b.fpsi());
353 }
354
355 inline bool operator==(const Default5Center& a, const Default5Center& b)
356 {
357     return std::tie(a.phi(), a.psi(), a.fphi(), a.fpsi())
358            == std::tie(b.phi(), b.psi(), b.fphi(), b.fpsi());
359 }
360
361
362 } // namespace nblib
363 #endif // NBLIB_LISTEDFORCES_BONDTYPES_H