5a1a01d1beb67d9336228d67106ce3f20f629a6b
[alexxy/gromacs.git] / src / gromacs / applied_forces / densityfittingoptions.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019, 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  * Implements force provider for density fitting
38  *
39  * \author Christian Blau <blau@kth.se>
40  * \ingroup module_applied_forces
41  */
42 #include "gmxpre.h"
43
44 #include "densityfittingoptions.h"
45
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/exceptions.h"
52 #include "gromacs/utility/keyvaluetreebuilder.h"
53 #include "gromacs/utility/keyvaluetreetransform.h"
54 #include "gromacs/utility/mdmodulenotification.h"
55 #include "gromacs/utility/strconvert.h"
56
57 #include "densityfittingamplitudelookup.h"
58
59 namespace gmx
60 {
61
62 namespace
63 {
64
65 /*! \brief Helper to declare mdp transform rules.
66  *
67  * Enforces uniform mdp options that are always prepended with the correct
68  * string for the densityfitting mdp options.
69  *
70  * \tparam ToType type to be transformed to
71  * \tparam TransformWithFunctionType type of transformation function to be used
72  *
73  * \param[in] rules KVT transformation rules
74  * \param[in] transformationFunction the function to transform the flat kvt tree
75  * \param[in] optionTag string tag that describes the mdp option, appended to the
76  *                      default string for the density guided simulation
77  */
78 template<class ToType, class TransformWithFunctionType>
79 void densityfittingMdpTransformFromString(IKeyValueTreeTransformRules* rules,
80                                           TransformWithFunctionType    transformationFunction,
81                                           const std::string&           optionTag)
82 {
83     rules->addRule()
84             .from<std::string>("/" + DensityFittingModuleInfo::name_ + "-" + optionTag)
85             .to<ToType>("/" + DensityFittingModuleInfo::name_ + "/" + optionTag)
86             .transformWith(transformationFunction);
87 }
88 /*! \brief Helper to declare mdp output.
89  *
90  * Enforces uniform mdp options output strings that are always prepended with the
91  * correct string for the densityfitting mdp options and are consistent with the
92  * options name and transformation type.
93  *
94  * \tparam OptionType the type of the mdp option
95  * \param[in] builder the KVT builder to generate the output
96  * \param[in] option the mdp option
97  * \param[in] optionTag string tag that describes the mdp option, appended to the
98  *                      default string for the density guided simulation
99  */
100 template<class OptionType>
101 void addDensityFittingMdpOutputValue(KeyValueTreeObjectBuilder* builder,
102                                      const OptionType&          option,
103                                      const std::string&         optionTag)
104 {
105     builder->addValue<OptionType>(DensityFittingModuleInfo::name_ + "-" + optionTag, option);
106 }
107
108 /*! \brief Helper to declare mdp output comments.
109  *
110  * Enforces uniform mdp options comment output strings that are always prepended
111  * with the correct string for the densityfitting mdp options and are consistent
112  * with the options name and transformation type.
113  *
114  * \param[in] builder the KVT builder to generate the output
115  * \param[in] comment on the mdp option
116  * \param[in] optionTag string tag that describes the mdp option
117  */
118 void addDensityFittingMdpOutputValueComment(KeyValueTreeObjectBuilder* builder,
119                                             const std::string&         comment,
120                                             const std::string&         optionTag)
121 {
122     builder->addValue<std::string>("comment-" + DensityFittingModuleInfo::name_ + "-" + optionTag, comment);
123 }
124
125 } // namespace
126
127 void DensityFittingOptions::initMdpTransform(IKeyValueTreeTransformRules* rules)
128 {
129     const auto& stringIdentityTransform = [](std::string s) { return s; };
130     densityfittingMdpTransformFromString<bool>(rules, &fromStdString<bool>, c_activeTag_);
131     densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_groupTag_);
132     densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform,
133                                                       c_similarityMeasureTag_);
134     densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_amplitudeMethodTag_);
135     densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>, c_forceConstantTag_);
136     densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>,
137                                                c_gaussianTransformSpreadingWidthTag_);
138     densityfittingMdpTransformFromString<real>(
139             rules, &fromStdString<real>, c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_);
140     densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform,
141                                                       c_referenceDensityFileNameTag_);
142     densityfittingMdpTransformFromString<std::int64_t>(rules, &fromStdString<std::int64_t>,
143                                                        c_everyNStepsTag_);
144     densityfittingMdpTransformFromString<bool>(rules, &fromStdString<bool>, c_normalizeDensitiesTag_);
145     densityfittingMdpTransformFromString<bool>(rules, &fromStdString<bool>, c_adaptiveForceScalingTag_);
146     densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>,
147                                                c_adaptiveForceScalingTimeConstantTag_);
148 }
149
150 void DensityFittingOptions::buildMdpOutput(KeyValueTreeObjectBuilder* builder) const
151 {
152
153     addDensityFittingMdpOutputValueComment(builder, "", "empty-line");
154     addDensityFittingMdpOutputValueComment(builder, "; Density guided simulation", "module");
155
156     addDensityFittingMdpOutputValue(builder, parameters_.active_, c_activeTag_);
157     if (parameters_.active_)
158     {
159         addDensityFittingMdpOutputValue(builder, groupString_, c_groupTag_);
160
161         addDensityFittingMdpOutputValueComment(
162                 builder,
163                 "; Similarity measure between densities: inner-product, relative-entropy, or "
164                 "cross-correlation",
165                 c_similarityMeasureTag_);
166         addDensityFittingMdpOutputValue<std::string>(
167                 builder, c_densitySimilarityMeasureMethodNames[parameters_.similarityMeasureMethod_],
168                 c_similarityMeasureTag_);
169
170         addDensityFittingMdpOutputValueComment(
171                 builder, "; Atom amplitude for spreading onto grid: unity, mass, or charge",
172                 c_amplitudeMethodTag_);
173         addDensityFittingMdpOutputValue<std::string>(
174                 builder, c_densityFittingAmplitudeMethodNames[parameters_.amplitudeLookupMethod_],
175                 c_amplitudeMethodTag_);
176
177         addDensityFittingMdpOutputValue(builder, parameters_.forceConstant_, c_forceConstantTag_);
178         addDensityFittingMdpOutputValue(builder, parameters_.gaussianTransformSpreadingWidth_,
179                                         c_gaussianTransformSpreadingWidthTag_);
180         addDensityFittingMdpOutputValue(builder, parameters_.gaussianTransformSpreadingRangeInMultiplesOfWidth_,
181                                         c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_);
182         addDensityFittingMdpOutputValueComment(builder,
183                                                "; Reference density file location as absolute path "
184                                                "or relative to the gmx mdrun calling location",
185                                                c_referenceDensityFileNameTag_);
186         addDensityFittingMdpOutputValue(builder, referenceDensityFileName_, c_referenceDensityFileNameTag_);
187         addDensityFittingMdpOutputValue(builder, parameters_.calculationIntervalInSteps_, c_everyNStepsTag_);
188         addDensityFittingMdpOutputValueComment(
189                 builder, "; Normalize the sum of density voxel values to one", c_normalizeDensitiesTag_);
190         addDensityFittingMdpOutputValue(builder, parameters_.normalizeDensities_, c_normalizeDensitiesTag_);
191         addDensityFittingMdpOutputValueComment(builder, "; Apply adaptive force scaling",
192                                                c_adaptiveForceScalingTag_);
193         addDensityFittingMdpOutputValue(builder, parameters_.adaptiveForceScaling_,
194                                         c_adaptiveForceScalingTag_);
195         addDensityFittingMdpOutputValueComment(builder,
196                                                "; Time constant for adaptive force scaling in ps",
197                                                c_adaptiveForceScalingTimeConstantTag_);
198         addDensityFittingMdpOutputValue(builder, parameters_.adaptiveForceScalingTimeConstant_,
199                                         c_adaptiveForceScalingTimeConstantTag_);
200     }
201 }
202
203 void DensityFittingOptions::initMdpOptions(IOptionsContainerWithSections* options)
204 {
205     auto section = options->addSection(OptionSection(DensityFittingModuleInfo::name_.c_str()));
206
207     section.addOption(BooleanOption(c_activeTag_.c_str()).store(&parameters_.active_));
208     section.addOption(StringOption(c_groupTag_.c_str()).store(&groupString_));
209
210     section.addOption(EnumOption<DensitySimilarityMeasureMethod>(c_similarityMeasureTag_.c_str())
211                               .enumValue(c_densitySimilarityMeasureMethodNames.m_elements)
212                               .store(&parameters_.similarityMeasureMethod_));
213
214     section.addOption(EnumOption<DensityFittingAmplitudeMethod>(c_amplitudeMethodTag_.c_str())
215                               .enumValue(c_densityFittingAmplitudeMethodNames.m_elements)
216                               .store(&parameters_.amplitudeLookupMethod_));
217
218     section.addOption(RealOption(c_forceConstantTag_.c_str()).store(&parameters_.forceConstant_));
219     section.addOption(RealOption(c_gaussianTransformSpreadingWidthTag_.c_str())
220                               .store(&parameters_.gaussianTransformSpreadingWidth_));
221     section.addOption(RealOption(c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_.c_str())
222                               .store(&parameters_.gaussianTransformSpreadingRangeInMultiplesOfWidth_));
223     section.addOption(StringOption(c_referenceDensityFileNameTag_.c_str()).store(&referenceDensityFileName_));
224     section.addOption(Int64Option(c_everyNStepsTag_.c_str()).store(&parameters_.calculationIntervalInSteps_));
225     section.addOption(BooleanOption(c_normalizeDensitiesTag_.c_str()).store(&parameters_.normalizeDensities_));
226     section.addOption(
227             BooleanOption(c_adaptiveForceScalingTag_.c_str()).store(&parameters_.adaptiveForceScaling_));
228     section.addOption(RealOption(c_adaptiveForceScalingTimeConstantTag_.c_str())
229                               .store(&parameters_.adaptiveForceScalingTimeConstant_));
230 }
231
232 bool DensityFittingOptions::active() const
233 {
234     return parameters_.active_;
235 }
236
237 const DensityFittingParameters& DensityFittingOptions::buildParameters()
238 {
239     // the options modules does not know unsigned integers so any input of this
240     // kind is rectified here
241     if (parameters_.calculationIntervalInSteps_ < 1)
242     {
243         parameters_.calculationIntervalInSteps_ = 1;
244     }
245     return parameters_;
246 }
247
248 void DensityFittingOptions::setFitGroupIndices(const IndexGroupsAndNames& indexGroupsAndNames)
249 {
250     if (!parameters_.active_)
251     {
252         return;
253     }
254     parameters_.indices_ = indexGroupsAndNames.indices(groupString_);
255 }
256
257 void DensityFittingOptions::writeInternalParametersToKvt(KeyValueTreeObjectBuilder treeBuilder)
258 {
259     auto groupIndexAdder = treeBuilder.addUniformArray<std::int64_t>(DensityFittingModuleInfo::name_
260                                                                      + "-" + c_groupTag_);
261     for (const auto& indexValue : parameters_.indices_)
262     {
263         groupIndexAdder.addValue(indexValue);
264     }
265 }
266
267 void DensityFittingOptions::readInternalParametersFromKvt(const KeyValueTreeObject& tree)
268 {
269     if (!parameters_.active_)
270     {
271         return;
272     }
273
274     if (!tree.keyExists(DensityFittingModuleInfo::name_ + "-" + c_groupTag_))
275     {
276         GMX_THROW(InconsistentInputError(
277                 "Cannot find atom index vector required for density guided simulation."));
278     }
279     auto kvtIndexArray = tree[DensityFittingModuleInfo::name_ + "-" + c_groupTag_].asArray().values();
280     parameters_.indices_.resize(kvtIndexArray.size());
281     std::transform(std::begin(kvtIndexArray), std::end(kvtIndexArray), std::begin(parameters_.indices_),
282                    [](const KeyValueTreeValue& val) { return val.cast<std::int64_t>(); });
283 }
284
285 void DensityFittingOptions::checkEnergyCaluclationFrequency(
286         EnergyCalculationFrequencyErrors* energyCalculationFrequencyErrors) const
287 {
288     if (energyCalculationFrequencyErrors->energyCalculationIntervalInSteps() % parameters_.calculationIntervalInSteps_
289         != 0)
290     {
291         energyCalculationFrequencyErrors->addError(
292                 "nstcalcenergy ("
293                 + toString(energyCalculationFrequencyErrors->energyCalculationIntervalInSteps())
294                 + ") is not a multiple of " + DensityFittingModuleInfo::name_ + "-" + c_everyNStepsTag_
295                 + " (" + toString(parameters_.calculationIntervalInSteps_) + ") .");
296     }
297 }
298
299 const std::string& DensityFittingOptions::referenceDensityFileName() const
300 {
301     return referenceDensityFileName_;
302 }
303 } // namespace gmx