Re-implemented mdrun -rerun tests
[gromacs.git] / src / programs / mdrun / tests / mdruncomparison.cpp
blob736b739afce06120e5e9bf679fa45afb3c16727a
1 /*
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.
35 /*! \internal \file
36 * \brief
37 * Implements declarations from in mdruncomparison.h.
39 * \author Mark Abraham <mark.j.abraham@gmail.com>
40 * \ingroup module_mdrun_integration_tests
42 #include "gmxpre.h"
44 #include "mdruncomparison.h"
46 #include <map>
47 #include <string>
48 #include <utility>
49 #include <vector>
51 #include "gromacs/utility/stringutil.h"
53 #include "testutils/testasserts.h"
55 namespace gmx
57 namespace test
60 namespace
63 //! Helper typedef
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
71 "argon12", { {
72 "ref-t", "80"
75 "compressibility", "5e-10"
78 "tau-p", "1000"
79 } }
81 // Simple system with 5 water molecules, fairly widely separated
83 "spc5", { {
84 "compressibility", "5e-10"
87 "tau-p", "1000"
88 } }
90 // Simple system with 5832 argon atoms, suitable for normal pressure coupling
92 "argon5832", { {
93 "ref-t", "80"
94 } }
96 // Simple system with 216 water molecules, condensed phase
98 "spc216", { }
100 // Capped alanine peptide in vacuo with virtual sites
102 "alanine_vsite_vacuo", { {
103 "constraints", "all-bonds"
106 "compressibility", "5e-10"
109 "tau-p", "1000"
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 testing FEP
119 // on KE, angles, dihedral restraints, coulomb and vdw
121 "nonanol_vacuo", { {
122 "nsteps", "16"
125 "compressibility", "5e-10"
128 "tau-p", "1000"
131 "constraints", "h-bonds"
134 "other",
135 R"(free-energy = yes
136 sc-alpha = 0.5
137 sc-r-power = 6
138 nstdhdl = 4
139 init-lambda-state = 3
140 mass-lambdas = 0.0 0.5 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
141 bonded-lambdas = 0.0 0.0 0.0 0.5 1.0 1.0 1.0 1.0 1.0 1.0 1.0
142 restraint-lambdas = 0.0 0.0 0.0 0.0 0.0 0.5 1.0 1.0 1.0 1.0 1.0
143 vdw-lambdas = 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.5 1.0 1.0 1.0
144 coul-lambdas = 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.5 1.0
145 ;couple-moltype = nonanol
146 ;couple-lambda0 = none
147 ;couple-lambda1 = vdw-q
148 ;couple-intramol = yes)"
153 /*! \brief Prepare default .mdp values
155 * Insert suitable .mdp defaults, so that \c mdpFileValueDatabase_g
156 * does not need to specify repetitive values. This works because
157 * std::map.insert() does not over-write elements that already exist.
159 * \todo ideally some of these default values should be the same as
160 * grompp uses, and sourced from the same place, but that code is a
161 * bit of a jungle until we transition to using IMdpOptions more.
163 * \throws std::bad_alloc if out of memory
164 * std::out_of_range if \c simulationName is not in the database */
165 MdpFieldValues prepareDefaultMdpFieldValues(const char *simulationName)
167 using MdpField = MdpFieldValues::value_type;
169 auto &mdpFieldValues = mdpFileValueDatabase_g.at(simulationName);
170 mdpFieldValues.insert(MdpField("nsteps", "16"));
171 mdpFieldValues.insert(MdpField("ref-t", "298"));
172 mdpFieldValues.insert(MdpField("tau-p", "1"));
173 mdpFieldValues.insert(MdpField("compressibility", "5e-5"));
174 mdpFieldValues.insert(MdpField("constraints", "none"));
175 mdpFieldValues.insert(MdpField("other", ""));
177 return mdpFieldValues;
180 } // namespace
182 MdpFieldValues
183 prepareMdpFieldValues(const char *simulationName,
184 const char *integrator,
185 const char *tcoupl,
186 const char *pcoupl)
188 using MdpField = MdpFieldValues::value_type;
190 auto mdpFieldValues = prepareDefaultMdpFieldValues(simulationName);
191 mdpFieldValues.insert(MdpField("integrator", integrator));
192 mdpFieldValues.insert(MdpField("tcoupl", tcoupl));
193 mdpFieldValues.insert(MdpField("pcoupl", pcoupl));
194 return mdpFieldValues;
197 std::string
198 prepareMdpFileContents(const MdpFieldValues &mdpFieldValues)
200 /* Set up an .mdp file that permits a highly reproducible
201 * simulation. The format string needs to be configured with
202 * values for various .mdp fields to suit the kind of system
203 * used and testing needed. It also
204 * - writes frames from different kinds of steps: starting, ending, intermediate NS, intermediate non-NS
205 * - has other steps between frame-writing steps
206 * - has enough buffer that e.g. a rerun will compute the same potential energy even though it does NS every frame
207 * - has very slow pressure coupling and weak compressibility (since otherwise the box will shrink too fast)
208 * - can have arbitrary chunks of .mdp content appended to it (but it is up to grompp how it deals with duplicate fields)
209 * - sets random seeds to known values
210 * - uses cutoffs that fit inside boxes even after GPU mdrun scales rlist
212 * Note that forces computed in the absence of energy computations
213 * generally follow a different code path from those computed with
214 * energies. However a rerun always computes energies, so we don't
215 * currently have a good way to compare forces at steps where
216 * energies were not computed with those from rerun on the same
217 * coordinates.
219 return formatString(R"(rcoulomb = 0.7
220 rvdw = 0.7
221 rlist = -1
222 bd-fric = 1000
223 verlet-buffer-tolerance = 0.000001
224 nsteps = %s
225 nstenergy = 4
226 nstlist = 8
227 nstxout = 4
228 nstvout = 4
229 nstfout = 4
230 integrator = %s
231 ld-seed = 234262
232 tcoupl = %s
233 ref-t = %s
234 tau-t = 1
235 tc-grps = System
236 pcoupl = %s
237 pcoupltype = isotropic
238 ref-p = 1
239 tau-p = %s
240 compressibility = %s
241 constraints = %s
242 constraint-algorithm = lincs
243 lincs-order = 2
244 lincs-iter = 5
245 %s)",
246 mdpFieldValues.at("nsteps").c_str(),
247 mdpFieldValues.at("integrator").c_str(),
248 mdpFieldValues.at("tcoupl").c_str(),
249 mdpFieldValues.at("ref-t").c_str(),
250 mdpFieldValues.at("pcoupl").c_str(),
251 mdpFieldValues.at("tau-p").c_str(),
252 mdpFieldValues.at("compressibility").c_str(),
253 mdpFieldValues.at("constraints").c_str(),
254 mdpFieldValues.at("other").c_str());
257 } // namespace test
258 } // namespace gmx