Add density fitting parameters
authorChristian Blau <cblau@gwdg.de>
Fri, 23 Aug 2019 12:24:13 +0000 (14:24 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Tue, 27 Aug 2019 15:03:37 +0000 (17:03 +0200)
Adds the mdp parameters that are required for running a
density fitting simulation.

refs #2282

Change-Id: Id22ceee0f0036450173522bd5e85af06542a2c39

src/gromacs/applied_forces/densityfittingoptions.cpp
src/gromacs/applied_forces/densityfittingoptions.h
src/gromacs/applied_forces/densityfittingparameters.cpp
src/gromacs/applied_forces/densityfittingparameters.h
src/gromacs/applied_forces/tests/densityfittingoptions.cpp
src/gromacs/math/densityfit.h

index 7d5400d5b71c170078362a8e5345ac023cee906d..68f1fb742b7eb4aed0672ef6478021ccfe2c1cbc 100644 (file)
@@ -44,6 +44,7 @@
 #include "densityfittingoptions.h"
 
 #include "gromacs/applied_forces/densityfitting.h"
+#include "gromacs/math/densityfit.h"
 #include "gromacs/options/basicoptions.h"
 #include "gromacs/options/optionsection.h"
 #include "gromacs/selection/indexutil.h"
@@ -52,6 +53,8 @@
 #include "gromacs/utility/keyvaluetreetransform.h"
 #include "gromacs/utility/strconvert.h"
 
+#include "densityfittingamplitudelookup.h"
+
 namespace gmx
 {
 
@@ -111,6 +114,11 @@ void DensityFittingOptions::initMdpTransform(IKeyValueTreeTransformRules * rules
         };
     densityfittingMdpTransformFromString<bool>(rules, &fromStdString<bool>, c_activeTag_);
     densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_groupTag_);
+    densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_similarityMeasureTag_);
+    densityfittingMdpTransformFromString<std::string>(rules, stringIdentityTransform, c_amplitudeMethodTag_);
+    densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>, c_forceConstantTag_);
+    densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>, c_gaussianTransformSpreadingWidthTag_);
+    densityfittingMdpTransformFromString<real>(rules, &fromStdString<real>, c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_);
 }
 
 void DensityFittingOptions::buildMdpOutput(KeyValueTreeObjectBuilder *builder) const
@@ -119,6 +127,17 @@ void DensityFittingOptions::buildMdpOutput(KeyValueTreeObjectBuilder *builder) c
     if (parameters_.active_)
     {
         addDensityFittingMdpOutputValue(builder, groupString_, c_groupTag_);
+        addDensityFittingMdpOutputValue<std::string>(builder,
+                                                     c_densitySimilarityMeasureMethodNames[parameters_.similarityMeasureMethod_],
+                                                     c_similarityMeasureTag_);
+
+        addDensityFittingMdpOutputValue<std::string>(builder,
+                                                     c_densityFittingAmplitudeMethodNames[parameters_.amplitudeLookupMethod_],
+                                                     c_amplitudeMethodTag_);
+
+        addDensityFittingMdpOutputValue(builder, parameters_.forceConstant_, c_forceConstantTag_);
+        addDensityFittingMdpOutputValue(builder, parameters_.gaussianTransformSpreadingWidth_, c_gaussianTransformSpreadingWidthTag_);
+        addDensityFittingMdpOutputValue(builder, parameters_.gaussianTransformSpreadingRangeInMultiplesOfWidth_, c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_);
     }
 }
 
@@ -128,6 +147,18 @@ void DensityFittingOptions::initMdpOptions(IOptionsContainerWithSections *option
 
     section.addOption(BooleanOption(c_activeTag_.c_str()).store(&parameters_.active_));
     section.addOption(StringOption(c_groupTag_.c_str()).store(&groupString_));
+
+    section.addOption(EnumOption<DensitySimilarityMeasureMethod>(c_similarityMeasureTag_.c_str())
+                          .enumValue(c_densitySimilarityMeasureMethodNames.m_elements)
+                          .store(&parameters_.similarityMeasureMethod_));
+
+    section.addOption(EnumOption<DensityFittingAmplitudeMethod>(c_amplitudeMethodTag_.c_str())
+                          .enumValue(c_densityFittingAmplitudeMethodNames.m_elements)
+                          .store(&parameters_.amplitudeLookupMethod_));
+
+    section.addOption(RealOption(c_forceConstantTag_.c_str()).store(&parameters_.forceConstant_));
+    section.addOption(RealOption(c_gaussianTransformSpreadingWidthTag_.c_str()).store(&parameters_.gaussianTransformSpreadingWidth_));
+    section.addOption(RealOption(c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_.c_str()).store(&parameters_.gaussianTransformSpreadingRangeInMultiplesOfWidth_));
 }
 
 bool DensityFittingOptions::active() const
index e0c73942692ffa385169a556ea1a485ad51d6fae..d57ab61b5e76156dfea26c90ae8724e38e96d4d9 100644 (file)
@@ -107,6 +107,16 @@ class DensityFittingOptions final : public IMdpOptionProvider
         const std::string        c_groupTag_  = "group";
         std::string              groupString_ = "protein";
 
+        const std::string        c_similarityMeasureTag_ = "similarity-measure";
+
+        const std::string        c_amplitudeMethodTag_ = "amplitude-method";
+
+        const std::string        c_forceConstantTag_ = "force-constant";
+
+        const std::string        c_gaussianTransformSpreadingWidthTag_ = "gaussian-transform-spreading-width";
+        const std::string        c_gaussianTransformSpreadingRangeInMultiplesOfWidthTag_
+            = "gaussian-transform-spreading-range-in-multiples-of-width";
+
         DensityFittingParameters parameters_;
 };
 
index f390d439e1366aee39df94bbc3b5f9f92da86d96..68af26e54ee4c48778f8e70025d57fbf1e29b383 100644 (file)
@@ -39,6 +39,7 @@
 
 namespace gmx
 {
+
 bool operator==(const DensityFittingParameters &lhs, const DensityFittingParameters &rhs)
 {
     if (lhs.active_ != rhs.active_)
@@ -49,6 +50,26 @@ bool operator==(const DensityFittingParameters &lhs, const DensityFittingParamet
     {
         return false;
     }
+    if (lhs.similarityMeasureMethod_ != rhs.similarityMeasureMethod_)
+    {
+        return false;
+    }
+    if (lhs.amplitudeLookupMethod_ != rhs.amplitudeLookupMethod_)
+    {
+        return false;
+    }
+    if (lhs.forceConstant_ != rhs.forceConstant_)
+    {
+        return false;
+    }
+    if (lhs.gaussianTransformSpreadingWidth_ != rhs.gaussianTransformSpreadingWidth_)
+    {
+        return false;
+    }
+    if (lhs.gaussianTransformSpreadingRangeInMultiplesOfWidth_ != rhs.gaussianTransformSpreadingRangeInMultiplesOfWidth_)
+    {
+        return false;
+    }
     return true;
 }
 
index 25fcc54b56220fdc1bf9ddf83075afc34c373e79..f585e26b465c2161b16965986a9a0adccc89331d 100644 (file)
 
 #include <vector>
 
+#include "gromacs/math/densityfit.h"
 #include "gromacs/utility/basedefinitions.h"
 
+#include "densityfittingamplitudelookup.h"
+
 namespace gmx
 {
 
@@ -57,12 +60,22 @@ namespace gmx
 struct DensityFittingParameters
 {
     //! Indicate if density fitting is active
-    bool               active_ = false;
+    bool                           active_ = false;
     //! Indices of the atoms that shall be fit to the density
-    std::vector<index> indices_;
+    std::vector<index>             indices_;
+    //! Determines how to measure similarity between simulated and reference density
+    DensitySimilarityMeasureMethod similarityMeasureMethod_ = DensitySimilarityMeasureMethod::innerProduct;
+    //! Determines with what weight atoms are spread
+    DensityFittingAmplitudeMethod  amplitudeLookupMethod_ = DensityFittingAmplitudeMethod::Unity;
+    //! The force constant to be used for the density fitting
+    real                           forceConstant_ = 1e9;
+    //! The spreading width used for the gauss transform of atoms onto the density grid
+    real                           gaussianTransformSpreadingWidth_ = 0.2;
+    //! The spreading range for spreading atoms onto the grid in multiples of the spreading width
+    real                           gaussianTransformSpreadingRangeInMultiplesOfWidth_ = 4.0;
 };
 
-/*! \brief Check if two structs holding density fitting parameters are equal.
+/*!\brief Check if two structs holding density fitting parameters are equal.
  *
  * \param[in] lhs left hand side to be compared
  * \param[in] rhs right hand side to be compared
@@ -70,7 +83,7 @@ struct DensityFittingParameters
  */
 bool operator==(const DensityFittingParameters &lhs, const DensityFittingParameters &rhs);
 
-/*! \brief Check if two structs holding density fitting parameters are not equal.
+/*!\brief Check if two structs holding density fitting parameters are not equal.
  *
  * \param[in] lhs left hand side to be compared
  * \param[in] rhs right hand side to be compared
index a6c7f57d1d35a5a8bf5cdb14ad01956635b826bc..7acd948da0ff2fd657046205f1cd7b3902ec6292 100644 (file)
@@ -137,7 +137,14 @@ class DensityFittingOptionsTest : public ::testing::Test
 
 TEST_F(DensityFittingOptionsTest, DefaultParameters)
 {
-    EXPECT_FALSE(densityFittingOptions_.buildParameters().active_);
+    const auto defaultParameters = densityFittingOptions_.buildParameters();
+    EXPECT_FALSE(defaultParameters.active_);
+    EXPECT_EQ(0, defaultParameters.indices_.size());
+    EXPECT_EQ(DensitySimilarityMeasureMethod::innerProduct, defaultParameters.similarityMeasureMethod_);
+    EXPECT_EQ(DensityFittingAmplitudeMethod::Unity, defaultParameters.amplitudeLookupMethod_);
+    EXPECT_REAL_EQ(1e9, defaultParameters.forceConstant_);
+    EXPECT_REAL_EQ(0.2, defaultParameters.gaussianTransformSpreadingWidth_);
+    EXPECT_REAL_EQ(4.0, defaultParameters.gaussianTransformSpreadingRangeInMultiplesOfWidth_);
 }
 
 TEST_F(DensityFittingOptionsTest, OptionSetsActive)
@@ -147,7 +154,7 @@ TEST_F(DensityFittingOptionsTest, OptionSetsActive)
     EXPECT_TRUE(densityFittingOptions_.buildParameters().active_);
 }
 
-TEST_F(DensityFittingOptionsTest, OutputDefaultValues)
+TEST_F(DensityFittingOptionsTest, OutputNoDefaultValuesWhenInactive)
 {
     // Transform module data into a flat key-value tree for output.
 
@@ -165,6 +172,36 @@ TEST_F(DensityFittingOptionsTest, OutputDefaultValues)
     EXPECT_EQ(stream.toString(), std::string("density-guided-simulation-active = false\n"));
 }
 
+TEST_F(DensityFittingOptionsTest, OutputDefaultValuesWhenActive)
+{
+    setFromMdpValues(densityFittingSetActiveAsMdpValues());
+    // Transform module data into a flat key-value tree for output.
+
+    StringOutputStream        stream;
+    KeyValueTreeBuilder       builder;
+    KeyValueTreeObjectBuilder builderObject = builder.rootObject();
+
+    densityFittingOptions_.buildMdpOutput(&builderObject);
+    {
+        TextWriter writer(&stream);
+        writeKeyValueTreeAsMdp(&writer, builder.build());
+    }
+    stream.close();
+    std::string expected
+        = {
+        "density-guided-simulation-active = true\n"
+        "density-guided-simulation-group = protein\n"
+        "density-guided-simulation-similarity-measure = inner-product\n"
+        "density-guided-simulation-amplitude-method = unity\n"
+        "density-guided-simulation-force-constant = 1e+09\n"
+        "density-guided-simulation-gaussian-transform-spreading-width = 0.2\n"
+        "density-guided-simulation-gaussian-transform-spreading-range-in-multiples-of-width = 4\n"
+        };
+
+    EXPECT_EQ(expected, stream.toString());
+}
+
+
 TEST_F(DensityFittingOptionsTest, CanConvertGroupStringToIndexGroup)
 {
     setFromMdpValues(densityFittingSetActiveAsMdpValues());
index f0d6758652f248e07041c959ab8dcdeb0bf5e9fd..3ad5f5c9ece14952e03195b98c8dbf901a21f121 100644 (file)
@@ -32,7 +32,7 @@
  * To help us fund GROMACS development, we humbly ask that you cite
  * the research papers on the package. Check out http://www.gromacs.org.
  */
-/*! \libinternal \file
+/*! \inlibraryapi \file
  * \brief
  * Declares density similarity measures and their derivatives.
  *
@@ -44,6 +44,7 @@
 
 #include "gromacs/mdspan/extensions.h"
 #include "gromacs/utility/classhelpers.h"
+#include "gromacs/utility/enumerationhelpers.h"
 
 namespace gmx
 {
@@ -61,12 +62,12 @@ enum class DensitySimilarityMeasureMethod
      * \f]
      */
     innerProduct,
-    count,
+    Count,
 };
 
 //! Name the methods that may be used to evaluate similarity between densities
-const char * const c_densitySimilarityMeasureMethodNames[static_cast<int>(DensitySimilarityMeasureMethod::count)]
-    = {"inner-product"};
+const EnumerationArray<DensitySimilarityMeasureMethod, const char *const>
+c_densitySimilarityMeasureMethodNames = {{ "inner-product" }};
 
 /* Forward declaration of implementation class outside class to allow
  * choose implementation class during construction of the DensitySimilarityMeasure*/