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