Grompp error for mismatching nst for energy calulation and densityfitting
[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,
106                                   option);
107 }
108
109 /*! \brief Helper to declare mdp output comments.
110  *
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.
114  *
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
118  */
119 void addDensityFittingMdpOutputValueComment(KeyValueTreeObjectBuilder *builder,
120                                             const std::string         &comment,
121                                             const std::string         &optionTag)
122 {
123     builder->addValue<std::string>("comment-" + DensityFittingModuleInfo::name_ + "-" + optionTag, comment);
124 }
125
126 }   // namespace
127
128 void DensityFittingOptions::initMdpTransform(IKeyValueTreeTransformRules * rules)
129 {
130     const auto &stringIdentityTransform = [](std::string s){
131             return s;
132         };
133     densityfittingMdpTransformFromString<bool>(rules, &fromStdString<bool>, c_activeTag_);
134     densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_groupTag_);
135     densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_similarityMeasureTag_);
136     densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_amplitudeMethodTag_);
137     densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>, c_forceConstantTag_);
138     densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>, c_gaussianTransformSpreadingWidthTag_);
139     densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>, c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_);
140     densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_referenceDensityFileNameTag_);
141     densityfittingMdpTransformFromString<std::int64_t>(rules, &fromStdString<std::int64_t>, c_everyNStepsTag_);
142     densityfittingMdpTransformFromString<bool>(rules, &fromStdString<bool>, c_normalizeDensitiesTag_);
143 }
144
145 void DensityFittingOptions::buildMdpOutput(KeyValueTreeObjectBuilder *builder) const
146 {
147
148     addDensityFittingMdpOutputValueComment(builder, "", "empty-line");
149     addDensityFittingMdpOutputValueComment(builder, "; Density guided simulation", "module");
150
151     addDensityFittingMdpOutputValue(builder, parameters_.active_, c_activeTag_);
152     if (parameters_.active_)
153     {
154         addDensityFittingMdpOutputValue(builder, groupString_, c_groupTag_);
155
156         addDensityFittingMdpOutputValueComment(builder, "; Similarity measure between densities: inner-product, relative-entropy, or cross-correlation", c_similarityMeasureTag_);
157         addDensityFittingMdpOutputValue<std::string>(builder,
158                                                      c_densitySimilarityMeasureMethodNames[parameters_.similarityMeasureMethod_],
159                                                      c_similarityMeasureTag_);
160
161         addDensityFittingMdpOutputValueComment(builder, "; Atom amplitude for spreading onto grid: unity, mass, or charge", c_amplitudeMethodTag_);
162         addDensityFittingMdpOutputValue<std::string>(builder,
163                                                      c_densityFittingAmplitudeMethodNames[parameters_.amplitudeLookupMethod_],
164                                                      c_amplitudeMethodTag_);
165
166         addDensityFittingMdpOutputValue(builder, parameters_.forceConstant_, c_forceConstantTag_);
167         addDensityFittingMdpOutputValue(builder, parameters_.gaussianTransformSpreadingWidth_, c_gaussianTransformSpreadingWidthTag_);
168         addDensityFittingMdpOutputValue(builder, parameters_.gaussianTransformSpreadingRangeInMultiplesOfWidth_, c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_);
169         addDensityFittingMdpOutputValueComment(builder, "; Reference density file location as absolute path or relative to the gmx mdrun calling location", c_referenceDensityFileNameTag_);
170         addDensityFittingMdpOutputValue(builder, referenceDensityFileName_, c_referenceDensityFileNameTag_);
171         addDensityFittingMdpOutputValue(builder, parameters_.calculationIntervalInSteps_, c_everyNStepsTag_);
172         addDensityFittingMdpOutputValueComment(builder, "; Normalize the sum of density voxel values to one", c_normalizeDensitiesTag_);
173         addDensityFittingMdpOutputValue(builder, parameters_.normalizeDensities_, c_normalizeDensitiesTag_);
174     }
175 }
176
177 void DensityFittingOptions::initMdpOptions(IOptionsContainerWithSections *options)
178 {
179     auto section = options->addSection(OptionSection(DensityFittingModuleInfo::name_.c_str()));
180
181     section.addOption(BooleanOption(c_activeTag_.c_str()).store(&parameters_.active_));
182     section.addOption(StringOption(c_groupTag_.c_str()).store(&groupString_));
183
184     section.addOption(EnumOption<DensitySimilarityMeasureMethod>(c_similarityMeasureTag_.c_str())
185                           .enumValue(c_densitySimilarityMeasureMethodNames.m_elements)
186                           .store(&parameters_.similarityMeasureMethod_));
187
188     section.addOption(EnumOption<DensityFittingAmplitudeMethod>(c_amplitudeMethodTag_.c_str())
189                           .enumValue(c_densityFittingAmplitudeMethodNames.m_elements)
190                           .store(&parameters_.amplitudeLookupMethod_));
191
192     section.addOption(RealOption(c_forceConstantTag_.c_str()).store(&parameters_.forceConstant_));
193     section.addOption(RealOption(c_gaussianTransformSpreadingWidthTag_.c_str()).store(&parameters_.gaussianTransformSpreadingWidth_));
194     section.addOption(RealOption(c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_.c_str()).store(&parameters_.gaussianTransformSpreadingRangeInMultiplesOfWidth_));
195     section.addOption(StringOption(c_referenceDensityFileNameTag_.c_str()).store(&referenceDensityFileName_));
196     section.addOption(Int64Option(c_everyNStepsTag_.c_str()).store(&parameters_.calculationIntervalInSteps_));
197     section.addOption(BooleanOption(c_normalizeDensitiesTag_.c_str()).store(&parameters_.normalizeDensities_));
198 }
199
200 bool DensityFittingOptions::active() const
201 {
202     return parameters_.active_;
203 }
204
205 const DensityFittingParameters &DensityFittingOptions::buildParameters()
206 {
207     // the options modules does not know unsigned integers so any input of this
208     // kind is rectified here
209     if (parameters_.calculationIntervalInSteps_ < 1)
210     {
211         parameters_.calculationIntervalInSteps_ = 1;
212     }
213     return parameters_;
214 }
215
216 void DensityFittingOptions::setFitGroupIndices(const IndexGroupsAndNames &indexGroupsAndNames)
217 {
218     if (!parameters_.active_)
219     {
220         return;
221     }
222     parameters_.indices_ = indexGroupsAndNames.indices(groupString_);
223 }
224
225 void DensityFittingOptions::writeInternalParametersToKvt(KeyValueTreeObjectBuilder treeBuilder)
226 {
227     auto groupIndexAdder = treeBuilder.addUniformArray<std::int64_t>(DensityFittingModuleInfo::name_ + "-" + c_groupTag_);
228     for (const auto &indexValue : parameters_.indices_)
229     {
230         groupIndexAdder.addValue(indexValue);
231     }
232 }
233
234 void DensityFittingOptions::readInternalParametersFromKvt(const KeyValueTreeObject &tree)
235 {
236     if (!parameters_.active_)
237     {
238         return;
239     }
240
241     if (!tree.keyExists(DensityFittingModuleInfo::name_ + "-" + c_groupTag_))
242     {
243         GMX_THROW(InconsistentInputError(
244                           "Cannot find atom index vector required for density guided simulation."));
245     }
246     auto kvtIndexArray = tree[DensityFittingModuleInfo::name_ + "-" + c_groupTag_].asArray().values();
247     parameters_.indices_.resize(kvtIndexArray.size());
248     std::transform(std::begin(kvtIndexArray), std::end(kvtIndexArray), std::begin(parameters_.indices_),
249                    [](const KeyValueTreeValue &val) { return val.cast<std::int64_t>(); });
250 }
251
252 void DensityFittingOptions::checkEnergyCaluclationFrequency(EnergyCalculationFrequencyErrors * energyCalculationFrequencyErrors) const
253 {
254     if (energyCalculationFrequencyErrors->energyCalculationIntervalInSteps() % parameters_.calculationIntervalInSteps_ != 0)
255     {
256         energyCalculationFrequencyErrors->addError("nstcalcenergy (" +
257                                                    toString(energyCalculationFrequencyErrors->energyCalculationIntervalInSteps())
258                                                    + ") is not a multiple of " + DensityFittingModuleInfo::name_ + "-" + c_everyNStepsTag_
259                                                    + " (" + toString(parameters_.calculationIntervalInSteps_) + ") .");
260     }
261 }
262
263 const std::string &DensityFittingOptions::referenceDensityFileName() const
264 {
265     return referenceDensityFileName_;
266 }
267 } // namespace gmx