2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2020, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 #ifndef GMX_MULTIPLETIMESTEPPING_H
36 #define GMX_MULTIPLETIMESTEPPING_H
40 #include "gromacs/utility/arrayref.h"
41 #include "gromacs/utility/enumerationhelpers.h"
48 //! Force group available for selection for multiple time step integration
49 enum class MtsForceGroups : int
51 LongrangeNonbonded, //!< PME-mesh or Ewald for electrostatics and/or LJ
52 Nonbonded, //!< Non-bonded pair interactions
53 Pair, //!< Bonded pair interactions
54 Dihedral, //!< Dihedrals, including cmap (not restraints)
55 Angle, //!< Bonded angle potentials (not restraints)
56 Pull, //!< COM pulling
57 Awh, //!< Accelerated weight histogram method
58 Count //!< The number of groups above
61 static const gmx::EnumerationArray<MtsForceGroups, std::string> mtsForceGroupNames = {
62 "longrange-nonbonded", "nonbonded", "pair", "dihedral", "angle", "pull", "awh"
65 //! Setting for a single level for multiple time step integration
68 //! The force group selection for this level;
69 std::bitset<static_cast<int>(MtsForceGroups::Count)> forceGroups;
70 //! The factor between the base, fastest, time step and the time step for this level
74 /*! \brief Returns the MTS level at which a force group is to be computed
76 * \param[in] mtsLevels List of force groups for each MTS level, can be empty without MTS
77 * \param[in] mtsForceGroup The force group to query the MTS level for
79 static inline int forceGroupMtsLevel(ArrayRef<const MtsLevel> mtsLevels, const MtsForceGroups mtsForceGroup)
81 GMX_ASSERT(mtsLevels.empty() || mtsLevels.size() == 2, "Only 0 or 2 MTS levels are supported");
83 return (mtsLevels.empty() || mtsLevels[0].forceGroups[static_cast<int>(mtsForceGroup)]) ? 0 : 1;
86 /*! \brief Returns the interval in steps at which the non-bonded pair forces are calculated
88 * Note: returns 1 when multiple time-stepping is not activated.
90 int nonbondedMtsFactor(const t_inputrec& ir);
92 //! (Release) Asserts that all multiple time-stepping requirements on \p ir are fulfilled
93 void assertMtsRequirements(const t_inputrec& ir);
97 #endif /* GMX_MULTIPLETIMESTEPPING_H */