Fix unsigned integer overflow bug in spline tables
[gromacs.git] / src / gromacs / applied_forces / densityfittingoptions.cpp
blob29b69cbbeced09d4fd6c42324c8a7247b8b4fd8d
1 /*
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.
35 /*! \internal \file
36 * \brief
37 * Implements force provider for density fitting
39 * \author Christian Blau <blau@kth.se>
40 * \ingroup module_applied_forces
42 #include "gmxpre.h"
44 #include "densityfittingoptions.h"
46 #include "gromacs/applied_forces/densityfitting.h"
47 #include "gromacs/math/densityfit.h"
48 #include "gromacs/options/basicoptions.h"
49 #include "gromacs/options/optionsection.h"
50 #include "gromacs/selection/indexutil.h"
51 #include "gromacs/utility/enumerationhelpers.h"
52 #include "gromacs/utility/exceptions.h"
53 #include "gromacs/utility/keyvaluetreebuilder.h"
54 #include "gromacs/utility/keyvaluetreetransform.h"
55 #include "gromacs/utility/mdmodulenotification.h"
56 #include "gromacs/utility/strconvert.h"
58 #include "densityfittingamplitudelookup.h"
60 namespace gmx
63 namespace
66 /*! \brief Helper to declare mdp transform rules.
68 * Enforces uniform mdp options that are always prepended with the correct
69 * string for the densityfitting mdp options.
71 * \tparam ToType type to be transformed to
72 * \tparam TransformWithFunctionType type of transformation function to be used
74 * \param[in] rules KVT transformation rules
75 * \param[in] transformationFunction the function to transform the flat kvt tree
76 * \param[in] optionTag string tag that describes the mdp option, appended to the
77 * default string for the density guided simulation
79 template<class ToType, class TransformWithFunctionType>
80 void densityfittingMdpTransformFromString(IKeyValueTreeTransformRules* rules,
81 TransformWithFunctionType transformationFunction,
82 const std::string& optionTag)
84 rules->addRule()
85 .from<std::string>("/" + DensityFittingModuleInfo::name_ + "-" + optionTag)
86 .to<ToType>("/" + DensityFittingModuleInfo::name_ + "/" + optionTag)
87 .transformWith(transformationFunction);
89 /*! \brief Helper to declare mdp output.
91 * Enforces uniform mdp options output strings that are always prepended with the
92 * correct string for the densityfitting mdp options and are consistent with the
93 * options name and transformation type.
95 * \tparam OptionType the type of the mdp option
96 * \param[in] builder the KVT builder to generate the output
97 * \param[in] option the mdp option
98 * \param[in] optionTag string tag that describes the mdp option, appended to the
99 * default string for the density guided simulation
101 template<class OptionType>
102 void addDensityFittingMdpOutputValue(KeyValueTreeObjectBuilder* builder,
103 const OptionType& option,
104 const std::string& optionTag)
106 builder->addValue<OptionType>(DensityFittingModuleInfo::name_ + "-" + optionTag, option);
109 /*! \brief Helper to declare mdp output comments.
111 * Enforces uniform mdp options comment output strings that are always prepended
112 * with the correct string for the densityfitting mdp options and are consistent
113 * with the options name and transformation type.
115 * \param[in] builder the KVT builder to generate the output
116 * \param[in] comment on the mdp option
117 * \param[in] optionTag string tag that describes the mdp option
119 void addDensityFittingMdpOutputValueComment(KeyValueTreeObjectBuilder* builder,
120 const std::string& comment,
121 const std::string& optionTag)
123 builder->addValue<std::string>("comment-" + DensityFittingModuleInfo::name_ + "-" + optionTag, comment);
126 } // namespace
128 void DensityFittingOptions::initMdpTransform(IKeyValueTreeTransformRules* rules)
130 const auto& stringIdentityTransform = [](std::string s) { return s; };
131 densityfittingMdpTransformFromString<bool>(rules, &fromStdString<bool>, c_activeTag_);
132 densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_groupTag_);
133 densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform,
134 c_similarityMeasureTag_);
135 densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_amplitudeMethodTag_);
136 densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>, c_forceConstantTag_);
137 densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>,
138 c_gaussianTransformSpreadingWidthTag_);
139 densityfittingMdpTransformFromString<real>(
140 rules, &fromStdString<real>, c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_);
141 densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform,
142 c_referenceDensityFileNameTag_);
143 densityfittingMdpTransformFromString<std::int64_t>(rules, &fromStdString<std::int64_t>,
144 c_everyNStepsTag_);
145 densityfittingMdpTransformFromString<bool>(rules, &fromStdString<bool>, c_normalizeDensitiesTag_);
146 densityfittingMdpTransformFromString<bool>(rules, &fromStdString<bool>, c_adaptiveForceScalingTag_);
147 densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>,
148 c_adaptiveForceScalingTimeConstantTag_);
151 //! Name the methods that may be used to evaluate similarity between densities
152 static const EnumerationArray<DensitySimilarityMeasureMethod, const char*> c_densitySimilarityMeasureMethodNames = {
153 { "inner-product", "relative-entropy", "cross-correlation" }
155 //! The names of the methods to determine the amplitude of the atoms to be spread on a grid
156 static const EnumerationArray<DensityFittingAmplitudeMethod, const char*> c_densityFittingAmplitudeMethodNames = {
157 { "unity", "mass", "charge" }
160 void DensityFittingOptions::buildMdpOutput(KeyValueTreeObjectBuilder* builder) const
163 addDensityFittingMdpOutputValueComment(builder, "", "empty-line");
164 addDensityFittingMdpOutputValueComment(builder, "; Density guided simulation", "module");
166 addDensityFittingMdpOutputValue(builder, parameters_.active_, c_activeTag_);
167 if (parameters_.active_)
169 addDensityFittingMdpOutputValue(builder, groupString_, c_groupTag_);
171 addDensityFittingMdpOutputValueComment(
172 builder,
173 "; Similarity measure between densities: inner-product, relative-entropy, or "
174 "cross-correlation",
175 c_similarityMeasureTag_);
176 addDensityFittingMdpOutputValue<std::string>(
177 builder, c_densitySimilarityMeasureMethodNames[parameters_.similarityMeasureMethod_],
178 c_similarityMeasureTag_);
180 addDensityFittingMdpOutputValueComment(
181 builder, "; Atom amplitude for spreading onto grid: unity, mass, or charge",
182 c_amplitudeMethodTag_);
183 addDensityFittingMdpOutputValue<std::string>(
184 builder, c_densityFittingAmplitudeMethodNames[parameters_.amplitudeLookupMethod_],
185 c_amplitudeMethodTag_);
187 addDensityFittingMdpOutputValue(builder, parameters_.forceConstant_, c_forceConstantTag_);
188 addDensityFittingMdpOutputValue(builder, parameters_.gaussianTransformSpreadingWidth_,
189 c_gaussianTransformSpreadingWidthTag_);
190 addDensityFittingMdpOutputValue(builder, parameters_.gaussianTransformSpreadingRangeInMultiplesOfWidth_,
191 c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_);
192 addDensityFittingMdpOutputValueComment(builder,
193 "; Reference density file location as absolute path "
194 "or relative to the gmx mdrun calling location",
195 c_referenceDensityFileNameTag_);
196 addDensityFittingMdpOutputValue(builder, referenceDensityFileName_, c_referenceDensityFileNameTag_);
197 addDensityFittingMdpOutputValue(builder, parameters_.calculationIntervalInSteps_, c_everyNStepsTag_);
198 addDensityFittingMdpOutputValueComment(
199 builder, "; Normalize the sum of density voxel values to one", c_normalizeDensitiesTag_);
200 addDensityFittingMdpOutputValue(builder, parameters_.normalizeDensities_, c_normalizeDensitiesTag_);
201 addDensityFittingMdpOutputValueComment(builder, "; Apply adaptive force scaling",
202 c_adaptiveForceScalingTag_);
203 addDensityFittingMdpOutputValue(builder, parameters_.adaptiveForceScaling_,
204 c_adaptiveForceScalingTag_);
205 addDensityFittingMdpOutputValueComment(builder,
206 "; Time constant for adaptive force scaling in ps",
207 c_adaptiveForceScalingTimeConstantTag_);
208 addDensityFittingMdpOutputValue(builder, parameters_.adaptiveForceScalingTimeConstant_,
209 c_adaptiveForceScalingTimeConstantTag_);
213 void DensityFittingOptions::initMdpOptions(IOptionsContainerWithSections* options)
215 auto section = options->addSection(OptionSection(DensityFittingModuleInfo::name_.c_str()));
217 section.addOption(BooleanOption(c_activeTag_.c_str()).store(&parameters_.active_));
218 section.addOption(StringOption(c_groupTag_.c_str()).store(&groupString_));
220 section.addOption(EnumOption<DensitySimilarityMeasureMethod>(c_similarityMeasureTag_.c_str())
221 .enumValue(c_densitySimilarityMeasureMethodNames)
222 .store(&parameters_.similarityMeasureMethod_));
224 section.addOption(EnumOption<DensityFittingAmplitudeMethod>(c_amplitudeMethodTag_.c_str())
225 .enumValue(c_densityFittingAmplitudeMethodNames)
226 .store(&parameters_.amplitudeLookupMethod_));
228 section.addOption(RealOption(c_forceConstantTag_.c_str()).store(&parameters_.forceConstant_));
229 section.addOption(RealOption(c_gaussianTransformSpreadingWidthTag_.c_str())
230 .store(&parameters_.gaussianTransformSpreadingWidth_));
231 section.addOption(RealOption(c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_.c_str())
232 .store(&parameters_.gaussianTransformSpreadingRangeInMultiplesOfWidth_));
233 section.addOption(StringOption(c_referenceDensityFileNameTag_.c_str()).store(&referenceDensityFileName_));
234 section.addOption(Int64Option(c_everyNStepsTag_.c_str()).store(&parameters_.calculationIntervalInSteps_));
235 section.addOption(BooleanOption(c_normalizeDensitiesTag_.c_str()).store(&parameters_.normalizeDensities_));
236 section.addOption(
237 BooleanOption(c_adaptiveForceScalingTag_.c_str()).store(&parameters_.adaptiveForceScaling_));
238 section.addOption(RealOption(c_adaptiveForceScalingTimeConstantTag_.c_str())
239 .store(&parameters_.adaptiveForceScalingTimeConstant_));
242 bool DensityFittingOptions::active() const
244 return parameters_.active_;
247 const DensityFittingParameters& DensityFittingOptions::buildParameters()
249 // the options modules does not know unsigned integers so any input of this
250 // kind is rectified here
251 if (parameters_.calculationIntervalInSteps_ < 1)
253 parameters_.calculationIntervalInSteps_ = 1;
255 return parameters_;
258 void DensityFittingOptions::setFitGroupIndices(const IndexGroupsAndNames& indexGroupsAndNames)
260 if (!parameters_.active_)
262 return;
264 parameters_.indices_ = indexGroupsAndNames.indices(groupString_);
267 void DensityFittingOptions::writeInternalParametersToKvt(KeyValueTreeObjectBuilder treeBuilder)
269 auto groupIndexAdder = treeBuilder.addUniformArray<std::int64_t>(DensityFittingModuleInfo::name_
270 + "-" + c_groupTag_);
271 for (const auto& indexValue : parameters_.indices_)
273 groupIndexAdder.addValue(indexValue);
277 void DensityFittingOptions::readInternalParametersFromKvt(const KeyValueTreeObject& tree)
279 if (!parameters_.active_)
281 return;
284 if (!tree.keyExists(DensityFittingModuleInfo::name_ + "-" + c_groupTag_))
286 GMX_THROW(InconsistentInputError(
287 "Cannot find atom index vector required for density guided simulation."));
289 auto kvtIndexArray = tree[DensityFittingModuleInfo::name_ + "-" + c_groupTag_].asArray().values();
290 parameters_.indices_.resize(kvtIndexArray.size());
291 std::transform(std::begin(kvtIndexArray), std::end(kvtIndexArray), std::begin(parameters_.indices_),
292 [](const KeyValueTreeValue& val) { return val.cast<std::int64_t>(); });
295 void DensityFittingOptions::checkEnergyCaluclationFrequency(
296 EnergyCalculationFrequencyErrors* energyCalculationFrequencyErrors) const
298 if (energyCalculationFrequencyErrors->energyCalculationIntervalInSteps() % parameters_.calculationIntervalInSteps_
299 != 0)
301 energyCalculationFrequencyErrors->addError(
302 "nstcalcenergy ("
303 + toString(energyCalculationFrequencyErrors->energyCalculationIntervalInSteps())
304 + ") is not a multiple of " + DensityFittingModuleInfo::name_ + "-" + c_everyNStepsTag_
305 + " (" + toString(parameters_.calculationIntervalInSteps_) + ") .");
309 const std::string& DensityFittingOptions::referenceDensityFileName() const
311 return referenceDensityFileName_;
313 } // namespace gmx