2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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.
38 * Test for Ewald 3DC and epsilon-surface terms.
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_mdrun_integration_tests
49 #include <gtest/gtest.h>
51 #include "gromacs/topology/idef.h"
52 #include "gromacs/topology/ifunc.h"
53 #include "gromacs/trajectory/trajectoryframe.h"
54 #include "gromacs/utility/basenetwork.h"
55 #include "gromacs/utility/filestream.h"
56 #include "gromacs/utility/strconvert.h"
57 #include "gromacs/utility/stringutil.h"
59 #include "testutils/mpitest.h"
60 #include "testutils/refdata.h"
61 #include "testutils/simulationdatabase.h"
62 #include "testutils/testasserts.h"
64 #include "energycomparison.h"
65 #include "moduletest.h"
66 #include "trajectoryreader.h"
76 using MdpField = MdpFieldValues::value_type;
78 /*! \brief Test fixture base for simple mdrun systems
80 * This test ensures mdrun can run a simulation, reaching reproducible
81 * energies and forces.
82 * The starting coordinates are set up such that both molecules are
83 * broken over PBC and the PBC treatment is tested.
85 class EwaldSurfaceTermTest : public MdrunTestFixture, public ::testing::WithParamInterface<std::string>
89 TEST_P(EwaldSurfaceTermTest, WithinTolerances)
91 auto simulationName = GetParam();
92 SCOPED_TRACE(formatString("Comparing simple mdrun for '%s'", simulationName.c_str()));
94 int numRanksAvailable = getNumberOfTestMpiRanks();
95 /* For epsilon-surface we need whole molecules.
96 * Without constraints we can make molecules whole on a sinlge rank.
97 * With constraints molecules are whole with update groups with DD.
99 * TODO: Remove the rank=1 check when DD also runs single rank.
101 if ((simulationName == "epsilon-surface" && numRanksAvailable > 1)
102 || (simulationName == "epsilon-surface-constraint" && numRanksAvailable == 1))
105 "Test system '%s' cannot run with %d ranks.\n"
106 "The supported numbers are %s 1.\n",
107 simulationName.c_str(),
109 numRanksAvailable == 1 ? ">" : "=");
113 std::string theMdpFile =
114 "coulombtype = PME\n"
115 "nstcalcenergy = 1\n"
121 "fourier-spacing = 0.2\n"
125 if (simulationName == "3DC")
128 "ewald-geometry = 3DC\n"
132 "wall-atomtype = C C\n"
133 "wall-ewald-zfac = 2\n";
137 theMdpFile += "epsilon-surface = 1\n";
138 if (simulationName == "epsilon-surface-constraint")
140 theMdpFile += "constraints = all-bonds";
144 // Prepare the .tpr file
147 runner_.useTopGroAndNdxFromDatabase("dipoles");
148 runner_.useStringAsMdpFile(theMdpFile);
149 EXPECT_EQ(0, runner_.callGrompp(caller));
153 CommandLine mdrunCaller;
154 ASSERT_EQ(0, runner_.callMdrun(mdrunCaller));
155 EnergyTermsToCompare energyTermsToCompare{
156 { { interaction_function[F_EPOT].longname, absoluteTolerance(1e-3) },
157 { interaction_function[F_ETOT].longname, absoluteTolerance(1e-3) } }
159 TestReferenceData refData;
160 auto checker = refData.rootChecker().checkCompound("Simulation", simulationName);
161 checkEnergiesAgainstReferenceData(runner_.edrFileName_, energyTermsToCompare, &checker);
162 // Now check the forces
163 TrajectoryFrameReader reader(runner_.fullPrecisionTrajectoryFileName_);
164 checker.setDefaultTolerance(relativeToleranceAsFloatingPoint(1, 1e-3));
167 auto frame = reader.frame();
168 auto force = frame.f();
170 for (auto& f : force)
172 std::string forceName = frame.frameName() + " F[" + toString(atom) + "]";
174 checker.checkVector(f, forceName.c_str());
177 } while (reader.readNextFrame());
181 //! Containers of systems to test.
183 std::vector<std::string> surfaceTerm = { "3DC", "epsilon-surface-constraint", "epsilon-surface" };
186 INSTANTIATE_TEST_CASE_P(EwaldSurfaceTerm, EwaldSurfaceTermTest, ::testing::ValuesIn(surfaceTerm));