2 * This file is part of the GROMACS molecular simulation package.
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.
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 * Declares data structure and utilities for density fitting
39 * \author Christian Blau <blau@kth.se>
40 * \ingroup module_applied_forces
44 #include "densityfitting.h"
48 #include "gromacs/domdec/localatomsetmanager.h"
49 #include "gromacs/fileio/mrcdensitymap.h"
50 #include "gromacs/math/multidimarray.h"
51 #include "gromacs/mdtypes/imdmodule.h"
52 #include "gromacs/utility/classhelpers.h"
53 #include "gromacs/utility/exceptions.h"
54 #include "gromacs/utility/keyvaluetreebuilder.h"
55 #include "gromacs/utility/mdmodulenotification.h"
57 #include "densityfittingforceprovider.h"
58 #include "densityfittingoptions.h"
59 #include "densityfittingoutputprovider.h"
64 class IMdpOptionProvider;
65 class DensityFittingForceProvider;
71 * \brief Collect density fitting parameters only available during simulation setup.
73 * \todo Implement builder pattern that will not use unqiue_ptr to check if
74 * parameters have been set or not.
76 * To build the density fitting force provider during simulation setup,
77 * the DensityFitting class needs access to parameters that become available
78 * only during simulation setup.
80 * This class collects these parameters via MdModuleNotifications in the
81 * simulation setup phase and provides a check if all necessary parameters have
84 class DensityFittingSimulationParameterSetup
87 DensityFittingSimulationParameterSetup() = default;
89 /*! \brief Set the local atom set for the density fitting.
90 * \param[in] localAtomSet of atoms to be fitted
92 void setLocalAtomSet(const LocalAtomSet &localAtomSet)
94 localAtomSet_ = std::make_unique<LocalAtomSet>(localAtomSet);
97 /*! \brief Return local atom set for density fitting.
98 * \throws InternalError if local atom set is not set
99 * \returns local atom set for density fitting
101 const LocalAtomSet &localAtomSet() const
103 if (localAtomSet_ == nullptr)
106 InternalError("Local atom set is not set for density "
107 "guided simulation."));
109 return *localAtomSet_;
112 /*! \brief Return transformation into density lattice.
113 * \throws InternalError if transformation into density lattice is not set
114 * \returns transformation into density lattice
116 const TranslateAndScale &transformationToDensityLattice() const
118 if (transformationToDensityLattice_ == nullptr)
121 InternalError("Transformation to reference density not set for density guided simulation."));
123 return *transformationToDensityLattice_;
125 /*! \brief Return reference density
126 * \throws InternalError if reference density is not set
127 * \returns the reference density
129 basic_mdspan<const float, dynamicExtents3D>
130 referenceDensity() const
132 if (referenceDensity_ == nullptr)
134 GMX_THROW(InternalError("Reference density not set for density guided simulation."));
136 return referenceDensity_->asConstView();
139 /*! \brief Reads the reference density from file.
141 * Reads and check file, then set and communicate the internal
142 * parameters related to the reference density with the file data.
144 * \throws FileIOError if reading from file was not successful
146 void readReferenceDensityFromFile(const std::string &referenceDensityFileName)
148 MrcDensityMapOfFloatFromFileReader reader(referenceDensityFileName);
149 referenceDensity_ = std::make_unique<MultiDimArray<std::vector<float>, dynamicExtents3D> >
150 (reader.densityDataCopy());
151 transformationToDensityLattice_
152 = std::make_unique<TranslateAndScale>(reader.transformationToDensityLattice());
154 /*! \brief Set the periodic boundary condition via MdModuleNotifier.
156 * The pbc type is wrapped in PeriodicBoundaryConditionType to
157 * allow the MdModuleNotifier to statically distinguish the callback
158 * function type from other 'int' function callbacks.
160 * \param[in] pbc MdModuleNotification class that contains a variable
161 * that enumerates the periodic boundary condition.
163 void setPeriodicBoundaryConditionType(PeriodicBoundaryConditionType pbc)
165 pbcType_ = std::make_unique<int>(pbc.pbcType);
168 //! Get the periodic boundary conditions
169 int periodicBoundaryConditionType()
171 if (pbcType_ == nullptr)
173 GMX_THROW(InternalError("Periodic boundary condition enum not set for density guided simulation."));
179 //! The reference density to fit to
180 std::unique_ptr<MultiDimArray<std::vector<float>, dynamicExtents3D> > referenceDensity_;
181 //! The coordinate transformation into the reference density
182 std::unique_ptr<TranslateAndScale> transformationToDensityLattice_;
183 //! The local atom set to act on
184 std::unique_ptr<LocalAtomSet> localAtomSet_;
185 //! The type of periodic boundary conditions in the simulation
186 std::unique_ptr<int> pbcType_;
188 GMX_DISALLOW_COPY_AND_ASSIGN(DensityFittingSimulationParameterSetup);
192 * \brief Density fitting
194 * Class that implements the density fitting forces and potential
195 * \note the virial calculation is not yet implemented
197 class DensityFitting final : public IMDModule
200 /*! \brief Construct the density fitting module.
202 * \param[in] notifier allows the module to subscribe to notifications from MdModules.
204 * The density fitting code subscribes to these notifications:
205 * - setting atom group indices in the densityFittingOptions_ by
206 * taking a parmeter const IndexGroupsAndNames &
207 * - storing its internal parameters in a tpr file by writing to a
208 * key-value-tree during pre-processing by a function taking a
209 * KeyValueTreeObjectBuilder as parameter
210 * - reading its internal parameters from a key-value-tree during
211 * simulation setup by taking a const KeyValueTreeObject & parameter
212 * - constructing local atom sets in the simulation parameter setup
213 * by taking a LocalAtomSetManager * as parameter
214 * - the type of periodic boundary conditions that are used
215 * by taking a PeriodicBoundaryConditionTypeEnum as parameter
217 explicit DensityFitting(MdModulesNotifier *notifier)
219 // Callbacks for several kinds of MdModuleNotification are created
220 // and subscribed, and will be dispatched correctly at run time
221 // based on the type of the parameter required by the lambda.
223 // Setting atom group indices
224 const auto setFitGroupIndicesFunction = [this](const IndexGroupsAndNames &indexGroupsAndNames) {
225 densityFittingOptions_.setFitGroupIndices(indexGroupsAndNames);
227 notifier->notifier_.subscribe(setFitGroupIndicesFunction);
229 // Writing internal parameters during pre-processing
230 const auto writeInternalParametersFunction = [this](KeyValueTreeObjectBuilder treeBuilder) {
231 densityFittingOptions_.writeInternalParametersToKvt(treeBuilder);
233 notifier->notifier_.subscribe(writeInternalParametersFunction);
235 // Reading internal parameters during simulation setup
236 const auto readInternalParametersFunction = [this](const KeyValueTreeObject &tree) {
237 densityFittingOptions_.readInternalParametersFromKvt(tree);
239 notifier->notifier_.subscribe(readInternalParametersFunction);
241 // constructing local atom sets during simulation setup
242 const auto setLocalAtomSetFunction = [this](LocalAtomSetManager *localAtomSetManager) {
243 this->constructLocalAtomSet(localAtomSetManager);
245 notifier->notifier_.subscribe(setLocalAtomSetFunction);
247 // constructing local atom sets during simulation setup
248 const auto setPeriodicBoundaryContionsFunction = [this](PeriodicBoundaryConditionType pbc) {
249 this->densityFittingSimulationParameters_.setPeriodicBoundaryConditionType(pbc);
251 notifier->notifier_.subscribe(setPeriodicBoundaryContionsFunction);
253 // adding output to energy file
254 const auto requestEnergyOutput
255 = [this](MdModulesEnergyOutputToDensityFittingRequestChecker *energyOutputRequest) {
256 this->setEnergyOutputRequest(energyOutputRequest);
258 notifier->notifier_.subscribe(requestEnergyOutput);
262 IMdpOptionProvider *mdpOptionProvider() override { return &densityFittingOptions_; }
264 //! Add this module to the force providers if active
265 void initForceProviders(ForceProviders *forceProviders) override
267 if (densityFittingOptions_.active())
269 const auto ¶meters = densityFittingOptions_.buildParameters();
270 densityFittingSimulationParameters_.readReferenceDensityFromFile(densityFittingOptions_.referenceDensityFileName());
271 forceProvider_ = std::make_unique<DensityFittingForceProvider>(
273 densityFittingSimulationParameters_.referenceDensity(),
274 densityFittingSimulationParameters_.transformationToDensityLattice(),
275 densityFittingSimulationParameters_.localAtomSet(),
276 densityFittingSimulationParameters_.periodicBoundaryConditionType());
277 forceProviders->addForceProvider(forceProvider_.get());
281 //! This MDModule provides its own output
282 IMDOutputProvider *outputProvider() override { return &densityFittingOutputProvider_; }
284 /*! \brief Set up the local atom sets that are used by this module.
286 * \note When density fitting is set up with MdModuleNotification in
287 * the constructor, this function is called back.
289 * \param[in] localAtomSetManager the manager to add local atom sets.
291 void constructLocalAtomSet(LocalAtomSetManager * localAtomSetManager)
293 if (densityFittingOptions_.active())
296 = localAtomSetManager->add(densityFittingOptions_.buildParameters().indices_);
297 densityFittingSimulationParameters_.setLocalAtomSet(atomSet);
301 /*! \brief Request energy output to energy file during simulation.
303 void setEnergyOutputRequest(MdModulesEnergyOutputToDensityFittingRequestChecker *energyOutputRequest)
305 energyOutputRequest->energyOutputToDensityFitting_ = densityFittingOptions_.active();
309 //! The output provider
310 DensityFittingOutputProvider densityFittingOutputProvider_;
311 //! The options provided for density fitting
312 DensityFittingOptions densityFittingOptions_;
313 //! Object that evaluates the forces
314 std::unique_ptr<DensityFittingForceProvider> forceProvider_;
315 /*! \brief Parameters for density fitting that become available at
316 * simulation setup time.
318 DensityFittingSimulationParameterSetup densityFittingSimulationParameters_;
320 GMX_DISALLOW_COPY_AND_ASSIGN(DensityFitting);
325 std::unique_ptr<IMDModule> DensityFittingModuleInfo::create(MdModulesNotifier * notifier)
327 return std::make_unique<DensityFitting>(notifier);
330 const std::string DensityFittingModuleInfo::name_ = "density-guided-simulation";