Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / src / gromacs / applied_forces / densityfitting / tests / densityfittingoptions.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 /*! \internal \file
36  * \brief
37  * Tests for density fitting module options.
38  *
39  * \author Christian Blau <blau@kth.se>
40  * \ingroup module_applied_forces
41  */
42 #include "gmxpre.h"
43
44 #include "gromacs/applied_forces/densityfitting/densityfittingoptions.h"
45
46 #include <string>
47 #include <vector>
48
49 #include <gtest/gtest.h>
50
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"
61
62 #include "testutils/testasserts.h"
63 #include "testutils/testmatchers.h"
64
65 namespace gmx
66 {
67
68 namespace
69 {
70
71 class DensityFittingOptionsTest : public ::testing::Test
72 {
73 public:
74     DensityFittingOptionsTest() { init_blocka(&defaultGroups_); }
75     ~DensityFittingOptionsTest() override { done_blocka(&defaultGroups_); }
76
77     void setFromMdpValues(const KeyValueTreeObject& densityFittingMdpValues)
78     {
79         // set up options
80         Options densityFittingModuleOptions;
81         densityFittingOptions_.initMdpOptions(&densityFittingModuleOptions);
82
83         // Add rules to transform mdp inputs to densityFittingModule data
84         KeyValueTreeTransformer transform;
85         transform.rules()->addRule().keyMatchType("/", StringCompareType::CaseAndDashInsensitive);
86
87         densityFittingOptions_.initMdpTransform(transform.rules());
88
89         // Execute the transform on the mdpValues
90         auto transformedMdpValues = transform.transform(densityFittingMdpValues, nullptr);
91         assignOptionsFromKeyValueTree(&densityFittingModuleOptions, transformedMdpValues.object(), nullptr);
92     }
93
94     static KeyValueTreeObject densityFittingSetActiveAsMdpValues()
95     {
96         // Prepare MDP inputs
97         KeyValueTreeBuilder mdpValueBuilder;
98         mdpValueBuilder.rootObject().addValue("density-guided-simulation-active", std::string("yes"));
99         return mdpValueBuilder.build();
100     }
101
102     IndexGroupsAndNames genericIndexGroupsAndNames()
103     {
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 };
111     }
112
113     IndexGroupsAndNames differingIndexGroupsAndNames()
114     {
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 };
122     }
123
124     void mangleInternalParameters()
125     {
126         densityFittingOptions_.setFitGroupIndices(differingIndexGroupsAndNames());
127     }
128
129 protected:
130     t_blocka              defaultGroups_;
131     DensityFittingOptions densityFittingOptions_;
132 };
133
134 TEST_F(DensityFittingOptionsTest, DefaultParameters)
135 {
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_);
144 }
145
146 TEST_F(DensityFittingOptionsTest, OptionSetsActive)
147 {
148     EXPECT_FALSE(densityFittingOptions_.buildParameters().active_);
149     setFromMdpValues(densityFittingSetActiveAsMdpValues());
150     EXPECT_TRUE(densityFittingOptions_.buildParameters().active_);
151 }
152
153 TEST_F(DensityFittingOptionsTest, OutputNoDefaultValuesWhenInactive)
154 {
155     // Transform module data into a flat key-value tree for output.
156
157     StringOutputStream        stream;
158     KeyValueTreeBuilder       builder;
159     KeyValueTreeObjectBuilder builderObject = builder.rootObject();
160
161     densityFittingOptions_.buildMdpOutput(&builderObject);
162     {
163         TextWriter writer(&stream);
164         writeKeyValueTreeAsMdp(&writer, builder.build());
165     }
166     stream.close();
167
168     EXPECT_EQ(stream.toString(),
169               std::string(
170                       "\n; Density guided simulation\ndensity-guided-simulation-active = false\n"));
171 }
172
173 TEST_F(DensityFittingOptionsTest, OutputDefaultValuesWhenActive)
174 {
175     setFromMdpValues(densityFittingSetActiveAsMdpValues());
176     // Transform module data into a flat key-value tree for output.
177
178     StringOutputStream        stream;
179     KeyValueTreeBuilder       builder;
180     KeyValueTreeObjectBuilder builderObject = builder.rootObject();
181
182     densityFittingOptions_.buildMdpOutput(&builderObject);
183     {
184         TextWriter writer(&stream);
185         writeKeyValueTreeAsMdp(&writer, builder.build());
186     }
187     stream.close();
188     std::string expectedString = {
189         "\n"
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 "
202         "location\n"
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"
211     };
212
213     EXPECT_EQ(expectedString, stream.toString());
214 }
215
216
217 TEST_F(DensityFittingOptionsTest, CanConvertGroupStringToIndexGroup)
218 {
219     setFromMdpValues(densityFittingSetActiveAsMdpValues());
220
221     const auto indexGroupAndNames = genericIndexGroupsAndNames();
222     densityFittingOptions_.setFitGroupIndices(indexGroupAndNames);
223
224     EXPECT_EQ(1, densityFittingOptions_.buildParameters().indices_.size());
225     EXPECT_EQ(1, densityFittingOptions_.buildParameters().indices_[0]);
226 }
227
228 TEST_F(DensityFittingOptionsTest, InternalsToKvt)
229 {
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();
238
239     EXPECT_EQ(0, storedIndex.size());
240 }
241
242 TEST_F(DensityFittingOptionsTest, KvtToInternal)
243 {
244     setFromMdpValues(densityFittingSetActiveAsMdpValues());
245
246     KeyValueTreeBuilder builder;
247     auto                addedArray =
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();
252
253     densityFittingOptions_.readInternalParametersFromKvt(tree);
254
255     EXPECT_EQ(2, densityFittingOptions_.buildParameters().indices_.size());
256     EXPECT_EQ(1, densityFittingOptions_.buildParameters().indices_[0]);
257     EXPECT_EQ(15, densityFittingOptions_.buildParameters().indices_[1]);
258 }
259
260 TEST_F(DensityFittingOptionsTest, RoundTripForInternalsIsIdempotent)
261 {
262     setFromMdpValues(densityFittingSetActiveAsMdpValues());
263     {
264         const IndexGroupsAndNames indexGroupAndNames = genericIndexGroupsAndNames();
265         densityFittingOptions_.setFitGroupIndices(indexGroupAndNames);
266     }
267
268     DensityFittingParameters parametersBefore = densityFittingOptions_.buildParameters();
269
270     KeyValueTreeBuilder builder;
271     densityFittingOptions_.writeInternalParametersToKvt(builder.rootObject());
272     const auto inputTree = builder.build();
273
274     mangleInternalParameters();
275
276     DensityFittingParameters parametersAfter = densityFittingOptions_.buildParameters();
277     EXPECT_NE(parametersBefore, parametersAfter);
278
279     densityFittingOptions_.readInternalParametersFromKvt(inputTree);
280
281     parametersAfter = densityFittingOptions_.buildParameters();
282     EXPECT_EQ(parametersBefore, parametersAfter);
283 }
284
285 } // namespace
286
287 } // namespace gmx