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,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 /*! \libinternal \file
39 * \brief This file contains declarations necessary for low-level
40 * functions for computing energies and forces for bonded interactions.
42 * \author Mark Abraham <mark.j.abraham@gmail.com>
44 * \ingroup module_listed_forces
47 #ifndef GMX_LISTED_FORCES_BONDED_H
48 #define GMX_LISTED_FORCES_BONDED_H
50 #include "gromacs/math/vectypes.h"
51 #include "gromacs/topology/ifunc.h"
52 #include "gromacs/utility/basedefinitions.h"
61 /*! \brief Calculate bond-angle. No PBC is taken into account (use mol-shift) */
62 real bond_angle(const rvec xi,
65 const struct t_pbc* pbc,
72 /*! \brief Calculate dihedral-angle. No PBC is taken into account (use mol-shift) */
73 real dih_angle(const rvec xi,
77 const struct t_pbc* pbc,
87 /*! \brief Do an update of the forces for dihedral potentials */
88 void do_dih_fup(int i,
100 const struct t_pbc* pbc,
101 const struct t_graph* g,
107 /*! \brief Make a dihedral fall in the range (-pi,pi) */
108 void make_dp_periodic(real* dp);
110 /*! \brief Compute CMAP dihedral energies and forces */
111 real cmap_dihs(int nbonds,
112 const t_iatom forceatoms[],
113 const t_iparams forceparams[],
114 const gmx_cmap_t* cmap_grid,
118 const struct t_pbc* pbc,
119 const struct t_graph* g,
120 real gmx_unused lambda,
121 real gmx_unused* dvdlambda,
122 const t_mdatoms gmx_unused* md,
123 t_fcdata gmx_unused* fcd,
124 int gmx_unused* global_atom_index);
126 /*! \brief For selecting which flavor of bonded kernel is used for simple bonded types */
127 enum class BondedKernelFlavor
129 ForcesSimdWhenAvailable, //!< Compute only forces, use SIMD when available; should not be used with perturbed parameters
130 ForcesNoSimd, //!< Compute only forces, do not use SIMD
131 ForcesAndVirialAndEnergy, //!< Compute forces, virial and energy (no SIMD)
132 ForcesAndEnergy, //!< Compute forces and energy (no SIMD)
133 Count //!< The number of flavors
136 /*! \brief Returns whether the energy should be computed */
137 static constexpr inline bool computeEnergy(const BondedKernelFlavor flavor)
139 return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy
140 || flavor == BondedKernelFlavor::ForcesAndEnergy);
143 /*! \brief Returns whether the virial should be computed */
144 static constexpr inline bool computeVirial(const BondedKernelFlavor flavor)
146 return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy);
149 /*! \brief Returns whether the energy and/or virial should be computed */
150 static constexpr inline bool computeEnergyOrVirial(const BondedKernelFlavor flavor)
152 return (flavor == BondedKernelFlavor::ForcesAndVirialAndEnergy
153 || flavor == BondedKernelFlavor::ForcesAndEnergy);
156 /*! \brief Calculates bonded interactions for simple bonded types
158 * Exits with an error when the bonded type is not simple
159 * All pointers should be non-null, except for pbc and g which can be nullptr.
160 * \returns the energy or 0 when \p bondedKernelFlavor did not request the energy.
162 real calculateSimpleBond(int ftype,
164 const t_iatom forceatoms[],
165 const t_iparams forceparams[],
169 const struct t_pbc* pbc,
170 const struct t_graph gmx_unused* g,
175 int gmx_unused* global_atom_index,
176 BondedKernelFlavor bondedKernelFlavor);
178 //! Getter for finding the flop count for an \c ftype interaction.
179 int nrnbIndex(int ftype);