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 * Implements force provider for density fitting
39 * \author Christian Blau <blau@kth.se>
40 * \ingroup module_applied_forces
44 #include "densityfittingoptions.h"
46 #include "gromacs/applied_forces/densityfitting/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"
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)
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);
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>(
134 rules, stringIdentityTransform, c_similarityMeasureTag_);
135 densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_amplitudeMethodTag_);
136 densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>, c_forceConstantTag_);
137 densityfittingMdpTransformFromString<real>(
138 rules, &fromStdString<real>, c_gaussianTransformSpreadingWidthTag_);
139 densityfittingMdpTransformFromString<real>(
140 rules, &fromStdString<real>, c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_);
141 densityfittingMdpTransformFromString<std::string>(
142 rules, stringIdentityTransform, c_referenceDensityFileNameTag_);
143 densityfittingMdpTransformFromString<std::int64_t>(
144 rules, &fromStdString<std::int64_t>, c_everyNStepsTag_);
145 densityfittingMdpTransformFromString<bool>(rules, &fromStdString<bool>, c_normalizeDensitiesTag_);
146 densityfittingMdpTransformFromString<bool>(rules, &fromStdString<bool>, c_adaptiveForceScalingTag_);
147 densityfittingMdpTransformFromString<real>(
148 rules, &fromStdString<real>, c_adaptiveForceScalingTimeConstantTag_);
149 const auto& stringRVecToStringRVecWithCheck = [](const std::string& str) {
150 return stringIdentityTransformWithArrayCheck<real, 3>(
152 "Reading three real values as vector while parsing the .mdp input failed in "
153 + DensityFittingModuleInfo::name_ + ".");
155 densityfittingMdpTransformFromString<std::string>(
156 rules, stringRVecToStringRVecWithCheck, c_translationTag_);
158 const auto& stringMatrixToStringMatrixWithCheck = [](const std::string& str) {
159 return stringIdentityTransformWithArrayCheck<real, 9>(
161 "Reading nine real values as vector while parsing the .mdp input failed in "
162 + DensityFittingModuleInfo::name_ + ".");
164 densityfittingMdpTransformFromString<std::string>(
165 rules, stringMatrixToStringMatrixWithCheck, c_transformationMatrixTag_);
168 //! Name the methods that may be used to evaluate similarity between densities
169 static const EnumerationArray<DensitySimilarityMeasureMethod, const char*> c_densitySimilarityMeasureMethodNames = {
170 { "inner-product", "relative-entropy", "cross-correlation" }
172 //! The names of the methods to determine the amplitude of the atoms to be spread on a grid
173 static const EnumerationArray<DensityFittingAmplitudeMethod, const char*> c_densityFittingAmplitudeMethodNames = {
174 { "unity", "mass", "charge" }
177 void DensityFittingOptions::buildMdpOutput(KeyValueTreeObjectBuilder* builder) const
180 addDensityFittingMdpOutputValueComment(builder, "", "empty-line");
181 addDensityFittingMdpOutputValueComment(builder, "; Density guided simulation", "module");
183 addDensityFittingMdpOutputValue(builder, parameters_.active_, c_activeTag_);
184 if (parameters_.active_)
186 addDensityFittingMdpOutputValue(builder, groupString_, c_groupTag_);
188 addDensityFittingMdpOutputValueComment(
190 "; Similarity measure between densities: inner-product, relative-entropy, or "
192 c_similarityMeasureTag_);
193 addDensityFittingMdpOutputValue<std::string>(
195 c_densitySimilarityMeasureMethodNames[parameters_.similarityMeasureMethod_],
196 c_similarityMeasureTag_);
198 addDensityFittingMdpOutputValueComment(
199 builder, "; Atom amplitude for spreading onto grid: unity, mass, or charge", c_amplitudeMethodTag_);
200 addDensityFittingMdpOutputValue<std::string>(
202 c_densityFittingAmplitudeMethodNames[parameters_.amplitudeLookupMethod_],
203 c_amplitudeMethodTag_);
205 addDensityFittingMdpOutputValue(builder, parameters_.forceConstant_, c_forceConstantTag_);
206 addDensityFittingMdpOutputValue(builder,
207 parameters_.gaussianTransformSpreadingWidth_,
208 c_gaussianTransformSpreadingWidthTag_);
209 addDensityFittingMdpOutputValue(builder,
210 parameters_.gaussianTransformSpreadingRangeInMultiplesOfWidth_,
211 c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_);
212 addDensityFittingMdpOutputValueComment(builder,
213 "; Reference density file location as absolute path "
214 "or relative to the gmx mdrun calling location",
215 c_referenceDensityFileNameTag_);
216 addDensityFittingMdpOutputValue(builder, referenceDensityFileName_, c_referenceDensityFileNameTag_);
217 addDensityFittingMdpOutputValue(builder, parameters_.calculationIntervalInSteps_, c_everyNStepsTag_);
218 addDensityFittingMdpOutputValueComment(
219 builder, "; Normalize the sum of density voxel values to one", c_normalizeDensitiesTag_);
220 addDensityFittingMdpOutputValue(builder, parameters_.normalizeDensities_, c_normalizeDensitiesTag_);
221 addDensityFittingMdpOutputValueComment(
222 builder, "; Apply adaptive force scaling", c_adaptiveForceScalingTag_);
223 addDensityFittingMdpOutputValue(
224 builder, parameters_.adaptiveForceScaling_, c_adaptiveForceScalingTag_);
225 addDensityFittingMdpOutputValueComment(builder,
226 "; Time constant for adaptive force scaling in ps",
227 c_adaptiveForceScalingTimeConstantTag_);
228 addDensityFittingMdpOutputValue(builder,
229 parameters_.adaptiveForceScalingTimeConstant_,
230 c_adaptiveForceScalingTimeConstantTag_);
234 void DensityFittingOptions::initMdpOptions(IOptionsContainerWithSections* options)
236 auto section = options->addSection(OptionSection(DensityFittingModuleInfo::name_.c_str()));
238 section.addOption(BooleanOption(c_activeTag_.c_str()).store(¶meters_.active_));
239 section.addOption(StringOption(c_groupTag_.c_str()).store(&groupString_));
241 section.addOption(EnumOption<DensitySimilarityMeasureMethod>(c_similarityMeasureTag_.c_str())
242 .enumValue(c_densitySimilarityMeasureMethodNames)
243 .store(¶meters_.similarityMeasureMethod_));
245 section.addOption(EnumOption<DensityFittingAmplitudeMethod>(c_amplitudeMethodTag_.c_str())
246 .enumValue(c_densityFittingAmplitudeMethodNames)
247 .store(¶meters_.amplitudeLookupMethod_));
249 section.addOption(RealOption(c_forceConstantTag_.c_str()).store(¶meters_.forceConstant_));
250 section.addOption(RealOption(c_gaussianTransformSpreadingWidthTag_.c_str())
251 .store(¶meters_.gaussianTransformSpreadingWidth_));
252 section.addOption(RealOption(c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_.c_str())
253 .store(¶meters_.gaussianTransformSpreadingRangeInMultiplesOfWidth_));
254 section.addOption(StringOption(c_referenceDensityFileNameTag_.c_str()).store(&referenceDensityFileName_));
255 section.addOption(Int64Option(c_everyNStepsTag_.c_str()).store(¶meters_.calculationIntervalInSteps_));
256 section.addOption(BooleanOption(c_normalizeDensitiesTag_.c_str()).store(¶meters_.normalizeDensities_));
258 BooleanOption(c_adaptiveForceScalingTag_.c_str()).store(¶meters_.adaptiveForceScaling_));
259 section.addOption(RealOption(c_adaptiveForceScalingTimeConstantTag_.c_str())
260 .store(¶meters_.adaptiveForceScalingTimeConstant_));
261 section.addOption(StringOption(c_translationTag_.c_str()).store(¶meters_.translationString_));
263 StringOption(c_transformationMatrixTag_.c_str()).store(¶meters_.transformationMatrixString_));
266 bool DensityFittingOptions::active() const
268 return parameters_.active_;
271 const DensityFittingParameters& DensityFittingOptions::buildParameters()
273 // the options modules does not know unsigned integers so any input of this
274 // kind is rectified here
275 if (parameters_.calculationIntervalInSteps_ < 1)
277 parameters_.calculationIntervalInSteps_ = 1;
282 void DensityFittingOptions::setFitGroupIndices(const IndexGroupsAndNames& indexGroupsAndNames)
284 if (!parameters_.active_)
288 parameters_.indices_ = indexGroupsAndNames.indices(groupString_);
291 void DensityFittingOptions::writeInternalParametersToKvt(KeyValueTreeObjectBuilder treeBuilder)
293 auto groupIndexAdder = treeBuilder.addUniformArray<std::int64_t>(DensityFittingModuleInfo::name_
294 + "-" + c_groupTag_);
295 for (const auto& indexValue : parameters_.indices_)
297 groupIndexAdder.addValue(indexValue);
301 void DensityFittingOptions::readInternalParametersFromKvt(const KeyValueTreeObject& tree)
303 if (!parameters_.active_)
308 if (!tree.keyExists(DensityFittingModuleInfo::name_ + "-" + c_groupTag_))
310 GMX_THROW(InconsistentInputError(
311 "Cannot find atom index vector required for density guided simulation."));
313 auto kvtIndexArray = tree[DensityFittingModuleInfo::name_ + "-" + c_groupTag_].asArray().values();
314 parameters_.indices_.resize(kvtIndexArray.size());
315 std::transform(std::begin(kvtIndexArray),
316 std::end(kvtIndexArray),
317 std::begin(parameters_.indices_),
318 [](const KeyValueTreeValue& val) { return val.cast<std::int64_t>(); });
321 void DensityFittingOptions::checkEnergyCaluclationFrequency(
322 EnergyCalculationFrequencyErrors* energyCalculationFrequencyErrors) const
324 if (energyCalculationFrequencyErrors->energyCalculationIntervalInSteps() % parameters_.calculationIntervalInSteps_
327 energyCalculationFrequencyErrors->addError(
329 + toString(energyCalculationFrequencyErrors->energyCalculationIntervalInSteps())
330 + ") is not a multiple of " + DensityFittingModuleInfo::name_ + "-" + c_everyNStepsTag_
331 + " (" + toString(parameters_.calculationIntervalInSteps_) + ") .");
335 const std::string& DensityFittingOptions::referenceDensityFileName() const
337 return referenceDensityFileName_;