2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2018, 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 * Implements declarations from in mdruncomparison.h.
39 * \author Mark Abraham <mark.j.abraham@gmail.com>
40 * \ingroup module_mdrun_integration_tests
44 #include "mdruncomparison.h"
51 #include "gromacs/utility/stringutil.h"
53 #include "testutils/testasserts.h"
64 using MdpFileValues
= std::map
<std::string
, MdpFieldValues
>;
66 //! Database of .mdp strings that supports prepareDefaultMdpValues()
67 MdpFileValues mdpFileValueDatabase_g
69 // Simple system with 12 argon atoms, fairly widely separated
75 "compressibility", "5e-10"
81 // Simple system with 5 water molecules, fairly widely separated
84 "compressibility", "5e-10"
90 // Simple system with 5832 argon atoms, suitable for normal pressure coupling
96 // Simple system with 216 water molecules, condensed phase
100 // Capped alanine peptide in vacuo with virtual sites
102 "alanine_vsite_vacuo", { {
103 "constraints", "all-bonds"
106 "compressibility", "5e-10"
112 // Capped alanine peptide in aqueous condensed phase, with virtual sites
114 "alanine_vsite_solvated", { {
115 "constraints", "all-bonds"
118 // Nonanol molecule in vacuo, topology suitable for FEP testing
124 "compressibility", "5e-10"
130 "constraints", "h-bonds"
138 init-lambda-state = 3
139 fep_lambdas = 0.00 0.50 1.00 1.00 1.00
140 vdw_lambdas = 0.00 0.00 0.00 0.50 1.00
141 couple-moltype = nonanol
142 couple-lambda0 = vdw-q
143 couple-lambda1 = none
144 couple-intramol = yes)"
149 /*! \brief Prepare default .mdp values
151 * Insert suitable .mdp defaults, so that \c mdpFileValueDatabase_g
152 * does not need to specify repetitive values. This works because
153 * std::map.insert() does not over-write elements that already exist.
155 * \todo ideally some of these default values should be the same as
156 * grompp uses, and sourced from the same place, but that code is a
157 * bit of a jungle until we transition to using IMdpOptions more.
159 * \throws std::bad_alloc if out of memory
160 * std::out_of_range if \c simulationName is not in the database */
161 MdpFieldValues
prepareDefaultMdpFieldValues(const char *simulationName
)
163 using MdpField
= MdpFieldValues::value_type
;
165 auto &mdpFieldValues
= mdpFileValueDatabase_g
.at(simulationName
);
166 mdpFieldValues
.insert(MdpField("nsteps", "16"));
167 mdpFieldValues
.insert(MdpField("ref-t", "298"));
168 mdpFieldValues
.insert(MdpField("tau-p", "1"));
169 mdpFieldValues
.insert(MdpField("compressibility", "5e-5"));
170 mdpFieldValues
.insert(MdpField("constraints", "none"));
171 mdpFieldValues
.insert(MdpField("other", ""));
173 return mdpFieldValues
;
179 prepareMdpFieldValues(const char *simulationName
,
180 const char *integrator
,
184 using MdpField
= MdpFieldValues::value_type
;
186 auto mdpFieldValues
= prepareDefaultMdpFieldValues(simulationName
);
187 mdpFieldValues
.insert(MdpField("integrator", integrator
));
188 mdpFieldValues
.insert(MdpField("tcoupl", tcoupl
));
189 mdpFieldValues
.insert(MdpField("pcoupl", pcoupl
));
190 return mdpFieldValues
;
194 prepareMdpFileContents(const MdpFieldValues
&mdpFieldValues
)
196 /* Set up an .mdp file that permits a highly reproducible
197 * simulation. The format string needs to be configured with
198 * values for various .mdp fields to suit the kind of system
199 * used and testing needed. It also
200 * - writes frames from different kinds of steps: starting, ending, intermediate NS, intermediate non-NS
201 * - has other steps between frame-writing steps
202 * - has enough buffer that e.g. a rerun will compute the same potential energy even though it does NS every frame
203 * - has very slow pressure coupling and weak compressibility (since otherwise the box will shrink too fast)
204 * - can have arbitrary chunks of .mdp content appended to it (but it is up to grompp how it deals with duplicate fields)
205 * - sets random seeds to known values
206 * - uses cutoffs that fit inside boxes even after GPU mdrun scales rlist
208 * Note that forces computed in the absence of energy computations
209 * generally follow a different code path from those computed with
210 * energies. However a rerun always computes energies, so we don't
211 * currently have a good way to compare forces at steps where
212 * energies were not computed with those from rerun on the same
215 return formatString(R
"(rcoulomb = 0.7
219 verlet-buffer-tolerance = 0.000001
233 pcoupltype = isotropic
238 constraint-algorithm = lincs
242 mdpFieldValues
.at("nsteps").c_str(),
243 mdpFieldValues
.at("integrator").c_str(),
244 mdpFieldValues
.at("tcoupl").c_str(),
245 mdpFieldValues
.at("ref-t").c_str(),
246 mdpFieldValues
.at("pcoupl").c_str(),
247 mdpFieldValues
.at("tau-p").c_str(),
248 mdpFieldValues
.at("compressibility").c_str(),
249 mdpFieldValues
.at("constraints").c_str(),
250 mdpFieldValues
.at("other").c_str());