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/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",
100 return mdpValueBuilder
.build();
103 IndexGroupsAndNames
genericIndexGroupsAndNames()
105 done_blocka(&defaultGroups_
);
106 stupid_fill_blocka(&defaultGroups_
, 3);
107 std::vector
<std::string
> groupNames
= { "A", "protein", "C" };
108 const char* const namesAsConstChar
[3] = { groupNames
[0].c_str(), 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(), groupNames
[1].c_str(),
119 groupNames
[2].c_str() };
120 return { defaultGroups_
, namesAsConstChar
};
123 void mangleInternalParameters()
125 densityFittingOptions_
.setFitGroupIndices(differingIndexGroupsAndNames());
129 t_blocka defaultGroups_
;
130 DensityFittingOptions densityFittingOptions_
;
133 TEST_F(DensityFittingOptionsTest
, DefaultParameters
)
135 const auto defaultParameters
= densityFittingOptions_
.buildParameters();
136 EXPECT_FALSE(defaultParameters
.active_
);
137 EXPECT_EQ(0, defaultParameters
.indices_
.size());
138 EXPECT_EQ(DensitySimilarityMeasureMethod::innerProduct
, defaultParameters
.similarityMeasureMethod_
);
139 EXPECT_EQ(DensityFittingAmplitudeMethod::Unity
, defaultParameters
.amplitudeLookupMethod_
);
140 EXPECT_REAL_EQ(1e9
, defaultParameters
.forceConstant_
);
141 EXPECT_REAL_EQ(0.2, defaultParameters
.gaussianTransformSpreadingWidth_
);
142 EXPECT_REAL_EQ(4.0, defaultParameters
.gaussianTransformSpreadingRangeInMultiplesOfWidth_
);
145 TEST_F(DensityFittingOptionsTest
, OptionSetsActive
)
147 EXPECT_FALSE(densityFittingOptions_
.buildParameters().active_
);
148 setFromMdpValues(densityFittingSetActiveAsMdpValues());
149 EXPECT_TRUE(densityFittingOptions_
.buildParameters().active_
);
152 TEST_F(DensityFittingOptionsTest
, OutputNoDefaultValuesWhenInactive
)
154 // Transform module data into a flat key-value tree for output.
156 StringOutputStream stream
;
157 KeyValueTreeBuilder builder
;
158 KeyValueTreeObjectBuilder builderObject
= builder
.rootObject();
160 densityFittingOptions_
.buildMdpOutput(&builderObject
);
162 TextWriter
writer(&stream
);
163 writeKeyValueTreeAsMdp(&writer
, builder
.build());
167 EXPECT_EQ(stream
.toString(),
169 "\n; Density guided simulation\ndensity-guided-simulation-active = false\n"));
172 TEST_F(DensityFittingOptionsTest
, OutputDefaultValuesWhenActive
)
174 setFromMdpValues(densityFittingSetActiveAsMdpValues());
175 // Transform module data into a flat key-value tree for output.
177 StringOutputStream stream
;
178 KeyValueTreeBuilder builder
;
179 KeyValueTreeObjectBuilder builderObject
= builder
.rootObject();
181 densityFittingOptions_
.buildMdpOutput(&builderObject
);
183 TextWriter
writer(&stream
);
184 writeKeyValueTreeAsMdp(&writer
, builder
.build());
187 std::string expectedString
= {
189 "; Density guided simulation\n"
190 "density-guided-simulation-active = true\n"
191 "density-guided-simulation-group = protein\n"
192 "; Similarity measure between densities: inner-product, relative-entropy, or "
193 "cross-correlation\n"
194 "density-guided-simulation-similarity-measure = inner-product\n"
195 "; Atom amplitude for spreading onto grid: unity, mass, or charge\n"
196 "density-guided-simulation-atom-spreading-weight = unity\n"
197 "density-guided-simulation-force-constant = 1e+09\n"
198 "density-guided-simulation-gaussian-transform-spreading-width = 0.2\n"
199 "density-guided-simulation-gaussian-transform-spreading-range-in-multiples-of-width = 4\n"
200 "; Reference density file location as absolute path or relative to the gmx mdrun calling "
202 "density-guided-simulation-reference-density-filename = reference.mrc\n"
203 "density-guided-simulation-nst = 1\n"
204 "; Normalize the sum of density voxel values to one\n"
205 "density-guided-simulation-normalize-densities = true\n"
206 "; Apply adaptive force scaling\n"
207 "density-guided-simulation-adaptive-force-scaling = false\n"
208 "; Time constant for adaptive force scaling in ps\n"
209 "density-guided-simulation-adaptive-force-scaling-time-constant = 4\n"
212 EXPECT_EQ(expectedString
, stream
.toString());
216 TEST_F(DensityFittingOptionsTest
, CanConvertGroupStringToIndexGroup
)
218 setFromMdpValues(densityFittingSetActiveAsMdpValues());
220 const auto indexGroupAndNames
= genericIndexGroupsAndNames();
221 densityFittingOptions_
.setFitGroupIndices(indexGroupAndNames
);
223 EXPECT_EQ(1, densityFittingOptions_
.buildParameters().indices_
.size());
224 EXPECT_EQ(1, densityFittingOptions_
.buildParameters().indices_
[0]);
227 TEST_F(DensityFittingOptionsTest
, InternalsToKvt
)
229 // stores the default internal options
230 DensityFittingOptions densityFittingOptions
;
231 KeyValueTreeBuilder builder
;
232 densityFittingOptions
.writeInternalParametersToKvt(builder
.rootObject());
233 const auto kvtTree
= builder
.build();
234 EXPECT_TRUE(kvtTree
.keyExists("density-guided-simulation-group"));
235 EXPECT_TRUE(kvtTree
["density-guided-simulation-group"].isArray());
236 auto storedIndex
= kvtTree
["density-guided-simulation-group"].asArray().values();
238 EXPECT_EQ(0, storedIndex
.size());
241 TEST_F(DensityFittingOptionsTest
, KvtToInternal
)
243 setFromMdpValues(densityFittingSetActiveAsMdpValues());
245 KeyValueTreeBuilder builder
;
247 builder
.rootObject().addUniformArray
<std::int64_t>("density-guided-simulation-group");
248 addedArray
.addValue(1);
249 addedArray
.addValue(15);
250 const auto tree
= builder
.build();
252 densityFittingOptions_
.readInternalParametersFromKvt(tree
);
254 EXPECT_EQ(2, densityFittingOptions_
.buildParameters().indices_
.size());
255 EXPECT_EQ(1, densityFittingOptions_
.buildParameters().indices_
[0]);
256 EXPECT_EQ(15, densityFittingOptions_
.buildParameters().indices_
[1]);
259 TEST_F(DensityFittingOptionsTest
, RoundTripForInternalsIsIdempotent
)
261 setFromMdpValues(densityFittingSetActiveAsMdpValues());
263 const IndexGroupsAndNames indexGroupAndNames
= genericIndexGroupsAndNames();
264 densityFittingOptions_
.setFitGroupIndices(indexGroupAndNames
);
267 DensityFittingParameters parametersBefore
= densityFittingOptions_
.buildParameters();
269 KeyValueTreeBuilder builder
;
270 densityFittingOptions_
.writeInternalParametersToKvt(builder
.rootObject());
271 const auto inputTree
= builder
.build();
273 mangleInternalParameters();
275 DensityFittingParameters parametersAfter
= densityFittingOptions_
.buildParameters();
276 EXPECT_NE(parametersBefore
, parametersAfter
);
278 densityFittingOptions_
.readInternalParametersFromKvt(inputTree
);
280 parametersAfter
= densityFittingOptions_
.buildParameters();
281 EXPECT_EQ(parametersBefore
, parametersAfter
);