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