beb95b55289bc9d108008e7985c4b7ff42ac406a
[alexxy/gromacs.git] / src / gromacs / listed_forces / bonded.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /*! \libinternal \file
39  *
40  * \brief This file contains declarations necessary for low-level
41  * functions for computing energies and forces for bonded interactions.
42  *
43  * \author Mark Abraham <mark.j.abraham@gmail.com>
44  * \inlibraryapi
45  * \ingroup module_listed_forces
46  */
47
48 #ifndef GMX_LISTED_FORCES_BONDED_H
49 #define GMX_LISTED_FORCES_BONDED_H
50
51 #include <string>
52
53 #include "gromacs/math/vectypes.h"
54 #include "gromacs/topology/ifunc.h"
55 #include "gromacs/utility/basedefinitions.h"
56
57 struct gmx_cmap_t;
58 struct t_fcdata;
59 struct t_mdatom;
60 struct t_nrnb;
61 struct t_pbc;
62
63 namespace gmx
64 {
65 template<typename EnumType, typename DataType, EnumType ArraySize>
66 struct EnumerationArray;
67 } // namespace gmx
68
69 /*! \brief Calculate bond-angle. No PBC is taken into account (use mol-shift) */
70 real bond_angle(const rvec          xi,
71                 const rvec          xj,
72                 const rvec          xk,
73                 const struct t_pbc* pbc,
74                 rvec                r_ij,
75                 rvec                r_kj,
76                 real*               costh,
77                 int*                t1,
78                 int*                t2); /* out */
79
80 /*! \brief Calculate dihedral-angle. No PBC is taken into account (use mol-shift) */
81 real dih_angle(const rvec          xi,
82                const rvec          xj,
83                const rvec          xk,
84                const rvec          xl,
85                const struct t_pbc* pbc,
86                rvec                r_ij,
87                rvec                r_kj,
88                rvec                r_kl,
89                rvec                m,
90                rvec                n, /* out */
91                int*                t1,
92                int*                t2,
93                int*                t3);
94
95 /*! \brief Do an update of the forces for dihedral potentials */
96 void do_dih_fup(int                 i,
97                 int                 j,
98                 int                 k,
99                 int                 l,
100                 real                ddphi,
101                 rvec                r_ij,
102                 rvec                r_kj,
103                 rvec                r_kl,
104                 rvec                m,
105                 rvec                n,
106                 rvec4               f[],
107                 rvec                fshift[],
108                 const struct t_pbc* pbc,
109                 const rvec*         x,
110                 int                 t1,
111                 int                 t2,
112                 int                 t3);
113
114 /*! \brief Make a dihedral fall in the range (-pi,pi) */
115 void make_dp_periodic(real* dp);
116
117 /*! \brief Compute CMAP dihedral energies and forces */
118 real cmap_dihs(int                 nbonds,
119                const t_iatom       forceatoms[],
120                const t_iparams     forceparams[],
121                const gmx_cmap_t*   cmap_grid,
122                const rvec          x[],
123                rvec4               f[],
124                rvec                fshift[],
125                const struct t_pbc* pbc,
126                real gmx_unused lambda,
127                real gmx_unused* dvdlambda,
128                const t_mdatoms gmx_unused* md,
129                t_fcdata gmx_unused* fcd,
130                int gmx_unused* global_atom_index);
131
132 /*! \brief For selecting which flavor of bonded kernel is used for simple bonded types */
133 enum class BondedKernelFlavor
134 {
135     ForcesSimdWhenAvailable, //!< Compute only forces, use SIMD when available; should not be used with perturbed parameters
136     ForcesNoSimd,             //!< Compute only forces, do not use SIMD
137     ForcesAndVirialAndEnergy, //!< Compute forces, virial and energy (no SIMD)
138     ForcesAndEnergy,          //!< Compute forces and energy (no SIMD)
139     Count                     //!< The number of flavors
140 };
141
142 //! Helper strings for human-readable messages
143 extern const gmx::EnumerationArray<BondedKernelFlavor, std::string, BondedKernelFlavor::Count> c_bondedKernelFlavorStrings;
144
145 /*! \brief Returns whether the energy should be computed */
146 static constexpr inline bool computeEnergy(const BondedKernelFlavor flavor)
147 {
148     return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy
149             || flavor == BondedKernelFlavor::ForcesAndEnergy);
150 }
151
152 /*! \brief Returns whether the virial should be computed */
153 static constexpr inline bool computeVirial(const BondedKernelFlavor flavor)
154 {
155     return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy);
156 }
157
158 /*! \brief Returns whether the energy and/or virial should be computed */
159 static constexpr inline bool computeEnergyOrVirial(const BondedKernelFlavor flavor)
160 {
161     return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy
162             || flavor == BondedKernelFlavor::ForcesAndEnergy);
163 }
164
165 /*! \brief Calculates bonded interactions for simple bonded types
166  *
167  * Exits with an error when the bonded type is not simple
168  * All pointers should be non-null, except for pbc and g which can be nullptr.
169  * \returns the energy or 0 when \p bondedKernelFlavor did not request the energy.
170  */
171 real calculateSimpleBond(int                 ftype,
172                          int                 numForceatoms,
173                          const t_iatom       forceatoms[],
174                          const t_iparams     forceparams[],
175                          const rvec          x[],
176                          rvec4               f[],
177                          rvec                fshift[],
178                          const struct t_pbc* pbc,
179                          real                lambda,
180                          real*               dvdlambda,
181                          const t_mdatoms*    md,
182                          t_fcdata*           fcd,
183                          int gmx_unused*    global_atom_index,
184                          BondedKernelFlavor bondedKernelFlavor);
185
186 //! Getter for finding the flop count for an \c ftype interaction.
187 int nrnbIndex(int ftype);
188
189 #endif