2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
38 /*! \libinternal \file
40 * \brief This file contains declarations necessary for low-level
41 * functions for computing energies and forces for bonded interactions.
43 * \author Mark Abraham <mark.j.abraham@gmail.com>
45 * \ingroup module_listed_forces
48 #ifndef GMX_LISTED_FORCES_BONDED_H
49 #define GMX_LISTED_FORCES_BONDED_H
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/topology/ifunc.h"
53 #include "gromacs/utility/basedefinitions.h"
62 /*! \brief Calculate bond-angle. No PBC is taken into account (use mol-shift) */
63 real bond_angle(const rvec xi,
66 const struct t_pbc* pbc,
73 /*! \brief Calculate dihedral-angle. No PBC is taken into account (use mol-shift) */
74 real dih_angle(const rvec xi,
78 const struct t_pbc* pbc,
88 /*! \brief Do an update of the forces for dihedral potentials */
89 void do_dih_fup(int i,
101 const struct t_pbc* pbc,
102 const struct t_graph* g,
108 /*! \brief Make a dihedral fall in the range (-pi,pi) */
109 void make_dp_periodic(real* dp);
111 /*! \brief Compute CMAP dihedral energies and forces */
112 real cmap_dihs(int nbonds,
113 const t_iatom forceatoms[],
114 const t_iparams forceparams[],
115 const gmx_cmap_t* cmap_grid,
119 const struct t_pbc* pbc,
120 const struct t_graph* g,
121 real gmx_unused lambda,
122 real gmx_unused* dvdlambda,
123 const t_mdatoms gmx_unused* md,
124 t_fcdata gmx_unused* fcd,
125 int gmx_unused* global_atom_index);
127 /*! \brief For selecting which flavor of bonded kernel is used for simple bonded types */
128 enum class BondedKernelFlavor
130 ForcesSimdWhenAvailable, //!< Compute only forces, use SIMD when available; should not be used with perturbed parameters
131 ForcesNoSimd, //!< Compute only forces, do not use SIMD
132 ForcesAndVirialAndEnergy, //!< Compute forces, virial and energy (no SIMD)
133 ForcesAndEnergy, //!< Compute forces and energy (no SIMD)
134 Count //!< The number of flavors
137 /*! \brief Returns whether the energy should be computed */
138 static constexpr inline bool computeEnergy(const BondedKernelFlavor flavor)
140 return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy
141 || flavor == BondedKernelFlavor::ForcesAndEnergy);
144 /*! \brief Returns whether the virial should be computed */
145 static constexpr inline bool computeVirial(const BondedKernelFlavor flavor)
147 return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy);
150 /*! \brief Returns whether the energy and/or virial should be computed */
151 static constexpr inline bool computeEnergyOrVirial(const BondedKernelFlavor flavor)
153 return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy
154 || flavor == BondedKernelFlavor::ForcesAndEnergy);
157 /*! \brief Calculates bonded interactions for simple bonded types
159 * Exits with an error when the bonded type is not simple
160 * All pointers should be non-null, except for pbc and g which can be nullptr.
161 * \returns the energy or 0 when \p bondedKernelFlavor did not request the energy.
163 real calculateSimpleBond(int ftype,
165 const t_iatom forceatoms[],
166 const t_iparams forceparams[],
170 const struct t_pbc* pbc,
171 const struct t_graph gmx_unused* g,
176 int gmx_unused* global_atom_index,
177 BondedKernelFlavor bondedKernelFlavor);
179 //! Getter for finding the flop count for an \c ftype interaction.
180 int nrnbIndex(int ftype);