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.
37 #include "multipletimestepping.h"
41 #include "gromacs/mdtypes/inputrec.h"
42 #include "gromacs/mdtypes/pull_params.h"
43 #include "gromacs/utility/arrayref.h"
44 #include "gromacs/utility/gmxassert.h"
45 #include "gromacs/utility/stringutil.h"
50 int nonbondedMtsFactor(const t_inputrec
& ir
)
52 GMX_RELEASE_ASSERT(!ir
.useMts
|| ir
.mtsLevels
.size() == 2, "Only 2 MTS levels supported here");
54 if (ir
.useMts
&& ir
.mtsLevels
[1].forceGroups
[static_cast<int>(MtsForceGroups::Nonbonded
)])
56 return ir
.mtsLevels
[1].stepFactor
;
64 std::vector
<MtsLevel
> setupMtsLevels(const GromppMtsOpts
& mtsOpts
, std::vector
<std::string
>* errorMessages
)
66 std::vector
<MtsLevel
> mtsLevels
;
68 if (mtsOpts
.numLevels
!= 2)
72 errorMessages
->push_back("Only mts-levels = 2 is supported");
79 const std::vector
<std::string
> inputForceGroups
= gmx::splitString(mtsOpts
.level2Forces
);
80 auto& forceGroups
= mtsLevels
[1].forceGroups
;
81 for (const auto& inputForceGroup
: inputForceGroups
)
85 for (const auto& forceGroupName
: gmx::mtsForceGroupNames
)
87 if (gmx::equalCaseInsensitive(inputForceGroup
, forceGroupName
))
89 forceGroups
.set(nameIndex
);
94 if (!found
&& errorMessages
)
96 errorMessages
->push_back(
97 gmx::formatString("Unknown MTS force group '%s'", inputForceGroup
.c_str()));
101 // Make the level 0 use the complement of the force groups of group 1
102 mtsLevels
[0].forceGroups
= ~mtsLevels
[1].forceGroups
;
103 mtsLevels
[0].stepFactor
= 1;
105 mtsLevels
[1].stepFactor
= mtsOpts
.level2Factor
;
107 if (errorMessages
&& mtsLevels
[1].stepFactor
<= 1)
109 errorMessages
->push_back("mts-factor should be larger than 1");
116 bool haveValidMtsSetup(const t_inputrec
& ir
)
118 return (ir
.useMts
&& ir
.mtsLevels
.size() == 2 && ir
.mtsLevels
[1].stepFactor
> 1);
124 //! Checks whether \p nstValue is a multiple of the largest MTS step, returns an error string for parameter \p param when this is not the case
125 std::optional
<std::string
> checkMtsInterval(ArrayRef
<const MtsLevel
> mtsLevels
, const char* param
, const int nstValue
)
127 GMX_RELEASE_ASSERT(mtsLevels
.size() >= 2, "Need at least two levels for MTS");
129 const int mtsFactor
= mtsLevels
.back().stepFactor
;
130 if (nstValue
% mtsFactor
== 0)
136 return gmx::formatString("With MTS, %s = %d should be a multiple of mts-factor = %d", param
,
137 nstValue
, mtsFactor
);
143 std::vector
<std::string
> checkMtsRequirements(const t_inputrec
& ir
)
145 std::vector
<std::string
> errorMessages
;
149 return errorMessages
;
152 GMX_RELEASE_ASSERT(haveValidMtsSetup(ir
), "MTS setup should be valid here");
154 ArrayRef
<const MtsLevel
> mtsLevels
= ir
.mtsLevels
;
156 if (!(ir
.eI
== eiMD
|| ir
.eI
== eiSD1
))
158 errorMessages
.push_back(gmx::formatString(
159 "Multiple time stepping is only supported with integrators %s and %s",
160 ei_names
[eiMD
], ei_names
[eiSD1
]));
163 if ((EEL_FULL(ir
.coulombtype
) || EVDW_PME(ir
.vdwtype
))
164 && forceGroupMtsLevel(ir
.mtsLevels
, MtsForceGroups::LongrangeNonbonded
) == 0)
166 errorMessages
.emplace_back(
167 "With long-range electrostatics and/or LJ treatment, the long-range part "
168 "has to be part of the mts-level2-forces");
171 std::optional
<std::string
> mesg
;
172 if (ir
.nstcalcenergy
> 0)
174 if ((mesg
= checkMtsInterval(mtsLevels
, "nstcalcenergy", ir
.nstcalcenergy
)))
176 errorMessages
.push_back(mesg
.value());
179 if ((mesg
= checkMtsInterval(mtsLevels
, "nstenergy", ir
.nstenergy
)))
181 errorMessages
.push_back(mesg
.value());
183 if ((mesg
= checkMtsInterval(mtsLevels
, "nstlog", ir
.nstlog
)))
185 errorMessages
.push_back(mesg
.value());
187 if ((mesg
= checkMtsInterval(mtsLevels
, "nstfout", ir
.nstfout
)))
189 errorMessages
.push_back(mesg
.value());
191 if (ir
.efep
!= efepNO
)
193 if ((mesg
= checkMtsInterval(mtsLevels
, "nstdhdl", ir
.fepvals
->nstdhdl
)))
195 errorMessages
.push_back(mesg
.value());
198 if (mtsLevels
.back().forceGroups
[static_cast<int>(gmx::MtsForceGroups::Nonbonded
)])
200 if ((mesg
= checkMtsInterval(mtsLevels
, "nstlist", ir
.nstlist
)))
202 errorMessages
.push_back(mesg
.value());
208 const int pullMtsLevel
= gmx::forceGroupMtsLevel(ir
.mtsLevels
, gmx::MtsForceGroups::Pull
);
209 const int mtsStepFactor
= ir
.mtsLevels
[pullMtsLevel
].stepFactor
;
210 if (ir
.pull
->nstxout
% mtsStepFactor
!= 0)
212 errorMessages
.emplace_back("pull-nstxout should be a multiple of mts-factor");
214 if (ir
.pull
->nstfout
% mtsStepFactor
!= 0)
216 errorMessages
.emplace_back("pull-nstfout should be a multiple of mts-factor");
220 return errorMessages
;