Simplify nblib listed forces type tests
[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 class G96BondType : public TwoParameterInteraction<struct G96BondTypeParameter>
124 {
125 public:
126     G96BondType() = default;
127     //! \brief Store square of equilibrium distance
128     G96BondType(ForceConstant f, EquilConstant equilConstant) :
129         TwoParameterInteraction<struct G96BondTypeParameter>{ f, equilConstant * equilConstant }
130     {
131     }
132 };
133
134
135 /*! \brief FENE bond type
136  *
137  * It represents the interaction of the form
138  * V(r, forceConstant, equilConstant) = - 0.5 * forceConstant * equilDistance^2 * log( 1 - (r / equilConstant)^2)
139  */
140 using FENEBondType = TwoParameterInteraction<struct FENEBondTypeParameter>;
141
142
143 /*! \brief Half-attractive quartic bond type
144  *
145  * It represents the interaction of the form
146  * V(r, forceConstant, equilConstant) = 0.5 * forceConstant * (r - equilConstant)^4
147  */
148 using HalfAttractiveQuarticBondType =
149         TwoParameterInteraction<struct HalfAttractiveQuarticBondTypeParameter>;
150
151
152 /*! \brief Cubic bond type
153  *
154  * It represents the interaction of the form
155  * V(r, quadraticForceConstant, cubicForceConstant, equilConstant) = quadraticForceConstant * (r -
156  * equilConstant)^2 + quadraticForceConstant * cubicForceConstant * (r - equilDistance)
157  */
158 struct CubicBondType
159 {
160     CubicBondType() = default;
161     CubicBondType(ForceConstant fq, ForceConstant fc, EquilConstant d) :
162         quadraticForceConstant_(fq), cubicForceConstant_(fc), equilDistance_(d)
163     {
164     }
165
166     [[nodiscard]] const ForceConstant& quadraticForceConstant() const
167     {
168         return quadraticForceConstant_;
169     }
170     [[nodiscard]] const ForceConstant& cubicForceConstant() const { return cubicForceConstant_; }
171     [[nodiscard]] const EquilConstant& equilDistance() const { return equilDistance_; }
172
173 private:
174     ForceConstant quadraticForceConstant_;
175     ForceConstant cubicForceConstant_;
176     EquilConstant equilDistance_;
177 };
178
179 inline bool operator<(const CubicBondType& a, const CubicBondType& b)
180 {
181     return std::tie(a.quadraticForceConstant(), a.cubicForceConstant(), a.equilDistance())
182            < std::tie(b.quadraticForceConstant(), b.cubicForceConstant(), b.equilDistance());
183 }
184
185 inline bool operator==(const CubicBondType& a, const CubicBondType& b)
186 {
187     return std::tie(a.quadraticForceConstant(), a.cubicForceConstant(), a.equilDistance())
188            == std::tie(b.quadraticForceConstant(), b.cubicForceConstant(), b.equilDistance());
189 }
190
191 /*! \brief Morse bond type
192  *
193  * It represents the interaction of the form
194  * V(r, forceConstant, exponent, equilDistance) = forceConstant * ( 1 - exp( -exponent * (r - equilConstant))
195  */
196 class MorseBondType
197 {
198 public:
199     MorseBondType() = default;
200     MorseBondType(ForceConstant f, Exponent e, EquilConstant d) :
201         forceConstant_(f), exponent_(e), equilDistance_(d)
202     {
203     }
204
205     [[nodiscard]] const ForceConstant& forceConstant() const { return forceConstant_; }
206     [[nodiscard]] const Exponent&      exponent() const { return exponent_; }
207     [[nodiscard]] const EquilConstant& equilDistance() const { return equilDistance_; }
208
209 private:
210     ForceConstant forceConstant_;
211     Exponent      exponent_;
212     EquilConstant equilDistance_;
213 };
214
215 inline bool operator<(const MorseBondType& a, const MorseBondType& b)
216 {
217     return std::tie(a.forceConstant(), a.exponent(), a.equilDistance())
218            < std::tie(b.forceConstant(), b.exponent(), b.equilDistance());
219 }
220
221 inline bool operator==(const MorseBondType& a, const MorseBondType& b)
222 {
223     return std::tie(a.forceConstant(), a.exponent(), a.equilDistance())
224            == std::tie(b.forceConstant(), b.exponent(), b.equilDistance());
225 }
226
227
228 /*! \brief Basic template for interactions with 2 parameters named forceConstant and equilAngle
229  *
230  * \tparam Phantom unused template parameter for type distinction
231  *
232  * Distinct angle types can be generated from this template with using declarations
233  * and declared, but undefined structs. For example:
234  * using HarmonicAngleType = AngleInteractionType<struct HarmonicAngleParameter>;
235  * HarmonicAngleParameter does not have to be defined.
236  *
237  * Note: the angle is always stored as radians internally
238  */
239 template<class Phantom>
240 class AngleInteractionType : public TwoParameterInteraction<Phantom>
241 {
242 public:
243     AngleInteractionType() = default;
244     //! \brief construct from angle given in radians
245     AngleInteractionType(ForceConstant f, Radians angle) :
246         TwoParameterInteraction<Phantom>{ f, angle }
247     {
248     }
249
250     //! \brief construct from angle given in degrees
251     AngleInteractionType(ForceConstant f, Degrees angle) :
252         TwoParameterInteraction<Phantom>{ f, angle * DEG2RAD }
253     {
254     }
255 };
256
257 /*! \brief Harmonic angle type
258  *
259  * It represents the interaction of the form
260  * V(theta, forceConstant, equilAngle) = 0.5 * forceConstant * (theta - equilAngle)^2
261  */
262 using HarmonicAngle = AngleInteractionType<struct HarmonicAngleParameter>;
263
264 /*! \brief Proper Dihedral Implementation
265  */
266 class ProperDihedral
267 {
268 public:
269     using Multiplicity = int;
270
271     ProperDihedral() = default;
272     ProperDihedral(Radians phi, ForceConstant f, Multiplicity m) :
273         phi_(phi), forceConstant_(f), multiplicity_(m)
274     {
275     }
276     ProperDihedral(Degrees phi, ForceConstant f, Multiplicity m) :
277         phi_(phi * DEG2RAD), forceConstant_(f), multiplicity_(m)
278     {
279     }
280
281     [[nodiscard]] const EquilConstant& equilDistance() const { return phi_; }
282     [[nodiscard]] const ForceConstant& forceConstant() const { return forceConstant_; }
283     [[nodiscard]] const Multiplicity&  multiplicity() const { return multiplicity_; }
284
285 private:
286     EquilConstant phi_;
287     ForceConstant forceConstant_;
288     Multiplicity  multiplicity_;
289 };
290
291 inline bool operator<(const ProperDihedral& a, const ProperDihedral& b)
292 {
293     return std::tie(a.equilDistance(), a.forceConstant(), a.multiplicity())
294            < std::tie(b.equilDistance(), b.forceConstant(), b.multiplicity());
295 }
296
297 inline bool operator==(const ProperDihedral& a, const ProperDihedral& b)
298 {
299     return std::tie(a.equilDistance(), a.forceConstant(), a.multiplicity())
300            == std::tie(b.equilDistance(), b.forceConstant(), b.multiplicity());
301 }
302
303
304 /*! \brief Improper Dihedral Implementation
305  */
306 struct ImproperDihedral : public TwoParameterInteraction<struct ImproperDihdedralParameter>
307 {
308     ImproperDihedral() = default;
309     ImproperDihedral(Radians phi, ForceConstant f) :
310         TwoParameterInteraction<struct ImproperDihdedralParameter>{ f, phi }
311     {
312     }
313     ImproperDihedral(Degrees phi, ForceConstant f) :
314         TwoParameterInteraction<struct ImproperDihdedralParameter>{ f, phi * DEG2RAD }
315     {
316     }
317 };
318
319 /*! \brief Ryckaert-Belleman Dihedral Implementation
320  */
321 class RyckaertBellemanDihedral
322 {
323 public:
324     RyckaertBellemanDihedral() = default;
325     RyckaertBellemanDihedral(real p1, real p2, real p3, real p4, real p5, real p6) :
326         parameters_{ p1, p2, p3, p4, p5, p6 }
327     {
328     }
329
330     const real& operator[](std::size_t i) const { return parameters_[i]; }
331
332     [[nodiscard]] const std::array<real, 6>& parameters() const { return parameters_; }
333
334     [[nodiscard]] std::size_t size() const { return parameters_.size(); }
335
336 private:
337     std::array<real, 6> parameters_;
338 };
339
340 inline bool operator<(const RyckaertBellemanDihedral& a, const RyckaertBellemanDihedral& b)
341 {
342     return a.parameters() < b.parameters();
343 }
344
345 inline bool operator==(const RyckaertBellemanDihedral& a, const RyckaertBellemanDihedral& b)
346 {
347     return a.parameters() == b.parameters();
348 }
349
350
351 /*! \brief Type for 5-center interaction (C-MAP)
352  *
353  *  Note: no kernels currently implemented
354  */
355 class Default5Center
356 {
357 public:
358     Default5Center() = default;
359     Default5Center(Radians phi, Radians psi, ForceConstant fphi, ForceConstant fpsi) :
360         phi_(phi), psi_(psi), fphi_(fphi), fpsi_(fpsi)
361     {
362     }
363
364     [[nodiscard]] const Radians&       phi() const { return phi_; }
365     [[nodiscard]] const Radians&       psi() const { return psi_; }
366     [[nodiscard]] const ForceConstant& fphi() const { return fphi_; }
367     [[nodiscard]] const ForceConstant& fpsi() const { return fpsi_; }
368
369 private:
370     Radians       phi_, psi_;
371     ForceConstant fphi_, fpsi_;
372 };
373
374 inline bool operator<(const Default5Center& a, const Default5Center& b)
375 {
376     return std::tie(a.phi(), a.psi(), a.fphi(), a.fpsi())
377            < std::tie(b.phi(), b.psi(), b.fphi(), b.fpsi());
378 }
379
380 inline bool operator==(const Default5Center& a, const Default5Center& b)
381 {
382     return std::tie(a.phi(), a.psi(), a.fphi(), a.fpsi())
383            == std::tie(b.phi(), b.psi(), b.fphi(), b.fpsi());
384 }
385
386
387 } // namespace nblib
388 #endif // NBLIB_LISTEDFORCES_BONDTYPES_H