2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2018,2019, 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 simulationdatabase.h
39 * \author Mark Abraham <mark.j.abraham@gmail.com>
40 * \ingroup module_testutils
44 #include "simulationdatabase.h"
50 #include "gromacs/utility/stringutil.h"
52 #include "testutils/testasserts.h"
64 MdpFieldValues mdpFieldValues
;
65 std::vector
<int> validPpRankCounts
;
69 using MdpFileValues
= std::map
<std::string
, DatabaseEntry
>;
71 //! Database of .mdp strings that supports prepareDefaultMdpValues()
72 const MdpFileValues mdpFileValueDatabase_g
74 // Simple system with 12 argon atoms, fairly widely separated
80 "compressibility", "5e-10"
89 // Simple system with 5 water molecules, fairly widely separated
92 "compressibility", "5e-10"
98 1, 2, 3, 4, 5, 6, 8, 9
101 // Simple system with 5832 argon atoms, suitable for normal pressure coupling
107 // TODO This test case is not currently used, so we
108 // have not tested which rank counts work.
109 1, 2, 3, 4, 5, 6, 7, 8, 9
112 // Simple system with 2 nearby water molecules
116 // TODO This test case is not currently used, so we
117 // have not tested which rank counts work.
118 1, 2, 3, 4, 5, 6, 7, 8, 9
121 // Simple system with 216 water molecules, condensed phase
125 // TODO This test case is not currently used, so we
126 // have not tested which rank counts work.
127 1, 2, 3, 4, 5, 6, 7, 8, 9 // TODO tpi test
130 // Capped alanine peptide in vacuo with virtual sites
132 "alanine_vsite_vacuo", { { {
133 "constraints", "all-bonds"
136 "compressibility", "5e-10"
145 // Capped alanine peptide in aqueous condensed phase, with virtual sites
147 "alanine_vsite_solvated", { { {
148 "constraints", "all-bonds"
151 "compressibility", "5e-10"
157 // TODO This test case is not currently used, so we
158 // have not tested which rank counts work.
159 1, 2, 3, 4, 5, 6, 7, 8, 9
162 // Zwitterionic glycine in vacuo
164 "glycine_vacuo", { { {
165 "constraints", "h-bonds"
168 1, 2, 3, 4, 5, 6, 7, 8, 9
171 // Zwitterionic glycine in vacuo, without constraints
173 "glycine_no_constraints_vacuo", { { {
174 "constraints", "none"
177 1, 2, 3, 4, 5, 6, 7, 8, 9
180 // Simple mdrun tests of energy
187 // Scaled water for NMA
189 "scaled-water", { { },
224 // ICE-Binding protein for NMA
226 "ice-binding", { { },
231 // Nonanol molecule in vacuo, topology suitable for testing FEP
232 // on KE, angles, dihedral restraints, coulomb and vdw
234 "nonanol_vacuo", { { {
238 "compressibility", "5e-10"
244 "constraints", "h-bonds"
251 mass-lambdas = 0.0 0.5 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
252 bonded-lambdas = 0.0 0.0 0.0 0.5 1.0 1.0 1.0 1.0 1.0 1.0 1.0
253 restraint-lambdas = 0.0 0.0 0.0 0.0 0.0 0.5 1.0 1.0 1.0 1.0 1.0
254 vdw-lambdas = 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.5 1.0 1.0 1.0
255 coul-lambdas = 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.5 1.0
256 ;couple-moltype = nonanol
257 ;couple-lambda0 = none
258 ;couple-lambda1 = vdw-q
259 ;couple-intramol = yes)"
262 1, 2, 3, 4, 5, 6, 8, 9
267 /*! \brief Prepare default .mdp values
269 * Insert suitable .mdp defaults, so that \c mdpFileValueDatabase_g
270 * does not need to specify repetitive values. This works because
271 * std::map.insert() does not over-write elements that already exist.
273 * \todo ideally some of these default values should be the same as
274 * grompp uses, and sourced from the same place, but that code is a
275 * bit of a jungle until we transition to using IMdpOptions more.
277 * \throws std::bad_alloc if out of memory
278 * std::out_of_range if \c simulationName is not in the database */
279 MdpFieldValues
prepareDefaultMdpFieldValues(const std::string
&simulationName
)
281 using MdpField
= MdpFieldValues::value_type
;
283 auto mdpFieldValues
= mdpFileValueDatabase_g
.at(simulationName
).mdpFieldValues
;
284 mdpFieldValues
.insert(MdpField("nsteps", "16"));
285 mdpFieldValues
.insert(MdpField("nstenergy", "4"));
286 mdpFieldValues
.insert(MdpField("nstxout", "4"));
287 mdpFieldValues
.insert(MdpField("nstvout", "4"));
288 mdpFieldValues
.insert(MdpField("nstfout", "4"));
289 mdpFieldValues
.insert(MdpField("nstxout-compressed", "0"));
290 mdpFieldValues
.insert(MdpField("nstdhdl", "4"));
291 mdpFieldValues
.insert(MdpField("comm-mode", "linear"));
292 mdpFieldValues
.insert(MdpField("nstcomm", "4"));
293 mdpFieldValues
.insert(MdpField("ref-t", "298"));
294 mdpFieldValues
.insert(MdpField("nsttcouple", "4"));
295 mdpFieldValues
.insert(MdpField("tau-p", "1"));
296 mdpFieldValues
.insert(MdpField("nstpcouple", "4"));
297 mdpFieldValues
.insert(MdpField("compressibility", "5e-5"));
298 mdpFieldValues
.insert(MdpField("constraints", "none"));
299 mdpFieldValues
.insert(MdpField("other", ""));
300 mdpFieldValues
.insert(MdpField("rcoulomb", "0.7"));
301 mdpFieldValues
.insert(MdpField("rvdw", "0.7"));
302 mdpFieldValues
.insert(MdpField("nstcalcenergy", "100"));
304 return mdpFieldValues
;
310 isNumberOfPpRanksSupported(const std::string
&simulationName
,
311 int possibleNumberOfPpRanks
)
313 const auto &possibleNumbers
= mdpFileValueDatabase_g
.at(simulationName
).validPpRankCounts
;
314 return (std::find(std::begin(possibleNumbers
), std::end(possibleNumbers
),
315 possibleNumberOfPpRanks
) != std::end(possibleNumbers
));
319 reportNumbersOfPpRanksSupported(const std::string
&simulationName
)
321 const auto &possibleNumbers
= mdpFileValueDatabase_g
.at(simulationName
).validPpRankCounts
;
322 return formatAndJoin(std::begin(possibleNumbers
), std::end(possibleNumbers
),
323 ",", StringFormatter("%d"));
327 prepareMdpFieldValues(const std::string
&simulationName
,
328 const std::string
&integrator
,
329 const std::string
&tcoupl
,
330 const std::string
&pcoupl
)
332 using MdpField
= MdpFieldValues::value_type
;
334 auto mdpFieldValues
= prepareDefaultMdpFieldValues(simulationName
);
335 mdpFieldValues
.insert(MdpField("integrator", integrator
));
336 mdpFieldValues
.insert(MdpField("tcoupl", tcoupl
));
337 mdpFieldValues
.insert(MdpField("pcoupl", pcoupl
));
338 return mdpFieldValues
;
342 prepareMdpFieldValues(const char *simulationName
,
343 const char *integrator
,
347 return prepareMdpFieldValues(std::string(simulationName
), integrator
, tcoupl
, pcoupl
);
350 prepareMdpFileContents(const MdpFieldValues
&mdpFieldValues
)
352 /* Set up an .mdp file that permits a highly reproducible
353 * simulation. The format string needs to be configured with
354 * values for various .mdp fields to suit the kind of system
355 * used and testing needed. It also
356 * - writes frames from different kinds of steps: starting, ending, intermediate NS, intermediate non-NS
357 * - has other steps between frame-writing steps
358 * - has enough buffer that e.g. a rerun will compute the same potential energy even though it does NS every frame
359 * - has very slow pressure coupling and weak compressibility (since otherwise the box will shrink too fast)
360 * - can have arbitrary chunks of .mdp content appended to it (but it is up to grompp how it deals with duplicate fields)
361 * - sets random seeds to known values
362 * - uses cutoffs that fit inside boxes even after GPU mdrun scales rlist
364 * Note that forces computed in the absence of energy computations
365 * generally follow a different code path from those computed with
366 * energies. However a rerun always computes energies, so we don't
367 * currently have a good way to compare forces at steps where
368 * energies were not computed with those from rerun on the same
371 return formatString(R
"(rcoulomb = %s
375 verlet-buffer-tolerance = 0.000001
381 nstxout-compressed = %s
393 pcoupltype = isotropic
398 constraint-algorithm = lincs
405 mdpFieldValues
.at("rcoulomb").c_str(),
406 mdpFieldValues
.at("rvdw").c_str(),
407 mdpFieldValues
.at("nsteps").c_str(),
408 mdpFieldValues
.at("nstenergy").c_str(),
409 mdpFieldValues
.at("nstxout").c_str(),
410 mdpFieldValues
.at("nstvout").c_str(),
411 mdpFieldValues
.at("nstfout").c_str(),
412 mdpFieldValues
.at("nstxout-compressed").c_str(),
413 mdpFieldValues
.at("nstdhdl").c_str(),
414 mdpFieldValues
.at("integrator").c_str(),
415 mdpFieldValues
.at("tcoupl").c_str(),
416 mdpFieldValues
.at("nsttcouple").c_str(),
417 mdpFieldValues
.at("ref-t").c_str(),
418 mdpFieldValues
.at("pcoupl").c_str(),
419 mdpFieldValues
.at("nstpcouple").c_str(),
420 mdpFieldValues
.at("tau-p").c_str(),
421 mdpFieldValues
.at("compressibility").c_str(),
422 mdpFieldValues
.at("constraints").c_str(),
423 mdpFieldValues
.at("nstcalcenergy").c_str(),
424 mdpFieldValues
.at("comm-mode").c_str(),
425 mdpFieldValues
.at("nstcomm").c_str(),
426 mdpFieldValues
.at("other").c_str());