14f63886a502b7c15aa0c739bc48daee9472669f
[alexxy/gromacs.git] / src / programs / mdrun / tests / mdruncomparisonfixture.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements classes in mdruncomparisonfixture.h.
38  *
39  * \author Mark Abraham <mark.j.abraham@gmail.com>
40  * \ingroup module_mdrun_integration_tests
41  */
42 #include "gmxpre.h"
43
44 #include "mdruncomparisonfixture.h"
45
46 #include <map>
47 #include <string>
48 #include <utility>
49 #include <vector>
50
51 #include "gromacs/utility/stringutil.h"
52
53 #include "testutils/testasserts.h"
54
55 namespace gmx
56 {
57 namespace test
58 {
59
60 MdrunComparisonFixture::~MdrunComparisonFixture()
61 {
62 }
63
64 namespace
65 {
66
67 //! Helper typedef
68 typedef std::map<std::string, MdrunComparisonFixture::MdpFieldValues> MdpFileValues;
69
70 //! Database of .mdp strings that supports MdrunComparisonFixture::prepareSimulation()
71 MdpFileValues mdpFileValueDatabase_g
72 {
73     // Simple system with 12 argon atoms, fairly widely separated
74     {
75         "argon12", { {
76                          "ref-t", "80"
77                      },
78                      {
79                          "compressibility", "5e-10"
80                      },
81                      {
82                          "tau-p", "1000"
83                      } }
84     },
85     // Simple system with 5 water molecules, fairly widely separated
86     {
87         "spc5", { {
88                       "compressibility", "5e-10"
89                   },
90                   {
91                       "tau-p", "1000"
92                   } }
93     },
94     // Simple system with 5832 argon atoms, suitable for normal pressure coupling
95     {
96         "argon5832", { {
97                            "ref-t", "80"
98                        } }
99     },
100     // Simple system with 216 water molecules, condensed phase
101     {
102         "spc216", { }
103     },
104     // Capped alanine peptide in vacuo with virtual sites
105     {
106         "alanine_vsite_vacuo", { {
107                                      "constraints", "all-bonds"
108                                  },
109                                  {
110                                      "compressibility", "5e-10"
111                                  },
112                                  {
113                                      "tau-p", "1000"
114                                  } }
115     },
116     // Capped alanine peptide in aqueous condensed phase, with virtual sites
117     {
118         "alanine_vsite_solvated", { {
119                                         "constraints", "all-bonds"
120                                     } }
121     },
122     // Nonanol molecule in vacuo, topology suitable for FEP testing
123     {
124         "nonanol_vacuo", { {
125                                "nsteps", "16"
126                            },
127                            {
128                                "compressibility", "5e-10"
129                            },
130                            {
131                                "tau-p", "1000"
132                            },
133                            {
134                                "constraints", "h-bonds"
135                            },
136                            {
137                                "other",
138                                "free-energy       = yes\n"
139                                "sc-alpha          = 0.5\n"
140                                "sc-r-power        = 6\n"
141                                "nstdhdl           = 4\n"
142                                "init-lambda-state = 3\n"
143                                "fep_lambdas       = 0.00 0.50 1.00 1.00 1.00\n"
144                                "vdw_lambdas       = 0.00 0.00 0.00 0.50 1.00\n"
145                                "couple-moltype    = nonanol\n"
146                                "couple-lambda0    = vdw-q\n"
147                                "couple-lambda1    = none\n"
148                                "couple-intramol   = yes\n"
149                            } }
150     }
151 };
152
153 }       // namespace
154
155 MdrunComparisonFixture::MdpFieldValues MdrunComparisonFixture::prepareMdpFieldValues(const char *simulationName)
156 {
157     /* Insert suitable .mdp defaults, so that the database set up
158      * above does not need to specify redundant values. This works
159      * because std::map.insert() does not over-write elements that
160      * already exist.
161      *
162      * TODO ideally some of these default values should be the same as
163      * grompp uses, and sourced from the same place, but that code is
164      * a bit of a jungle. */
165
166     typedef MdpFieldValues::value_type MdpField;
167
168     auto &mdpFieldValues = mdpFileValueDatabase_g.at(simulationName);
169     mdpFieldValues.insert(MdpField("nsteps", "16"));
170     mdpFieldValues.insert(MdpField("ref-t", "298"));
171     mdpFieldValues.insert(MdpField("tau-p", "1"));
172     mdpFieldValues.insert(MdpField("compressibility", "5e-5"));
173     mdpFieldValues.insert(MdpField("constraints", "none"));
174     mdpFieldValues.insert(MdpField("other", ""));
175
176     return mdpFieldValues;
177 }
178
179 void MdrunComparisonFixture::prepareMdpFile(const MdpFieldValues &mdpFieldValues,
180                                             const char           *integrator,
181                                             const char           *tcoupl,
182                                             const char           *pcoupl)
183 {
184     /* Set up an .mdp file that permits a highly reproducible
185      * simulation. The format string needs to be configured with
186      * values for various .mdp fields to suit the kind of system
187      * used and testing needed. It also
188      * - writes frames from different kinds of steps: starting, ending, intermediate NS, intermediate non-NS
189      * - has other steps between frame-writing steps
190      * - has enough buffer that e.g. a rerun will compute the same potential energy even though it does NS every frame
191      * - has very slow pressure coupling and weak compressibility (since otherwise the box will shrink too fast)
192      * - can have arbitrary chunks of .mdp content appended to it (but it is up to grompp how it deals with duplicate fields)
193      * - sets random seeds to known values
194      * - uses cutoffs that fit inside boxes even after GPU mdrun scales rlist
195      *
196      * Note that forces computed in the absence of energy computations
197      * generally follow a different code path from those computed with
198      * energies. However a rerun always computes energies, so we don't
199      * currently have a good way to compare forces at steps where
200      * energies were not computed with those from rerun on the same
201      * coordinates.
202      */
203     runner_.useStringAsMdpFile(formatString("rcoulomb                = 0.7\n"
204                                             "rvdw                    = 0.7\n"
205                                             "rlist                   = -1\n"
206                                             "bd-fric                 = 1000\n"
207                                             "verlet-buffer-tolerance = 0.000001\n"
208                                             "nsteps                  = %s\n"
209                                             "nstenergy               = 4\n"
210                                             "nstlist                 = 8\n"
211                                             "nstxout                 = 4\n"
212                                             "nstvout                 = 4\n"
213                                             "nstfout                 = 4\n"
214                                             "integrator              = %s\n"
215                                             "ld-seed                 = 234262\n"
216                                             "tcoupl                  = %s\n"
217                                             "ref-t                   = %s\n"
218                                             "tau-t                   = 1\n"
219                                             "tc-grps                 = System\n"
220                                             "pcoupl                  = %s\n"
221                                             "pcoupltype              = isotropic\n"
222                                             "ref-p                   = 1\n"
223                                             "tau-p                   = %s\n"
224                                             "compressibility         = %s\n"
225                                             "constraints             = %s\n"
226                                             "constraint-algorithm    = lincs\n"
227                                             "lincs-order             = 2\n"
228                                             "lincs-iter              = 5\n"
229                                             "%s",
230                                             mdpFieldValues.at("nsteps").c_str(),
231                                             integrator, tcoupl,
232                                             mdpFieldValues.at("ref-t").c_str(),
233                                             pcoupl,
234                                             mdpFieldValues.at("tau-p").c_str(),
235                                             mdpFieldValues.at("compressibility").c_str(),
236                                             mdpFieldValues.at("constraints").c_str(),
237                                             mdpFieldValues.at("other").c_str()));
238 }
239
240 void MdrunComparisonFixture::runTest(const char            *simulationName,
241                                      const char            *integrator,
242                                      const char            *tcoupl,
243                                      const char            *pcoupl,
244                                      FloatingPointTolerance tolerance)
245 {
246     CommandLine caller;
247     caller.append("grompp");
248     runTest(caller, simulationName, integrator, tcoupl, pcoupl, tolerance);
249 }
250
251 } // namespace test
252 } // namespace gmx