Apply clang-format to source tree
[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, 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 <string>
46
47 #include <gmock/gmock.h>
48 #include <gtest/gtest.h>
49
50 #include "gromacs/topology/ifunc.h"
51 #include "gromacs/utility/stringutil.h"
52
53 #include "testutils/refdata.h"
54 #include "testutils/testfilemanager.h"
55
56 #include "energycomparison.h"
57 #include "energyreader.h"
58 #include "moduletest.h"
59
60 namespace gmx
61 {
62 namespace test
63 {
64
65 /*! \internal
66  * \brief End to end test for the density fitting functionality in mdrun.
67  *
68  * This test checks whether the density fitting related .mdp input options are
69  * understood by grompp and that the density fitting module in mdrun computes
70  * the correct density fitting energies and stores them in the .edr output file.
71  */
72 class DensityFittingTest : public MdrunTestFixture
73 {
74 public:
75     DensityFittingTest()
76     {
77         runner_.useTopGroAndNdxFromDatabase("argon12");
78         runner_.edrFileName_ = fileManager_.getTemporaryFilePath(".edr");
79     };
80
81     //! Check the output of mdrun
82     void checkMdrun(real energyTermMagnitude)
83     {
84         const FloatingPointTolerance energyTermTolerance =
85                 relativeToleranceAsFloatingPoint(energyTermMagnitude, 1e-4);
86
87         EnergyTermsToCompare energyTermsToCompare{
88             { { interaction_function[F_DENSITYFITTING].longname, energyTermTolerance },
89               { interaction_function[F_EPOT].longname, energyTermTolerance } }
90         };
91
92         TestReferenceData refData;
93         auto              checker = refData.rootChecker();
94         checkEnergiesAgainstReferenceData(runner_.edrFileName_, energyTermsToCompare, &checker);
95     }
96
97     //! Mdp values for steepest-decent energy minimization with default density fitting parameters.
98     const std::string mdpEminDensfitYesUnsetValues = formatString(
99             "integrator                       = steep\n"
100             "nsteps                           = 2\n"
101             "cutoff-scheme                    = verlet\n"
102             "density-guided-simulation-active = yes\n"
103             "density-guided-simulation-group  = FirstThreeOfTwelve\n"
104             "density-guided-simulation-reference-density-filename = %s\n",
105             TestFileManager::getInputFilePath("ellipsoid-density.mrc").c_str());
106
107     //! Mdp values for md integrator with default density fitting parameters.
108     const std::string mdpMdDensfitYesUnsetValues = formatString(
109             "integrator                       = md\n"
110             "nsteps                           = 2\n"
111             "cutoff-scheme                    = verlet\n"
112             "density-guided-simulation-active = yes\n"
113             "density-guided-simulation-group  = FirstThreeOfTwelve\n"
114             "density-guided-simulation-reference-density-filename = %s\n",
115             TestFileManager::getInputFilePath("ellipsoid-density.mrc").c_str());
116
117     //! Mdp values for steepest-decent energy minimization with density fitting values set to non-defaults.
118     const std::string mdpDensiftAllDefaultsChanged_ = formatString(
119             "density-guided-simulation-similarity-measure = relative-entropy\n"
120             "density-guided-simulation-atom-spreading-weight = mass\n"
121             "density-guided-simulation-force-constant = -1\n"
122             "density-guided-simulation-gaussian-transform-spreading-width = 0.8\n"
123             "density-guided-simulation-gaussian-transform-spreading-range-in-multiples-of-width = "
124             "6\n"
125             "density-guided-simulation-normalize-densities = false\n");
126     //! Set mdp values so that energy calculation interval and density guided simulation interval mismatch.
127     const std::string mdpEnergyAndDensityfittingIntervalMismatch_ = formatString(
128             "nstcalcenergy = 7\n"
129             "density-guided-simulation-nst = 3\n");
130     //! The command line to call mdrun
131     CommandLine commandLineForMdrun_;
132 };
133
134 /* Fit a subset of three of twelve argon atoms into a reference density
135  * whose origin is offset from the simulation box origin.
136  *
137  * All density fitting mdp parameters are set to defaults
138  */
139 TEST_F(DensityFittingTest, EnergyMinimizationEnergyCorrectInnerProduct)
140 {
141     runner_.useStringAsMdpFile(mdpEminDensfitYesUnsetValues);
142
143     ASSERT_EQ(0, runner_.callGrompp());
144     ASSERT_EQ(0, runner_.callMdrun(commandLineForMdrun_));
145
146     const real expectedEnergyTermMagnitude = -3378.825928;
147     checkMdrun(expectedEnergyTermMagnitude);
148 }
149
150 /* Like above, but with as many parameters reversed as possible
151  */
152 TEST_F(DensityFittingTest, EnergyMinimizationEnergyCorrectForRelativeEntropy)
153 {
154     runner_.useStringAsMdpFile(mdpEminDensfitYesUnsetValues + mdpDensiftAllDefaultsChanged_);
155
156     ASSERT_EQ(0, runner_.callGrompp());
157     ASSERT_EQ(0, runner_.callMdrun(commandLineForMdrun_));
158
159     const real expectedEnergyTermMagnitude = -27000;
160     checkMdrun(expectedEnergyTermMagnitude);
161 }
162
163 /* Test that grompp exits with error message if energy evaluation frequencies
164  * do not match.
165  */
166 TEST_F(DensityFittingTest, GromppErrorWhenEnergyEvaluationFrequencyMismatch)
167 {
168     runner_.useStringAsMdpFile(mdpMdDensfitYesUnsetValues + mdpEnergyAndDensityfittingIntervalMismatch_);
169
170     EXPECT_DEATH_IF_SUPPORTED(runner_.callGrompp(),
171                               ".*is not a multiple of density-guided-simulation-nst.*");
172 }
173
174 } // namespace test
175 } // namespace gmx