7149eee08b3d87a291fdd7143e6900c34bcc08a2
[alexxy/gromacs.git] / src / programs / mdrun / tests / densityfittingmodule.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020, 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
36 /*! \internal \file
37  * \brief
38  * Tests utilities for "Densityfitting" setups.
39  *
40  * \author Christian Blau <blau@kth.se>
41  * \ingroup module_mdrun_integration_tests
42  */
43 #include "gmxpre.h"
44
45 #include "config.h"
46
47 #include <string>
48
49 #include <gmock/gmock.h>
50 #include <gtest/gtest.h>
51
52 #include "gromacs/topology/ifunc.h"
53 #include "gromacs/utility/stringutil.h"
54
55 #include "testutils/refdata.h"
56 #include "testutils/testfilemanager.h"
57
58 #include "energycomparison.h"
59 #include "energyreader.h"
60 #include "moduletest.h"
61
62 namespace gmx
63 {
64 namespace test
65 {
66
67 /*! \internal
68  * \brief End to end test for the density fitting functionality in mdrun.
69  *
70  * This test checks whether the density fitting related .mdp input options are
71  * understood by grompp and that the density fitting module in mdrun computes
72  * the correct density fitting energies and stores them in the .edr output file.
73  */
74 class DensityFittingTest : public MdrunTestFixture
75 {
76 public:
77     DensityFittingTest()
78     {
79         runner_.useTopGroAndNdxFromDatabase("argon12");
80         runner_.edrFileName_ = fileManager_.getTemporaryFilePath(".edr");
81     }
82
83     //! Check the output of mdrun
84     void checkMdrun(real energyTermMagnitude)
85     {
86         const FloatingPointTolerance energyTermTolerance =
87                 relativeToleranceAsFloatingPoint(energyTermMagnitude, 1e-4);
88
89         EnergyTermsToCompare energyTermsToCompare{
90             { { interaction_function[F_DENSITYFITTING].longname, energyTermTolerance },
91               { interaction_function[F_EPOT].longname, energyTermTolerance } }
92         };
93
94         TestReferenceData refData;
95         auto              checker = refData.rootChecker();
96         checkEnergiesAgainstReferenceData(runner_.edrFileName_, energyTermsToCompare, &checker);
97     }
98
99     //! Mdp values for steepest-decent energy minimization with default density fitting parameters.
100     const std::string mdpEminDensfitYesUnsetValues = formatString(
101             "integrator                       = steep\n"
102             "nsteps                           = 2\n"
103             "cutoff-scheme                    = verlet\n"
104             "density-guided-simulation-active = yes\n"
105             "density-guided-simulation-group  = FirstThreeOfTwelve\n"
106             "density-guided-simulation-reference-density-filename = %s\n",
107             TestFileManager::getInputFilePath("ellipsoid-density.mrc").c_str());
108
109     //! Mdp values for md integrator with default density fitting parameters.
110     const std::string mdpMdDensfitYesUnsetValues = formatString(
111             "integrator                       = md\n"
112             "nsteps                           = 2\n"
113             "cutoff-scheme                    = verlet\n"
114             "density-guided-simulation-active = yes\n"
115             "density-guided-simulation-group  = FirstThreeOfTwelve\n"
116             "density-guided-simulation-reference-density-filename = %s\n",
117             TestFileManager::getInputFilePath("ellipsoid-density.mrc").c_str());
118
119     //! Mdp values for steepest-decent energy minimization with density fitting values set to non-defaults.
120     const std::string mdpDensiftAllDefaultsChanged_ = formatString(
121             "density-guided-simulation-similarity-measure = relative-entropy\n"
122             "density-guided-simulation-atom-spreading-weight = mass\n"
123             "density-guided-simulation-force-constant = -1\n"
124             "density-guided-simulation-gaussian-transform-spreading-width = 0.8\n"
125             "density-guided-simulation-gaussian-transform-spreading-range-in-multiples-of-width = "
126             "6\n"
127             "density-guided-simulation-normalize-densities = false\n");
128     //! Set mdp values so that energy calculation interval and density guided simulation interval mismatch.
129     const std::string mdpEnergyAndDensityfittingIntervalMismatch_ = formatString(
130             "nstcalcenergy = 7\n"
131             "density-guided-simulation-nst = 3\n");
132     //! Set mdp values so that we skip every other step
133     const std::string mdpSkipDensityfittingEveryOtherStep_ = formatString(
134             "nstenergy = 2\n"
135             "density-guided-simulation-nst = 2\n");
136     //! The command line to call mdrun
137     CommandLine commandLineForMdrun_;
138 };
139
140 /* Fit a subset of three of twelve argon atoms into a reference density
141  * whose origin is offset from the simulation box origin.
142  *
143  * All density fitting mdp parameters are set to defaults
144  */
145 TEST_F(DensityFittingTest, EnergyMinimizationEnergyCorrectInnerProduct)
146 {
147     runner_.useStringAsMdpFile(mdpEminDensfitYesUnsetValues);
148
149     ASSERT_EQ(0, runner_.callGrompp());
150     ASSERT_EQ(0, runner_.callMdrun(commandLineForMdrun_));
151
152     const real expectedEnergyTermMagnitude = -3378.825928;
153     checkMdrun(expectedEnergyTermMagnitude);
154 }
155
156 /* Like above, but with as many parameters reversed as possible
157  */
158 TEST_F(DensityFittingTest, EnergyMinimizationEnergyCorrectForRelativeEntropy)
159 {
160     runner_.useStringAsMdpFile(mdpEminDensfitYesUnsetValues + mdpDensiftAllDefaultsChanged_);
161
162     ASSERT_EQ(0, runner_.callGrompp());
163     ASSERT_EQ(0, runner_.callMdrun(commandLineForMdrun_));
164
165     const real expectedEnergyTermMagnitude = -27000;
166     checkMdrun(expectedEnergyTermMagnitude);
167 }
168
169 /* Test that grompp exits with error message if energy evaluation frequencies
170  * do not match.
171  */
172 TEST_F(DensityFittingTest, GromppErrorWhenEnergyEvaluationFrequencyMismatch)
173 {
174     runner_.useStringAsMdpFile(mdpMdDensfitYesUnsetValues + mdpEnergyAndDensityfittingIntervalMismatch_);
175
176     GMX_EXPECT_DEATH_IF_SUPPORTED(runner_.callGrompp(),
177                                   ".*is not a multiple of density-guided-simulation-nst.*");
178 }
179
180 /* Fit a subset of three of twelve argon atoms into a reference density
181  * whose origin is offset from the simulation box origin. Stop the simulation,
182  * then restart.
183  *
184  * All density fitting mdp parameters are set to defaults
185  */
186 TEST_F(DensityFittingTest, CheckpointWorks)
187 {
188     runner_.useStringAsMdpFile(mdpMdDensfitYesUnsetValues + mdpSkipDensityfittingEveryOtherStep_);
189     runner_.cptFileName_ = fileManager_.getTemporaryFilePath(".cpt");
190     commandLineForMdrun_.addOption("-cpo", runner_.cptFileName_);
191
192     ASSERT_EQ(0, runner_.callGrompp());
193     ASSERT_EQ(0, runner_.callMdrun(commandLineForMdrun_));
194
195     // checkMdrun(expectedEnergyTermMagnitude);
196
197     CommandLine commandLineForRestart;
198     commandLineForRestart.addOption("-cpi", runner_.cptFileName_);
199     commandLineForRestart.addOption("-noappend");
200     runner_.nsteps_ = 4;
201     ASSERT_EQ(0, runner_.callMdrun(commandLineForRestart));
202
203     const real expectedEnergyTermMagnitude = -3378.825928;
204     checkMdrun(expectedEnergyTermMagnitude);
205 }
206
207
208 } // namespace test
209 } // namespace gmx