73bdfe2ee96a2d0bbaee51217727247014fbe979
[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 EquilDistance = 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 equilDistance
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, EquilDistance d) : forceConstant_(f), equilDistance_(d)
85     {
86     }
87
88     [[nodiscard]] const ForceConstant& forceConstant() const { return forceConstant_; }
89     [[nodiscard]] const EquilDistance& equilDistance() const { return equilDistance_; }
90
91 private:
92     ForceConstant forceConstant_;
93     EquilDistance equilDistance_;
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.equilDistance())
100            < std::tie(b.forceConstant(), b.equilDistance());
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.equilDistance())
107            == std::tie(b.forceConstant(), b.equilDistance());
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 - equilDistance)^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 - equilDistance^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, equilDistance) = - 0.5 * forceConstant * equilDistance^2 * log( 1 - (r / equilDistance)^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, equilDistance) = 0.5 * forceConstant * (r - equilDistance)^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, equilDistance) = quadraticForceConstant * (r -
147  * equilDistance)^2 + quadraticForceConstant * cubicForceConstant * (r - equilDistance)
148  */
149 struct CubicBondType
150 {
151     CubicBondType() = default;
152     CubicBondType(ForceConstant fq, ForceConstant fc, EquilDistance d) :
153         quadraticForceConstant_(fq),
154         cubicForceConstant_(fc),
155         equilDistance_(d)
156     {
157     }
158
159     [[nodiscard]] const ForceConstant& quadraticForceConstant() const
160     {
161         return quadraticForceConstant_;
162     }
163     [[nodiscard]] const ForceConstant& cubicForceConstant() const { return cubicForceConstant_; }
164     [[nodiscard]] const EquilDistance& equilDistance() const { return equilDistance_; }
165
166 private:
167     ForceConstant quadraticForceConstant_;
168     ForceConstant cubicForceConstant_;
169     EquilDistance equilDistance_;
170 };
171
172 inline bool operator<(const CubicBondType& a, const CubicBondType& b)
173 {
174     return std::tie(a.quadraticForceConstant(), a.cubicForceConstant(), a.equilDistance())
175            < std::tie(b.quadraticForceConstant(), b.cubicForceConstant(), b.equilDistance());
176 }
177
178 inline bool operator==(const CubicBondType& a, const CubicBondType& b)
179 {
180     return std::tie(a.quadraticForceConstant(), a.cubicForceConstant(), a.equilDistance())
181            == std::tie(b.quadraticForceConstant(), b.cubicForceConstant(), b.equilDistance());
182 }
183
184 /*! \brief Morse bond type
185  *
186  * It represents the interaction of the form
187  * V(r; forceConstant, exponent, equilDistance) = forceConstant * ( 1 - exp( -exponent * (r - equilDistance))
188  */
189 class MorseBondType
190 {
191 public:
192     MorseBondType() = default;
193     MorseBondType(ForceConstant f, Exponent e, EquilDistance d) :
194         forceConstant_(f),
195         exponent_(e),
196         equilDistance_(d)
197     {
198     }
199
200     [[nodiscard]] const ForceConstant& forceConstant() const { return forceConstant_; }
201     [[nodiscard]] const Exponent&      exponent() const { return exponent_; }
202     [[nodiscard]] const EquilDistance& equilDistance() const { return equilDistance_; }
203
204 private:
205     ForceConstant forceConstant_;
206     Exponent      exponent_;
207     EquilDistance equilDistance_;
208 };
209
210 inline bool operator<(const MorseBondType& a, const MorseBondType& b)
211 {
212     return std::tie(a.forceConstant(), a.exponent(), a.equilDistance())
213            < std::tie(b.forceConstant(), b.exponent(), b.equilDistance());
214 }
215
216 inline bool operator==(const MorseBondType& a, const MorseBondType& b)
217 {
218     return std::tie(a.forceConstant(), a.exponent(), a.equilDistance())
219            == std::tie(b.forceConstant(), b.exponent(), b.equilDistance());
220 }
221
222
223 /*! \brief Harmonic angle type
224  *
225  * Note: the angle is always stored as radians internally
226  */
227 struct HarmonicAngleType : public TwoParameterInteraction<struct HarmonicAngleTypeParameter>
228 {
229     HarmonicAngleType() = default;
230     //! \brief construct from angle given in radians
231     HarmonicAngleType(Radians angle, ForceConstant f) :
232         TwoParameterInteraction<struct HarmonicAngleTypeParameter>{ f, angle }
233     {
234     }
235
236     //! \brief construct from angle given in degrees
237     HarmonicAngleType(Degrees angle, ForceConstant f) :
238         TwoParameterInteraction<struct HarmonicAngleTypeParameter>{ f, angle * DEG2RAD }
239     {
240     }
241 };
242
243 /*! \brief Proper Dihedral Implementation
244  */
245 class ProperDihedral
246 {
247 public:
248     using Multiplicity = int;
249
250     ProperDihedral() = default;
251     ProperDihedral(Radians phi, ForceConstant f, Multiplicity m) :
252         phi_(phi),
253         forceConstant_(f),
254         multiplicity_(m)
255     {
256     }
257     ProperDihedral(Degrees phi, ForceConstant f, Multiplicity m) :
258         phi_(phi * DEG2RAD),
259         forceConstant_(f),
260         multiplicity_(m)
261     {
262     }
263
264     [[nodiscard]] const EquilDistance& equilDistance() const { return phi_; }
265     [[nodiscard]] const ForceConstant& forceConstant() const { return forceConstant_; }
266     [[nodiscard]] const Multiplicity&  multiplicity() const { return multiplicity_; }
267
268 private:
269     EquilDistance phi_;
270     ForceConstant forceConstant_;
271     Multiplicity  multiplicity_;
272 };
273
274 inline bool operator<(const ProperDihedral& a, const ProperDihedral& b)
275 {
276     return std::tie(a.equilDistance(), a.forceConstant(), a.multiplicity())
277            < std::tie(b.equilDistance(), b.forceConstant(), b.multiplicity());
278 }
279
280 inline bool operator==(const ProperDihedral& a, const ProperDihedral& b)
281 {
282     return std::tie(a.equilDistance(), a.forceConstant(), a.multiplicity())
283            == std::tie(b.equilDistance(), b.forceConstant(), b.multiplicity());
284 }
285
286
287 /*! \brief Improper Dihedral Implementation
288  */
289 struct ImproperDihedral : public TwoParameterInteraction<struct ImproperDihdedralParameter>
290 {
291     ImproperDihedral() = default;
292     ImproperDihedral(Radians phi, ForceConstant f) :
293         TwoParameterInteraction<struct ImproperDihdedralParameter>{ f, phi }
294     {
295     }
296     ImproperDihedral(Degrees phi, ForceConstant f) :
297         TwoParameterInteraction<struct ImproperDihdedralParameter>{ f, phi * DEG2RAD }
298     {
299     }
300 };
301
302 /*! \brief Ryckaert-Belleman Dihedral Implementation
303  */
304 class RyckaertBellemanDihedral
305 {
306 public:
307     RyckaertBellemanDihedral() = default;
308     RyckaertBellemanDihedral(real p1, real p2, real p3, real p4, real p5, real p6) :
309         parameters_{ p1, p2, p3, p4, p5, p6 }
310     {
311     }
312
313     const real& operator[](std::size_t i) const { return parameters_[i]; }
314
315     [[nodiscard]] const std::array<real, 6>& parameters() const { return parameters_; }
316
317     [[nodiscard]] std::size_t size() const { return parameters_.size(); }
318
319 private:
320     std::array<real, 6> parameters_;
321 };
322
323 inline bool operator<(const RyckaertBellemanDihedral& a, const RyckaertBellemanDihedral& b)
324 {
325     return a.parameters() < b.parameters();
326 }
327
328 inline bool operator==(const RyckaertBellemanDihedral& a, const RyckaertBellemanDihedral& b)
329 {
330     return a.parameters() == b.parameters();
331 }
332
333
334 /*! \brief Type for 5-center interaction (C-MAP)
335  *
336  *  Note: no kernels currently implemented
337  */
338 class Default5Center
339 {
340 public:
341     Default5Center() = default;
342     Default5Center(Radians phi, Radians psi, ForceConstant fphi, ForceConstant fpsi) :
343         phi_(phi),
344         psi_(psi),
345         fphi_(fphi),
346         fpsi_(fpsi)
347     {
348     }
349
350     [[nodiscard]] const Radians&       phi() const { return phi_; }
351     [[nodiscard]] const Radians&       psi() const { return psi_; }
352     [[nodiscard]] const ForceConstant& fphi() const { return fphi_; }
353     [[nodiscard]] const ForceConstant& fpsi() const { return fpsi_; }
354
355 private:
356     Radians       phi_, psi_;
357     ForceConstant fphi_, fpsi_;
358 };
359
360 inline bool operator<(const Default5Center& a, const Default5Center& b)
361 {
362     return std::tie(a.phi(), a.psi(), a.fphi(), a.fpsi())
363            < std::tie(b.phi(), b.psi(), b.fphi(), b.fpsi());
364 }
365
366 inline bool operator==(const Default5Center& a, const Default5Center& b)
367 {
368     return std::tie(a.phi(), a.psi(), a.fphi(), a.fpsi())
369            == std::tie(b.phi(), b.psi(), b.fphi(), b.fpsi());
370 }
371
372
373 } // namespace nblib
374 #endif // NBLIB_LISTEDFORCES_BONDTYPES_H