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
43 #include "gromacs/utility/arrayref.h"
44 #include "gromacs/utility/enumerationhelpers.h"
51 //! Force group available for selection for multiple time step integration
52 enum class MtsForceGroups
: int
54 LongrangeNonbonded
, //!< PME-mesh or Ewald for electrostatics and/or LJ
55 Nonbonded
, //!< Non-bonded pair interactions
56 Pair
, //!< Bonded pair interactions
57 Dihedral
, //!< Dihedrals, including cmap (not restraints)
58 Angle
, //!< Bonded angle potentials (not restraints)
59 Pull
, //!< COM pulling
60 Awh
, //!< Accelerated weight histogram method
61 Count
//!< The number of groups above
64 //! Names for the MTS force groups
65 static const gmx::EnumerationArray
<MtsForceGroups
, std::string
> mtsForceGroupNames
= {
66 "longrange-nonbonded", "nonbonded", "pair", "dihedral", "angle", "pull", "awh"
69 //! Setting for a single level for multiple time step integration
72 //! The force group selection for this level;
73 std::bitset
<static_cast<int>(MtsForceGroups::Count
)> forceGroups
;
74 //! The factor between the base, fastest, time step and the time step for this level
78 /*! \brief Returns the MTS level at which a force group is to be computed
80 * \param[in] mtsLevels List of force groups for each MTS level, can be empty without MTS
81 * \param[in] mtsForceGroup The force group to query the MTS level for
83 static inline int forceGroupMtsLevel(ArrayRef
<const MtsLevel
> mtsLevels
, const MtsForceGroups mtsForceGroup
)
85 GMX_ASSERT(mtsLevels
.empty() || mtsLevels
.size() == 2, "Only 0 or 2 MTS levels are supported");
87 return (mtsLevels
.empty() || mtsLevels
[0].forceGroups
[static_cast<int>(mtsForceGroup
)]) ? 0 : 1;
90 /*! \brief Returns the interval in steps at which the non-bonded pair forces are calculated
92 * Note: returns 1 when multiple time-stepping is not activated.
94 int nonbondedMtsFactor(const t_inputrec
& ir
);
96 //! Struct for passing the MTS mdp options to setupMtsLevels()
99 //! The number of MTS levels
101 //! The names of the force groups assigned by the user to level 2, internal index 1
102 std::string level2Forces
;
103 //! The step factor assigned by the user to level 2, internal index 1
104 int level2Factor
= 0;
107 /*! \brief Sets up and returns the MTS levels and checks requirements of MTS
109 * Appends errors about allowed input values ir to errorMessages, when not nullptr.
111 * \param[in] mtsOpts Options for setting the MTS levels
112 * \param[in,out] errorMessages List of error messages, can be nullptr
114 std::vector
<MtsLevel
> setupMtsLevels(const GromppMtsOpts
& mtsOpts
, std::vector
<std::string
>* errorMessages
);
116 /*! \brief Returns whether we use MTS and the MTS setup is internally valid
118 * Note that setupMtsLevels would have returned at least one error message
119 * when this function returns false
121 bool haveValidMtsSetup(const t_inputrec
& ir
);
123 /*! \brief Checks whether the MTS requirements on other algorithms and output frequencies are met
125 * Note: exits with an assertion failure when
126 * ir.useMts == true && haveValidMtsSetup(ir) == false
128 * \param[in] ir Complete input record
129 * \returns list of error messages, empty when all MTS requirements are met
131 std::vector
<std::string
> checkMtsRequirements(const t_inputrec
& ir
);
135 #endif /* GMX_MULTIPLETIMESTEPPING_H */