2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2021, 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.
36 * \brief Implements position restraint tests.
38 * \author Kevin Boyd <kevin44boyd@gmail.com>
39 * \ingroup module_listed_forces
43 #include "gromacs/listed_forces/position_restraints.h"
49 #include <gtest/gtest.h>
51 #include "gromacs/gmxlib/nrnb.h"
52 #include "gromacs/listed_forces/listed_forces.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/math/vectypes.h"
55 #include "gromacs/mdtypes/enerdata.h"
56 #include "gromacs/mdtypes/forcerec.h"
57 #include "gromacs/mdtypes/forceoutput.h"
58 #include "gromacs/mdtypes/inputrec.h"
59 #include "gromacs/mdtypes/interaction_const.h"
60 #include "gromacs/mdtypes/mdatom.h"
61 #include "gromacs/mdtypes/nblist.h"
62 #include "gromacs/mdtypes/simulation_workload.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/tables/forcetable.h"
65 #include "gromacs/topology/forcefieldparameters.h"
66 #include "gromacs/topology/idef.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/stringstream.h"
70 #include "testutils/refdata.h"
71 #include "testutils/testasserts.h"
80 //! Tolerance for float evaluation
81 constexpr float c_precisionTolerance = 1e-6;
84 class PositionRestraintsTest : public ::testing::TestWithParam<std::tuple<RefCoordScaling, PbcType>>
96 InteractionDefinitions idef_;
97 gmx_enerdata_t enerd_;
98 std::unique_ptr<ForceWithVirial> forceWithVirial_;
99 RefCoordScaling refCoordScaling_;
101 TestReferenceData refData_;
102 TestReferenceChecker checker_;
104 PositionRestraintsTest() : idef_({}), enerd_(1, 0), checker_(refData_.rootChecker())
106 refCoordScaling_ = std::get<0>(GetParam());
107 pbcType_ = std::get<1>(GetParam());
113 set_pbc(&pbc_, pbcType_, box_);
115 FloatingPointTolerance tolerance = relativeToleranceAsFloatingPoint(1.0, c_precisionTolerance);
116 checker_.setDefaultTolerance(tolerance);
118 fr_.rc_scaling = refCoordScaling_;
119 fr_.pbcType = pbcType_;
120 fr_.posres_com[1] = 0.5;
123 //! Prepares the test with the coordinate and force constant input.
124 void setValues(gmx::ArrayRef<const RVec> positions,
125 gmx::ArrayRef<const RVec> referencePositions,
126 gmx::ArrayRef<const RVec> forceConstants)
128 x_.resize(positions.size());
129 std::copy(positions.begin(), positions.end(), x_.begin());
131 for (size_t i = 0; i < positions.size(); i++)
133 // First item is "type" - each atom will have a different forceparam type
134 // Second item is index - we'll just go from 0.
135 idef_.il[F_POSRES].iatoms.push_back(i);
136 idef_.il[F_POSRES].iatoms.push_back(i);
138 auto& entry = idef_.iparams_posres.emplace_back();
139 copy_rvec(referencePositions[i], entry.posres.pos0A);
140 copy_rvec(forceConstants[i], entry.posres.fcA);
141 clear_rvec(entry.posres.pos0B);
142 clear_rvec(entry.posres.fcB);
144 f_.resize(x_.size(), { 0, 0, 0 });
145 forceWithVirial_ = std::make_unique<ForceWithVirial>(f_, /*computeVirial=*/true);
149 std::array<real, static_cast<size_t>(FreeEnergyPerturbationCouplingType::Count)> c_emptyLambdas = { { 0 } };
151 TEST_P(PositionRestraintsTest, BasicPosResNoFreeEnergy)
153 SCOPED_TRACE(formatString("Testing PBC type: %s, refcoord type: %s",
154 c_pbcTypeNames[pbcType_].c_str(),
155 enumValueToString(refCoordScaling_)));
156 const std::vector<RVec> positions = { { 0.0, 0.0, 0.0 }, { 0.4, 0.5, 0.6 } };
157 const std::vector<RVec> referencePositions = { { 0.0, 0.0, 0.0 }, { 0.5, 0.6, 0.0 } };
158 const std::vector<RVec> forceConstants = { { 1000, 500, 250 }, { 0, 200, 400 } };
159 setValues(positions, referencePositions, forceConstants);
160 posres_wrapper(&nrnb_,
163 as_rvec_array(x_.data()),
167 forceWithVirial_.get());
168 checker_.checkSequence(
169 std::begin(forceWithVirial_->force_), std::end(forceWithVirial_->force_), "Forces");
170 checker_.checkSequenceArray(3, forceWithVirial_->getVirial(), "Virial contribution");
171 checker_.checkReal(enerd_.term[F_POSRES], "Potential energy");
174 //! PBC values for testing
175 std::vector<PbcType> c_pbcForTests = { PbcType::No, PbcType::XY, PbcType::Xyz };
176 //! Reference Coordinate Scaling values for testing
177 std::vector<RefCoordScaling> c_refCoordScalingForTests = { RefCoordScaling::No,
178 RefCoordScaling::Com,
179 RefCoordScaling::All };
181 INSTANTIATE_TEST_SUITE_P(PosResBasicTest,
182 PositionRestraintsTest,
183 ::testing::Combine(::testing::ValuesIn(c_refCoordScalingForTests),
184 ::testing::ValuesIn(c_pbcForTests)));