2 * This file is part of the GROMACS molecular simulation package.
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.
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.
37 * Tests for density fitting module options.
39 * \author Christian Blau <blau@kth.se>
40 * \ingroup module_applied_forces
44 #include "gromacs/applied_forces/densityfitting/densityfittingoptions.h"
49 #include <gtest/gtest.h>
51 #include "gromacs/options/options.h"
52 #include "gromacs/options/treesupport.h"
53 #include "gromacs/selection/indexutil.h"
54 #include "gromacs/utility/keyvaluetreebuilder.h"
55 #include "gromacs/utility/keyvaluetreemdpwriter.h"
56 #include "gromacs/utility/keyvaluetreetransform.h"
57 #include "gromacs/utility/smalloc.h"
58 #include "gromacs/utility/stringcompare.h"
59 #include "gromacs/utility/stringstream.h"
60 #include "gromacs/utility/textwriter.h"
62 #include "testutils/testasserts.h"
63 #include "testutils/testmatchers.h"
71 class DensityFittingOptionsTest : public ::testing::Test
74 DensityFittingOptionsTest() { init_blocka(&defaultGroups_); }
75 ~DensityFittingOptionsTest() override { done_blocka(&defaultGroups_); }
77 void setFromMdpValues(const KeyValueTreeObject& densityFittingMdpValues)
80 Options densityFittingModuleOptions;
81 densityFittingOptions_.initMdpOptions(&densityFittingModuleOptions);
83 // Add rules to transform mdp inputs to densityFittingModule data
84 KeyValueTreeTransformer transform;
85 transform.rules()->addRule().keyMatchType("/", StringCompareType::CaseAndDashInsensitive);
87 densityFittingOptions_.initMdpTransform(transform.rules());
89 // Execute the transform on the mdpValues
90 auto transformedMdpValues = transform.transform(densityFittingMdpValues, nullptr);
91 assignOptionsFromKeyValueTree(&densityFittingModuleOptions, transformedMdpValues.object(), nullptr);
94 static KeyValueTreeObject densityFittingSetActiveAsMdpValues()
97 KeyValueTreeBuilder mdpValueBuilder;
98 mdpValueBuilder.rootObject().addValue("density-guided-simulation-active", std::string("yes"));
99 return mdpValueBuilder.build();
102 IndexGroupsAndNames genericIndexGroupsAndNames()
104 done_blocka(&defaultGroups_);
105 stupid_fill_blocka(&defaultGroups_, 3);
106 std::vector<std::string> groupNames = { "A", "protein", "C" };
107 const char* const namesAsConstChar[3] = { groupNames[0].c_str(),
108 groupNames[1].c_str(),
109 groupNames[2].c_str() };
110 return { defaultGroups_, namesAsConstChar };
113 IndexGroupsAndNames differingIndexGroupsAndNames()
115 done_blocka(&defaultGroups_);
116 stupid_fill_blocka(&defaultGroups_, 3);
117 std::vector<std::string> groupNames = { "protein", "C", "A" };
118 const char* const namesAsConstChar[3] = { groupNames[0].c_str(),
119 groupNames[1].c_str(),
120 groupNames[2].c_str() };
121 return { defaultGroups_, namesAsConstChar };
124 void mangleInternalParameters()
126 densityFittingOptions_.setFitGroupIndices(differingIndexGroupsAndNames());
130 t_blocka defaultGroups_;
131 DensityFittingOptions densityFittingOptions_;
134 TEST_F(DensityFittingOptionsTest, DefaultParameters)
136 const auto defaultParameters = densityFittingOptions_.buildParameters();
137 EXPECT_FALSE(defaultParameters.active_);
138 EXPECT_EQ(0, defaultParameters.indices_.size());
139 EXPECT_EQ(DensitySimilarityMeasureMethod::innerProduct, defaultParameters.similarityMeasureMethod_);
140 EXPECT_EQ(DensityFittingAmplitudeMethod::Unity, defaultParameters.amplitudeLookupMethod_);
141 EXPECT_REAL_EQ(1e9, defaultParameters.forceConstant_);
142 EXPECT_REAL_EQ(0.2, defaultParameters.gaussianTransformSpreadingWidth_);
143 EXPECT_REAL_EQ(4.0, defaultParameters.gaussianTransformSpreadingRangeInMultiplesOfWidth_);
146 TEST_F(DensityFittingOptionsTest, OptionSetsActive)
148 EXPECT_FALSE(densityFittingOptions_.buildParameters().active_);
149 setFromMdpValues(densityFittingSetActiveAsMdpValues());
150 EXPECT_TRUE(densityFittingOptions_.buildParameters().active_);
153 TEST_F(DensityFittingOptionsTest, OutputNoDefaultValuesWhenInactive)
155 // Transform module data into a flat key-value tree for output.
157 StringOutputStream stream;
158 KeyValueTreeBuilder builder;
159 KeyValueTreeObjectBuilder builderObject = builder.rootObject();
161 densityFittingOptions_.buildMdpOutput(&builderObject);
163 TextWriter writer(&stream);
164 writeKeyValueTreeAsMdp(&writer, builder.build());
168 EXPECT_EQ(stream.toString(),
170 "\n; Density guided simulation\ndensity-guided-simulation-active = false\n"));
173 TEST_F(DensityFittingOptionsTest, OutputDefaultValuesWhenActive)
175 setFromMdpValues(densityFittingSetActiveAsMdpValues());
176 // Transform module data into a flat key-value tree for output.
178 StringOutputStream stream;
179 KeyValueTreeBuilder builder;
180 KeyValueTreeObjectBuilder builderObject = builder.rootObject();
182 densityFittingOptions_.buildMdpOutput(&builderObject);
184 TextWriter writer(&stream);
185 writeKeyValueTreeAsMdp(&writer, builder.build());
188 std::string expectedString = {
190 "; Density guided simulation\n"
191 "density-guided-simulation-active = true\n"
192 "density-guided-simulation-group = protein\n"
193 "; Similarity measure between densities: inner-product, relative-entropy, or "
194 "cross-correlation\n"
195 "density-guided-simulation-similarity-measure = inner-product\n"
196 "; Atom amplitude for spreading onto grid: unity, mass, or charge\n"
197 "density-guided-simulation-atom-spreading-weight = unity\n"
198 "density-guided-simulation-force-constant = 1e+09\n"
199 "density-guided-simulation-gaussian-transform-spreading-width = 0.2\n"
200 "density-guided-simulation-gaussian-transform-spreading-range-in-multiples-of-width = 4\n"
201 "; Reference density file location as absolute path or relative to the gmx mdrun calling "
203 "density-guided-simulation-reference-density-filename = reference.mrc\n"
204 "density-guided-simulation-nst = 1\n"
205 "; Normalize the sum of density voxel values to one\n"
206 "density-guided-simulation-normalize-densities = true\n"
207 "; Apply adaptive force scaling\n"
208 "density-guided-simulation-adaptive-force-scaling = false\n"
209 "; Time constant for adaptive force scaling in ps\n"
210 "density-guided-simulation-adaptive-force-scaling-time-constant = 4\n"
213 EXPECT_EQ(expectedString, stream.toString());
217 TEST_F(DensityFittingOptionsTest, CanConvertGroupStringToIndexGroup)
219 setFromMdpValues(densityFittingSetActiveAsMdpValues());
221 const auto indexGroupAndNames = genericIndexGroupsAndNames();
222 densityFittingOptions_.setFitGroupIndices(indexGroupAndNames);
224 EXPECT_EQ(1, densityFittingOptions_.buildParameters().indices_.size());
225 EXPECT_EQ(1, densityFittingOptions_.buildParameters().indices_[0]);
228 TEST_F(DensityFittingOptionsTest, InternalsToKvt)
230 // stores the default internal options
231 DensityFittingOptions densityFittingOptions;
232 KeyValueTreeBuilder builder;
233 densityFittingOptions.writeInternalParametersToKvt(builder.rootObject());
234 const auto kvtTree = builder.build();
235 EXPECT_TRUE(kvtTree.keyExists("density-guided-simulation-group"));
236 EXPECT_TRUE(kvtTree["density-guided-simulation-group"].isArray());
237 auto storedIndex = kvtTree["density-guided-simulation-group"].asArray().values();
239 EXPECT_EQ(0, storedIndex.size());
242 TEST_F(DensityFittingOptionsTest, KvtToInternal)
244 setFromMdpValues(densityFittingSetActiveAsMdpValues());
246 KeyValueTreeBuilder builder;
248 builder.rootObject().addUniformArray<std::int64_t>("density-guided-simulation-group");
249 addedArray.addValue(1);
250 addedArray.addValue(15);
251 const auto tree = builder.build();
253 densityFittingOptions_.readInternalParametersFromKvt(tree);
255 EXPECT_EQ(2, densityFittingOptions_.buildParameters().indices_.size());
256 EXPECT_EQ(1, densityFittingOptions_.buildParameters().indices_[0]);
257 EXPECT_EQ(15, densityFittingOptions_.buildParameters().indices_[1]);
260 TEST_F(DensityFittingOptionsTest, RoundTripForInternalsIsIdempotent)
262 setFromMdpValues(densityFittingSetActiveAsMdpValues());
264 const IndexGroupsAndNames indexGroupAndNames = genericIndexGroupsAndNames();
265 densityFittingOptions_.setFitGroupIndices(indexGroupAndNames);
268 DensityFittingParameters parametersBefore = densityFittingOptions_.buildParameters();
270 KeyValueTreeBuilder builder;
271 densityFittingOptions_.writeInternalParametersToKvt(builder.rootObject());
272 const auto inputTree = builder.build();
274 mangleInternalParameters();
276 DensityFittingParameters parametersAfter = densityFittingOptions_.buildParameters();
277 EXPECT_NE(parametersBefore, parametersAfter);
279 densityFittingOptions_.readInternalParametersFromKvt(inputTree);
281 parametersAfter = densityFittingOptions_.buildParameters();
282 EXPECT_EQ(parametersBefore, parametersAfter);