Introduce expanded ensemble in modular simulator
authorPascal Merz <pascal.merz@me.com>
Wed, 22 Sep 2021 07:53:11 +0000 (07:53 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 22 Sep 2021 07:53:11 +0000 (07:53 +0000)
src/gromacs/fileio/checkpoint.cpp
src/gromacs/mdlib/md_support.cpp
src/gromacs/modularsimulator/expandedensembleelement.cpp [new file with mode: 0644]
src/gromacs/modularsimulator/expandedensembleelement.h [new file with mode: 0644]
src/gromacs/modularsimulator/modularsimulator.cpp

index 673ae2c689ac74ccae37c29913e0e12979e59d76..e8d2409c2d47687e4098f1c2e7184653d2103d9b 100644 (file)
@@ -1904,9 +1904,11 @@ static int do_cpt_df_hist(XDR* xd, int fflags, int nlambda, df_history_t** dfhis
         return 0;
     }
 
+    std::unique_ptr<df_history_t> localDFHistory = nullptr;
     if (*dfhistPtr == nullptr)
     {
-        snew(*dfhistPtr, 1);
+        localDFHistory        = std::make_unique<df_history_t>();
+        *dfhistPtr            = localDFHistory.get();
         (*dfhistPtr)->nlambda = nlambda;
         init_df_history(*dfhistPtr, nlambda);
     }
index 46ebb11617187e02736905a240e3b500fbe864e6..297a41bf5aadfdec56f6c8ccb9ecd8e434342e13 100644 (file)
@@ -584,7 +584,7 @@ void set_state_entries(t_state* state, const t_inputrec* ir, bool useModularSimu
     init_gtc_state(state, state->ngtc, state->nnhpres, ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
     init_ekinstate(&state->ekinstate, ir);
 
-    if (ir->bExpanded)
+    if (ir->bExpanded && !useModularSimulator)
     {
         snew(state->dfhist, 1);
         init_df_history(state->dfhist, ir->fepvals->n_lambda);
diff --git a/src/gromacs/modularsimulator/expandedensembleelement.cpp b/src/gromacs/modularsimulator/expandedensembleelement.cpp
new file mode 100644 (file)
index 0000000..17ec878
--- /dev/null
@@ -0,0 +1,223 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2021, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \internal \file
+ * \brief Defines the expanded ensemble element for the modular simulator
+ *
+ * \author Pascal Merz <pascal.merz@me.com>
+ * \ingroup module_modularsimulator
+ */
+
+#include "gmxpre.h"
+
+#include "expandedensembleelement.h"
+
+#include "gromacs/domdec/distribute.h"
+#include "gromacs/mdlib/expanded.h"
+#include "gromacs/mdlib/stat.h"
+#include "gromacs/mdtypes/checkpointdata.h"
+#include "gromacs/mdtypes/commrec.h"
+#include "gromacs/mdtypes/df_history.h"
+#include "gromacs/mdtypes/inputrec.h"
+#include "gromacs/utility/fatalerror.h"
+
+#include "energydata.h"
+#include "simulatoralgorithm.h"
+
+namespace gmx
+{
+
+void ExpandedEnsembleElement::apply(Step step, bool doLambdaStep, bool doLog)
+{
+    if (doLambdaStep)
+    {
+        const int newFepState =
+                expandedEnsembleUpdateLambdaState(fplog_,
+                                                  inputrec_,
+                                                  energyData_->enerdata(),
+                                                  freeEnergyPerturbationData_->currentFEPState(),
+                                                  dfhist_.get(),
+                                                  step);
+        // Set new state at next step
+        fepStateSetting_->setNewState(newFepState, step + 1);
+    }
+    if (doLog)
+    {
+        /* only needed if doing expanded ensemble */
+        PrintFreeEnergyInfoToFile(fplog_,
+                                  inputrec_->fepvals.get(),
+                                  inputrec_->expandedvals.get(),
+                                  inputrec_->bSimTemp ? inputrec_->simtempvals.get() : nullptr,
+                                  dfhist_.get(),
+                                  freeEnergyPerturbationData_->currentFEPState(),
+                                  inputrec_->nstlog,
+                                  step);
+    }
+}
+
+void ExpandedEnsembleElement::elementSetup()
+{
+    // Check nstexpanded here, because the grompp check was broken (#2714)
+    if (inputrec_->expandedvals->nstexpanded % inputrec_->nstcalcenergy != 0)
+    {
+        gmx_fatal(FARGS,
+                  "With expanded ensemble, nstexpanded should be a multiple of nstcalcenergy");
+    }
+    init_expanded_ensemble(restoredFromCheckpoint_, inputrec_, dfhist_.get());
+}
+
+void ExpandedEnsembleElement::scheduleTask(Step step, Time /*unused*/, const RegisterRunFunction& registerRunFunction)
+{
+    const bool isFirstStep  = (step == initialStep_);
+    const bool doLambdaStep = (do_per_step(step, frequency_) && !isFirstStep);
+    const bool doLog        = (isMasterRank_ && step == nextLogWritingStep_ && (fplog_ != nullptr));
+
+    if (doLambdaStep || doLog)
+    {
+        registerRunFunction([this, step, doLambdaStep, doLog]() { apply(step, doLambdaStep, doLog); });
+    }
+    if (doLambdaStep)
+    {
+        // We'll compute a new lambda state and want it applied for next step
+        fepStateSetting_->signalSettingStep(step + 1);
+    }
+}
+
+namespace
+{
+/*!
+ * \brief Enum describing the contents FreeEnergyPerturbationData::Element writes to modular checkpoint
+ *
+ * When changing the checkpoint content, add a new element just above Count, and adjust the
+ * checkpoint functionality.
+ */
+enum class CheckpointVersion
+{
+    Base, //!< First version of modular checkpointing
+    Count //!< Number of entries. Add new versions right above this!
+};
+constexpr auto c_currentVersion = CheckpointVersion(int(CheckpointVersion::Count) - 1);
+} // namespace
+
+template<CheckpointDataOperation operation>
+void ExpandedEnsembleElement::doCheckpointData(CheckpointData<operation>* checkpointData)
+{
+    checkpointVersion(checkpointData, "ExpandedEnsembleElement version", c_currentVersion);
+
+    dfhist_->doCheckpoint<operation>(checkpointData->subCheckpointData("dfhist"),
+                                     inputrec_->expandedvals->elamstats);
+}
+
+void ExpandedEnsembleElement::saveCheckpointState(std::optional<WriteCheckpointData> checkpointData,
+                                                  const t_commrec*                   cr)
+{
+    if (MASTER(cr))
+    {
+        doCheckpointData<CheckpointDataOperation::Write>(&checkpointData.value());
+    }
+}
+
+void ExpandedEnsembleElement::restoreCheckpointState(std::optional<ReadCheckpointData> checkpointData,
+                                                     const t_commrec*                  cr)
+{
+    if (MASTER(cr))
+    {
+        doCheckpointData<CheckpointDataOperation::Read>(&checkpointData.value());
+    }
+    if (DOMAINDECOMP(cr))
+    {
+        dd_distribute_dfhist(cr->dd, dfhist_.get());
+    }
+    restoredFromCheckpoint_ = true;
+}
+
+const std::string& ExpandedEnsembleElement::clientID()
+{
+    return identifier_;
+}
+
+std::optional<SignallerCallback> ExpandedEnsembleElement::registerLoggingCallback()
+{
+    if (isMasterRank_)
+    {
+        return [this](Step step, Time /*unused*/) { nextLogWritingStep_ = step; };
+    }
+    else
+    {
+        return std::nullopt;
+    }
+}
+
+ExpandedEnsembleElement::ExpandedEnsembleElement(bool                              isMasterRank,
+                                                 Step                              initialStep,
+                                                 int                               frequency,
+                                                 const EnergyData*                 energyData,
+                                                 const FreeEnergyPerturbationData* freeEnergyPerturbationData,
+                                                 FILE*                             fplog,
+                                                 const t_inputrec*                 inputrec) :
+    fepStateSetting_(freeEnergyPerturbationData->enableExternalFepStateSetting()),
+    isMasterRank_(isMasterRank),
+    initialStep_(initialStep),
+    frequency_(frequency),
+    nextLogWritingStep_(-1),
+    dfhist_(std::make_unique<df_history_t>()),
+    restoredFromCheckpoint_(false),
+    energyData_(energyData),
+    freeEnergyPerturbationData_(freeEnergyPerturbationData),
+    fplog_(fplog),
+    inputrec_(inputrec)
+{
+    init_df_history(dfhist_.get(), inputrec_->fepvals->n_lambda);
+}
+
+ISimulatorElement* ExpandedEnsembleElement::getElementPointerImpl(
+        LegacySimulatorData*                    legacySimulatorData,
+        ModularSimulatorAlgorithmBuilderHelper* builderHelper,
+        StatePropagatorData gmx_unused* statePropagatorData,
+        EnergyData*                     energyData,
+        FreeEnergyPerturbationData*     freeEnergyPerturbationData,
+        GlobalCommunicationHelper gmx_unused* globalCommunicationHelper,
+        ObservablesReducer* /*observablesReducer*/)
+{
+    return builderHelper->storeElement(std::make_unique<ExpandedEnsembleElement>(
+            MASTER(legacySimulatorData->cr),
+            legacySimulatorData->inputrec->init_step,
+            legacySimulatorData->inputrec->expandedvals->nstexpanded,
+            energyData,
+            freeEnergyPerturbationData,
+            legacySimulatorData->fplog,
+            legacySimulatorData->inputrec));
+}
+
+} // namespace gmx
diff --git a/src/gromacs/modularsimulator/expandedensembleelement.h b/src/gromacs/modularsimulator/expandedensembleelement.h
new file mode 100644 (file)
index 0000000..bf08567
--- /dev/null
@@ -0,0 +1,161 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2021, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \internal \file
+ * \brief Declares the expanded ensemble element for the modular simulator
+ *
+ * \author Pascal Merz <pascal.merz@me.com>
+ * \ingroup module_modularsimulator
+ *
+ * This header is only used within the modular simulator module
+ */
+
+#ifndef GMX_MODULARSIMULATOR_EXPANDEDENSEMBLEELEMENT_H
+#define GMX_MODULARSIMULATOR_EXPANDEDENSEMBLEELEMENT_H
+
+#include "freeenergyperturbationdata.h"
+
+struct df_history_t;
+struct t_inputrec;
+
+namespace gmx
+{
+enum class CheckpointDataOperation;
+class EnergyData;
+class GlobalCommunicationHelper;
+class LegacySimulatorData;
+class MDAtoms;
+class ModularSimulatorAlgorithmBuilderHelper;
+class ObservablesReducer;
+class StatePropagatorData;
+
+/*! \internal
+ * \ingroup module_modularsimulator
+ * \brief The expanded ensemble element
+ *
+ * This element periodically attempts Monte Carlo moves in lambda
+ * space and sets the new lambda state in FreeEnergyPerturbationData::Element.
+ */
+class ExpandedEnsembleElement final : public ISimulatorElement, public ICheckpointHelperClient, public ILoggingSignallerClient
+{
+public:
+    //! Constructor
+    explicit ExpandedEnsembleElement(bool                              isMasterRank,
+                                     Step                              initialStep,
+                                     int                               frequency,
+                                     const EnergyData*                 energyData,
+                                     const FreeEnergyPerturbationData* freeEnergyPerturbationData,
+                                     FILE*                             fplog,
+                                     const t_inputrec*                 inputrec);
+
+    //! Attempt lambda MC step and write log
+    void scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction) override;
+    //! Set up FEP history object
+    void elementSetup() override;
+    //! No teardown needed
+    void elementTeardown() override{};
+
+    //! ICheckpointHelperClient write checkpoint implementation
+    void saveCheckpointState(std::optional<WriteCheckpointData> checkpointData, const t_commrec* cr) override;
+    //! ICheckpointHelperClient read checkpoint implementation
+    void restoreCheckpointState(std::optional<ReadCheckpointData> checkpointData, const t_commrec* cr) override;
+    //! ICheckpointHelperClient key implementation
+    const std::string& clientID() override;
+
+    /*! \brief Factory method implementation
+     *
+     * \param legacySimulatorData  Pointer allowing access to simulator level data
+     * \param builderHelper  ModularSimulatorAlgorithmBuilder helper object
+     * \param statePropagatorData  Pointer to the \c StatePropagatorData object
+     * \param energyData  Pointer to the \c EnergyData object
+     * \param freeEnergyPerturbationData  Pointer to the \c FreeEnergyPerturbationData object
+     * \param globalCommunicationHelper  Pointer to the \c GlobalCommunicationHelper object
+     * \param observablesReducer          Pointer to the \c ObservablesReducer object
+     *
+     * \return  Pointer to the element to be added. Element needs to have been stored using \c storeElement
+     */
+    static ISimulatorElement* getElementPointerImpl(LegacySimulatorData* legacySimulatorData,
+                                                    ModularSimulatorAlgorithmBuilderHelper* builderHelper,
+                                                    StatePropagatorData*        statePropagatorData,
+                                                    EnergyData*                 energyData,
+                                                    FreeEnergyPerturbationData* freeEnergyPerturbationData,
+                                                    GlobalCommunicationHelper* globalCommunicationHelper,
+                                                    ObservablesReducer*        observablesReducer);
+
+private:
+    //! Use expanded ensemble to determine new FEP state or write log
+    void apply(Step step, bool doLambdaStep, bool doLog);
+
+    //! Helper object allowing to set new FEP state
+    FepStateSetting* fepStateSetting_;
+
+    //! Whether this runs on master
+    const bool isMasterRank_;
+    //! The initial Step
+    const Step initialStep_;
+    //! The frequency of lambda MC steps
+    const int frequency_;
+
+    //! ILoggingSignallerClient implementation
+    std::optional<SignallerCallback> registerLoggingCallback() override;
+    //! The next logging step
+    Step nextLogWritingStep_;
+
+    //! The free energy sampling history
+    std::unique_ptr<df_history_t> dfhist_;
+
+    //! CheckpointHelper identifier
+    const std::string identifier_ = "ExpandedEnsembleElement";
+    //! Helper function to read from / write to CheckpointData
+    template<CheckpointDataOperation operation>
+    void doCheckpointData(CheckpointData<operation>* checkpointData);
+    //! Whether this object was restored from checkpoint
+    bool restoredFromCheckpoint_;
+
+    // TODO: Clarify relationship to data objects and find a more robust alternative to raw pointers (#3583)
+    //! Pointer to the energy data
+    const EnergyData* energyData_;
+    //! Pointer to the free energy perturbation data
+    const FreeEnergyPerturbationData* freeEnergyPerturbationData_;
+
+    // Access to ISimulator data
+    //! Handles logging.
+    FILE* fplog_;
+    //! Contains user input mdp options.
+    const t_inputrec* inputrec_;
+};
+
+} // namespace gmx
+
+#endif // GMX_MODULARSIMULATOR_EXPANDEDENSEMBLEELEMENT_H
index 4e91aa95d633aaf613b979739f3afc8df147f1ec..d3fce2c0b5563bf8bce6df5c8821aaf171f55e58 100644 (file)
@@ -78,6 +78,7 @@
 #include "andersentemperaturecoupling.h"
 #include "computeglobalselement.h"
 #include "constraintelement.h"
+#include "expandedensembleelement.h"
 #include "firstorderpressurecoupling.h"
 #include "forceelement.h"
 #include "mttk.h"
@@ -154,7 +155,12 @@ void ModularSimulator::addIntegrationElements(ModularSimulatorAlgorithmBuilder*
             builder->add<ConstraintsElement<ConstraintVariable::Velocities>>();
         }
         builder->add<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>();
+        // Here, we have x / v / f at the full time step
         builder->add<StatePropagatorData::Element>();
+        if (legacySimulatorData_->inputrec->bExpanded)
+        {
+            builder->add<ExpandedEnsembleElement>();
+        }
         if (legacySimulatorData_->inputrec->etc == TemperatureCoupling::VRescale
             || legacySimulatorData_->inputrec->etc == TemperatureCoupling::Berendsen)
         {
@@ -246,6 +252,10 @@ void ModularSimulator::addIntegrationElements(ModularSimulatorAlgorithmBuilder*
         }
         // We have a full state at time t here
         builder->add<StatePropagatorData::Element>();
+        if (legacySimulatorData_->inputrec->bExpanded)
+        {
+            builder->add<ExpandedEnsembleElement>();
+        }
 
         // Propagate extended system variables from t to t+dt/2
         if (legacySimulatorData_->inputrec->epc == PressureCoupling::Mttk)
@@ -370,13 +380,6 @@ bool ModularSimulator::isInputCompatible(bool                             exitOn
     isInputCompatible =
             isInputCompatible
             && conditionalAssert(!doRerun, "Rerun is not supported by the modular simulator.");
-    isInputCompatible =
-            isInputCompatible
-            && conditionalAssert(inputrec->efep == FreeEnergyPerturbationType::No
-                                         || inputrec->efep == FreeEnergyPerturbationType::Yes
-                                         || inputrec->efep == FreeEnergyPerturbationType::SlowGrowth,
-                                 "Expanded ensemble free energy calculation is not "
-                                 "supported by the modular simulator.");
     isInputCompatible = isInputCompatible
                         && conditionalAssert(!inputrec->bPull,
                                              "Pulling is not supported by the modular simulator.");
@@ -447,10 +450,6 @@ bool ModularSimulator::isInputCompatible(bool                             exitOn
             isInputCompatible
             && conditionalAssert(!inputrec->bSimTemp,
                                  "Simulated tempering is not supported by the modular simulator.");
-    isInputCompatible = isInputCompatible
-                        && conditionalAssert(!inputrec->bExpanded,
-                                             "Expanded ensemble simulations are not supported by "
-                                             "the modular simulator.");
     isInputCompatible =
             isInputCompatible
             && conditionalAssert(!doEssentialDynamics,