4a7182c88df2b848c0c5b1ab8ee1bd84ce6f1669
[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 "gromacs/math/vectypes.h"
52 #include "gromacs/topology/ifunc.h"
53 #include "gromacs/utility/basedefinitions.h"
54
55 struct gmx_cmap_t;
56 struct t_fcdata;
57 struct t_mdatom;
58 struct t_nrnb;
59 struct t_pbc;
60
61 /*! \brief Calculate bond-angle. No PBC is taken into account (use mol-shift) */
62 real bond_angle(const rvec          xi,
63                 const rvec          xj,
64                 const rvec          xk,
65                 const struct t_pbc* pbc,
66                 rvec                r_ij,
67                 rvec                r_kj,
68                 real*               costh,
69                 int*                t1,
70                 int*                t2); /* out */
71
72 /*! \brief Calculate dihedral-angle. No PBC is taken into account (use mol-shift) */
73 real dih_angle(const rvec          xi,
74                const rvec          xj,
75                const rvec          xk,
76                const rvec          xl,
77                const struct t_pbc* pbc,
78                rvec                r_ij,
79                rvec                r_kj,
80                rvec                r_kl,
81                rvec                m,
82                rvec                n, /* out */
83                int*                t1,
84                int*                t2,
85                int*                t3);
86
87 /*! \brief Do an update of the forces for dihedral potentials */
88 void do_dih_fup(int                 i,
89                 int                 j,
90                 int                 k,
91                 int                 l,
92                 real                ddphi,
93                 rvec                r_ij,
94                 rvec                r_kj,
95                 rvec                r_kl,
96                 rvec                m,
97                 rvec                n,
98                 rvec4               f[],
99                 rvec                fshift[],
100                 const struct t_pbc* pbc,
101                 const rvec*         x,
102                 int                 t1,
103                 int                 t2,
104                 int                 t3);
105
106 /*! \brief Make a dihedral fall in the range (-pi,pi) */
107 void make_dp_periodic(real* dp);
108
109 /*! \brief Compute CMAP dihedral energies and forces */
110 real cmap_dihs(int                 nbonds,
111                const t_iatom       forceatoms[],
112                const t_iparams     forceparams[],
113                const gmx_cmap_t*   cmap_grid,
114                const rvec          x[],
115                rvec4               f[],
116                rvec                fshift[],
117                const struct t_pbc* pbc,
118                real gmx_unused lambda,
119                real gmx_unused* dvdlambda,
120                const t_mdatoms gmx_unused* md,
121                t_fcdata gmx_unused* fcd,
122                int gmx_unused* global_atom_index);
123
124 /*! \brief For selecting which flavor of bonded kernel is used for simple bonded types */
125 enum class BondedKernelFlavor
126 {
127     ForcesSimdWhenAvailable, //!< Compute only forces, use SIMD when available; should not be used with perturbed parameters
128     ForcesNoSimd,             //!< Compute only forces, do not use SIMD
129     ForcesAndVirialAndEnergy, //!< Compute forces, virial and energy (no SIMD)
130     ForcesAndEnergy,          //!< Compute forces and energy (no SIMD)
131     Count                     //!< The number of flavors
132 };
133
134 /*! \brief Returns whether the energy should be computed */
135 static constexpr inline bool computeEnergy(const BondedKernelFlavor flavor)
136 {
137     return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy
138             || flavor == BondedKernelFlavor::ForcesAndEnergy);
139 }
140
141 /*! \brief Returns whether the virial should be computed */
142 static constexpr inline bool computeVirial(const BondedKernelFlavor flavor)
143 {
144     return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy);
145 }
146
147 /*! \brief Returns whether the energy and/or virial should be computed */
148 static constexpr inline bool computeEnergyOrVirial(const BondedKernelFlavor flavor)
149 {
150     return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy
151             || flavor == BondedKernelFlavor::ForcesAndEnergy);
152 }
153
154 /*! \brief Calculates bonded interactions for simple bonded types
155  *
156  * Exits with an error when the bonded type is not simple
157  * All pointers should be non-null, except for pbc and g which can be nullptr.
158  * \returns the energy or 0 when \p bondedKernelFlavor did not request the energy.
159  */
160 real calculateSimpleBond(int                 ftype,
161                          int                 numForceatoms,
162                          const t_iatom       forceatoms[],
163                          const t_iparams     forceparams[],
164                          const rvec          x[],
165                          rvec4               f[],
166                          rvec                fshift[],
167                          const struct t_pbc* pbc,
168                          real                lambda,
169                          real*               dvdlambda,
170                          const t_mdatoms*    md,
171                          t_fcdata*           fcd,
172                          int gmx_unused*    global_atom_index,
173                          BondedKernelFlavor bondedKernelFlavor);
174
175 //! Getter for finding the flop count for an \c ftype interaction.
176 int nrnbIndex(int ftype);
177
178 #endif