Free topology in non-DD simulations
[gromacs.git] / src / testutils / simulationdatabase.cpp
blob78ad2f8e5d673501d01c9f800308c19c8fbd5f81
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 simulationdatabase.h
39 * \author Mark Abraham <mark.j.abraham@gmail.com>
40 * \ingroup module_testutils
42 #include "gmxpre.h"
44 #include "simulationdatabase.h"
46 #include <algorithm>
47 #include <map>
48 #include <string>
50 #include "gromacs/utility/stringutil.h"
52 #include "testutils/testasserts.h"
54 namespace gmx
56 namespace test
59 namespace
62 struct DatabaseEntry
64 MdpFieldValues mdpFieldValues;
65 std::vector<int> validPpRankCounts;
68 //! Helper typedef
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
76 "argon12", { { {
77 "ref-t", "80"
80 "compressibility", "5e-10"
83 "tau-p", "1000"
84 } },
86 1, 2, 3, 4
87 } }
89 // Simple system with 5 water molecules, fairly widely separated
91 "spc5", { { {
92 "compressibility", "5e-10"
95 "tau-p", "1000"
96 } },
98 1, 2, 3, 4, 5, 6, 8, 9
99 } }
101 // Simple system with 5832 argon atoms, suitable for normal pressure coupling
103 "argon5832", { { {
104 "ref-t", "80"
105 } },
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 216 water molecules, condensed phase
114 "spc216", { { },
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 // TODO tpi test
121 // Capped alanine peptide in vacuo with virtual sites
123 "alanine_vsite_vacuo", { { {
124 "constraints", "all-bonds"
127 "compressibility", "5e-10"
130 "tau-p", "1000"
131 } },
133 1, 2, 3, 4, 6, 9
136 // Capped alanine peptide in aqueous condensed phase, with virtual sites
138 "alanine_vsite_solvated", { { {
139 "constraints", "all-bonds"
142 "compressibility", "5e-10"
145 "tau-p", "1000"
146 } },
148 // TODO This test case is not currently used, so we
149 // have not tested which rank counts work.
150 1, 2, 3, 4, 5, 6, 7, 8, 9
153 // Zwitterionic glycine in vacuo
155 "glycine_vacuo", { { {
156 "constraints", "h-bonds"
157 } },
159 1, 2, 3, 4, 5, 6, 7, 8, 9
162 // Zwitterionic glycine in vacuo, without constraints
164 "glycine_no_constraints_vacuo", { { {
165 "constraints", "none"
166 } },
168 1, 2, 3, 4, 5, 6, 7, 8, 9
171 // Scaled water for NMA
173 "scaled-water", { { },
175 1, 2, 3, 4, 5, 6
178 // Villin for NMA
180 "villin", { { },
182 1, 2, 3, 4, 5, 6
185 // SPC-Dimer for NMA
187 "spc-dimer", { { },
189 1, 2, 3, 4, 5, 6
192 // SW-Dimer for NMA
194 "sw-dimer", { { {
195 "nstcalcenergy", "1"
196 } },
198 1, 2, 3, 4, 5, 6
201 // TIP5P for NMA
203 "one-tip5p", { { },
205 1, 2, 3, 4, 5, 6
208 // ICE-Binding protein for NMA
210 "ice-binding", { { },
212 1, 2, 3, 4, 5, 6
215 // Nonanol molecule in vacuo, topology suitable for testing FEP
216 // on KE, angles, dihedral restraints, coulomb and vdw
218 "nonanol_vacuo", { { {
219 "nsteps", "16"
222 "compressibility", "5e-10"
225 "tau-p", "1000"
228 "constraints", "h-bonds"
231 "other",
232 R"(free-energy = yes
233 sc-alpha = 0.5
234 sc-r-power = 6
235 nstdhdl = 4
236 init-lambda-state = 3
237 mass-lambdas = 0.0 0.5 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
238 bonded-lambdas = 0.0 0.0 0.0 0.5 1.0 1.0 1.0 1.0 1.0 1.0 1.0
239 restraint-lambdas = 0.0 0.0 0.0 0.0 0.0 0.5 1.0 1.0 1.0 1.0 1.0
240 vdw-lambdas = 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.5 1.0 1.0 1.0
241 coul-lambdas = 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.5 1.0
242 ;couple-moltype = nonanol
243 ;couple-lambda0 = none
244 ;couple-lambda1 = vdw-q
245 ;couple-intramol = yes)"
246 } },
248 1, 2, 3, 4, 5, 6, 8, 9
253 /*! \brief Prepare default .mdp values
255 * Insert suitable .mdp defaults, so that \c mdpFileValueDatabase_g
256 * does not need to specify repetitive values. This works because
257 * std::map.insert() does not over-write elements that already exist.
259 * \todo ideally some of these default values should be the same as
260 * grompp uses, and sourced from the same place, but that code is a
261 * bit of a jungle until we transition to using IMdpOptions more.
263 * \throws std::bad_alloc if out of memory
264 * std::out_of_range if \c simulationName is not in the database */
265 MdpFieldValues prepareDefaultMdpFieldValues(const char *simulationName)
267 using MdpField = MdpFieldValues::value_type;
269 auto mdpFieldValues = mdpFileValueDatabase_g.at(simulationName).mdpFieldValues;
270 mdpFieldValues.insert(MdpField("nsteps", "16"));
271 mdpFieldValues.insert(MdpField("ref-t", "298"));
272 mdpFieldValues.insert(MdpField("tau-p", "1"));
273 mdpFieldValues.insert(MdpField("compressibility", "5e-5"));
274 mdpFieldValues.insert(MdpField("constraints", "none"));
275 mdpFieldValues.insert(MdpField("other", ""));
276 mdpFieldValues.insert(MdpField("rcoulomb", "0.7"));
277 mdpFieldValues.insert(MdpField("rvdw", "0.7"));
278 mdpFieldValues.insert(MdpField("nstcalcenergy", "100"));
280 return mdpFieldValues;
283 } // namespace
285 bool
286 isNumberOfPpRanksSupported(const std::string &simulationName,
287 int possibleNumberOfPpRanks)
289 const auto &possibleNumbers = mdpFileValueDatabase_g.at(simulationName).validPpRankCounts;
290 return (std::find(std::begin(possibleNumbers), std::end(possibleNumbers),
291 possibleNumberOfPpRanks) != std::end(possibleNumbers));
294 std::string
295 reportNumbersOfPpRanksSupported(const std::string &simulationName)
297 const auto &possibleNumbers = mdpFileValueDatabase_g.at(simulationName).validPpRankCounts;
298 return formatAndJoin(std::begin(possibleNumbers), std::end(possibleNumbers),
299 ",", StringFormatter("%d"));
302 MdpFieldValues
303 prepareMdpFieldValues(const char *simulationName,
304 const char *integrator,
305 const char *tcoupl,
306 const char *pcoupl)
308 using MdpField = MdpFieldValues::value_type;
310 auto mdpFieldValues = prepareDefaultMdpFieldValues(simulationName);
311 mdpFieldValues.insert(MdpField("integrator", integrator));
312 mdpFieldValues.insert(MdpField("tcoupl", tcoupl));
313 mdpFieldValues.insert(MdpField("pcoupl", pcoupl));
314 return mdpFieldValues;
317 std::string
318 prepareMdpFileContents(const MdpFieldValues &mdpFieldValues)
320 /* Set up an .mdp file that permits a highly reproducible
321 * simulation. The format string needs to be configured with
322 * values for various .mdp fields to suit the kind of system
323 * used and testing needed. It also
324 * - writes frames from different kinds of steps: starting, ending, intermediate NS, intermediate non-NS
325 * - has other steps between frame-writing steps
326 * - has enough buffer that e.g. a rerun will compute the same potential energy even though it does NS every frame
327 * - has very slow pressure coupling and weak compressibility (since otherwise the box will shrink too fast)
328 * - can have arbitrary chunks of .mdp content appended to it (but it is up to grompp how it deals with duplicate fields)
329 * - sets random seeds to known values
330 * - uses cutoffs that fit inside boxes even after GPU mdrun scales rlist
332 * Note that forces computed in the absence of energy computations
333 * generally follow a different code path from those computed with
334 * energies. However a rerun always computes energies, so we don't
335 * currently have a good way to compare forces at steps where
336 * energies were not computed with those from rerun on the same
337 * coordinates.
339 return formatString(R"(rcoulomb = %s
340 rvdw = %s
341 rlist = -1
342 bd-fric = 1000
343 verlet-buffer-tolerance = 0.000001
344 nsteps = %s
345 nstenergy = 4
346 nstlist = 8
347 nstxout = 4
348 nstvout = 4
349 nstfout = 4
350 integrator = %s
351 ld-seed = 234262
352 tcoupl = %s
353 ref-t = %s
354 tau-t = 1
355 tc-grps = System
356 pcoupl = %s
357 pcoupltype = isotropic
358 ref-p = 1
359 tau-p = %s
360 compressibility = %s
361 constraints = %s
362 constraint-algorithm = lincs
363 lincs-order = 2
364 lincs-iter = 5
365 nstcalcenergy = %s
366 %s)",
367 mdpFieldValues.at("rcoulomb").c_str(),
368 mdpFieldValues.at("rvdw").c_str(),
369 mdpFieldValues.at("nsteps").c_str(),
370 mdpFieldValues.at("integrator").c_str(),
371 mdpFieldValues.at("tcoupl").c_str(),
372 mdpFieldValues.at("ref-t").c_str(),
373 mdpFieldValues.at("pcoupl").c_str(),
374 mdpFieldValues.at("tau-p").c_str(),
375 mdpFieldValues.at("compressibility").c_str(),
376 mdpFieldValues.at("constraints").c_str(),
377 mdpFieldValues.at("nstcalcenergy").c_str(),
378 mdpFieldValues.at("other").c_str());
381 } // namespace test
382 } // namespace gmx