2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020,2021, 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"
49 #include "gromacs/domdec/localatomset.h"
50 #include "gromacs/domdec/localatomsetmanager.h"
51 #include "gromacs/fileio/checkpoint.h"
52 #include "gromacs/fileio/mrcdensitymap.h"
53 #include "gromacs/math/coordinatetransformation.h"
54 #include "gromacs/math/multidimarray.h"
55 #include "gromacs/mdtypes/imdmodule.h"
56 #include "gromacs/utility/classhelpers.h"
57 #include "gromacs/utility/exceptions.h"
58 #include "gromacs/utility/keyvaluetreebuilder.h"
59 #include "gromacs/utility/mdmodulesnotifiers.h"
61 #include "densityfittingforceprovider.h"
62 #include "densityfittingoptions.h"
63 #include "densityfittingoutputprovider.h"
68 class IMdpOptionProvider;
69 class DensityFittingForceProvider;
75 * \brief Collect density fitting parameters only available during simulation setup.
77 * \todo Implement builder pattern that will not use unqiue_ptr to check if
78 * parameters have been set or not.
80 * To build the density fitting force provider during simulation setup,
81 * the DensityFitting class needs access to parameters that become available
82 * only during simulation setup.
84 * This class collects these parameters via MDModulesNotifiers in the
85 * simulation setup phase and provides a check if all necessary parameters have
88 class DensityFittingSimulationParameterSetup
91 DensityFittingSimulationParameterSetup() = default;
93 /*! \brief Set the local atom set for the density fitting.
94 * \param[in] localAtomSet of atoms to be fitted
96 void setLocalAtomSet(const LocalAtomSet& localAtomSet)
98 localAtomSet_ = std::make_unique<LocalAtomSet>(localAtomSet);
101 /*! \brief Return local atom set for density fitting.
102 * \throws InternalError if local atom set is not set
103 * \returns local atom set for density fitting
105 const LocalAtomSet& localAtomSet() const
107 if (localAtomSet_ == nullptr)
110 InternalError("Local atom set is not set for density "
111 "guided simulation."));
113 return *localAtomSet_;
116 /*! \brief Return transformation into density lattice.
117 * \throws InternalError if transformation into density lattice is not set
118 * \returns transformation into density lattice
120 const TranslateAndScale& transformationToDensityLattice() const
122 if (transformationToDensityLattice_ == nullptr)
124 GMX_THROW(InternalError(
125 "Transformation to reference density not set for density guided simulation."));
127 return *transformationToDensityLattice_;
129 /*! \brief Return reference density
130 * \throws InternalError if reference density is not set
131 * \returns the reference density
133 basic_mdspan<const float, dynamicExtents3D> referenceDensity() const
135 if (referenceDensity_ == nullptr)
137 GMX_THROW(InternalError("Reference density not set for density guided simulation."));
139 return referenceDensity_->asConstView();
142 /*! \brief Reads the reference density from file.
144 * Reads and check file, then set and communicate the internal
145 * parameters related to the reference density with the file data.
147 * \throws FileIOError if reading from file was not successful
149 void readReferenceDensityFromFile(const std::string& referenceDensityFileName)
151 MrcDensityMapOfFloatFromFileReader reader(referenceDensityFileName);
152 referenceDensity_ = std::make_unique<MultiDimArray<std::vector<float>, dynamicExtents3D>>(
153 reader.densityDataCopy());
154 transformationToDensityLattice_ =
155 std::make_unique<TranslateAndScale>(reader.transformationToDensityLattice());
158 //! Normalize the reference density so that the sum over all voxels is unity
159 void normalizeReferenceDensity()
161 if (referenceDensity_ == nullptr)
163 GMX_THROW(InternalError("Need to set reference density before normalizing it."));
166 const real sumOfDensityData = std::accumulate(
167 begin(referenceDensity_->asView()), end(referenceDensity_->asView()), 0.);
168 for (float& referenceDensityVoxel : referenceDensity_->asView())
170 referenceDensityVoxel /= sumOfDensityData;
173 /*! \brief Set the periodic boundary condition via MDModuleNotifier.
175 * The pbc type is wrapped in PeriodicBoundaryConditionType to
176 * allow the MDModuleNotifier to statically distinguish the callback
177 * function type from other 'int' function callbacks.
179 * \param[in] pbcType enumerates the periodic boundary condition.
181 void setPeriodicBoundaryConditionType(const PbcType& pbcType)
183 pbcType_ = std::make_unique<PbcType>(pbcType);
186 //! Get the periodic boundary conditions
187 PbcType periodicBoundaryConditionType()
189 if (pbcType_ == nullptr)
191 GMX_THROW(InternalError(
192 "Periodic boundary condition enum not set for density guided simulation."));
197 //! Set the simulation time step
198 void setSimulationTimeStep(double timeStep) { simulationTimeStep_ = timeStep; }
200 //! Return the simulation time step
201 double simulationTimeStep() { return simulationTimeStep_; }
204 //! The reference density to fit to
205 std::unique_ptr<MultiDimArray<std::vector<float>, dynamicExtents3D>> referenceDensity_;
206 //! The coordinate transformation into the reference density
207 std::unique_ptr<TranslateAndScale> transformationToDensityLattice_;
208 //! The local atom set to act on
209 std::unique_ptr<LocalAtomSet> localAtomSet_;
210 //! The type of periodic boundary conditions in the simulation
211 std::unique_ptr<PbcType> pbcType_;
212 //! The simulation time step
213 double simulationTimeStep_ = 1;
215 GMX_DISALLOW_COPY_AND_ASSIGN(DensityFittingSimulationParameterSetup);
219 * \brief Density fitting
221 * Class that implements the density fitting forces and potential
222 * \note the virial calculation is not yet implemented
224 class DensityFitting final : public IMDModule
228 //! Construct the density fitting module.
229 explicit DensityFitting() = default;
231 /*! \brief Request to be notified during pre-processing.
233 * \param[in] notifiers allows the module to subscribe to notifications from MDModules.
235 * The density fitting code subscribes to these notifications:
236 * - setting atom group indices in the densityFittingOptions_ from an
237 * index group string by taking a parmeter const IndexGroupsAndNames &
238 * - storing its internal parameters in a tpr file by writing to a
239 * key-value-tree during pre-processing by a function taking a
240 * KeyValueTreeObjectBuilder as parameter
242 void subscribeToPreProcessingNotifications(MDModulesNotifiers* notifiers) override
244 // Callbacks for several kinds of MDModulesNotifier are created
245 // and subscribed, and will be dispatched correctly at run time
246 // based on the type of the parameter required by the lambda.
248 if (!densityFittingOptions_.active())
253 // Setting the atom group indices from index group string
254 const auto setFitGroupIndicesFunction = [this](const IndexGroupsAndNames& indexGroupsAndNames) {
255 densityFittingOptions_.setFitGroupIndices(indexGroupsAndNames);
257 notifiers->preProcessingNotifier_.subscribe(setFitGroupIndicesFunction);
259 // Writing internal parameters during pre-processing
260 const auto writeInternalParametersFunction = [this](KeyValueTreeObjectBuilder treeBuilder) {
261 densityFittingOptions_.writeInternalParametersToKvt(treeBuilder);
263 notifiers->preProcessingNotifier_.subscribe(writeInternalParametersFunction);
265 // Checking for consistency with all .mdp options
266 const auto checkEnergyCaluclationFrequencyFunction =
267 [this](EnergyCalculationFrequencyErrors* energyCalculationFrequencyErrors) {
268 densityFittingOptions_.checkEnergyCaluclationFrequency(energyCalculationFrequencyErrors);
270 notifiers->preProcessingNotifier_.subscribe(checkEnergyCaluclationFrequencyFunction);
273 /*! \brief Request to be notified.
274 * The density fitting code subscribes to these notifications:
275 * - reading its internal parameters from a key-value-tree during
276 * simulation setup by taking a const KeyValueTreeObject & parameter
277 * - constructing local atom sets in the simulation parameter setup
278 * by taking a LocalAtomSetManager * as parameter
279 * - the type of periodic boundary conditions that are used
280 * by taking a PeriodicBoundaryConditionType as parameter
281 * - the writing of checkpoint data
282 * by taking a MDModulesWriteCheckpointData as parameter
283 * - the reading of checkpoint data
284 * by taking a MDModulesCheckpointReadingDataOnMaster as parameter
285 * - the broadcasting of checkpoint data
286 * by taking MDModulesCheckpointReadingBroadcast as parameter
288 void subscribeToSimulationSetupNotifications(MDModulesNotifiers* notifiers) override
290 if (!densityFittingOptions_.active())
295 // Reading internal parameters during simulation setup
296 const auto readInternalParametersFunction = [this](const KeyValueTreeObject& tree) {
297 densityFittingOptions_.readInternalParametersFromKvt(tree);
299 notifiers->simulationSetupNotifier_.subscribe(readInternalParametersFunction);
301 // constructing local atom sets during simulation setup
302 const auto setLocalAtomSetFunction = [this](LocalAtomSetManager* localAtomSetManager) {
303 this->constructLocalAtomSet(localAtomSetManager);
305 notifiers->simulationSetupNotifier_.subscribe(setLocalAtomSetFunction);
307 // constructing local atom sets during simulation setup
308 const auto setPeriodicBoundaryContionsFunction = [this](const PbcType& pbc) {
309 this->densityFittingSimulationParameters_.setPeriodicBoundaryConditionType(pbc);
311 notifiers->simulationSetupNotifier_.subscribe(setPeriodicBoundaryContionsFunction);
313 // setting the simulation time step
314 const auto setSimulationTimeStepFunction = [this](const SimulationTimeStep& simulationTimeStep) {
315 this->densityFittingSimulationParameters_.setSimulationTimeStep(simulationTimeStep.delta_t);
317 notifiers->simulationSetupNotifier_.subscribe(setSimulationTimeStepFunction);
319 // adding output to energy file
320 const auto requestEnergyOutput =
321 [this](MDModulesEnergyOutputToDensityFittingRequestChecker* energyOutputRequest) {
322 this->setEnergyOutputRequest(energyOutputRequest);
324 notifiers->simulationSetupNotifier_.subscribe(requestEnergyOutput);
326 // writing checkpoint data
327 const auto checkpointDataWriting = [this](MDModulesWriteCheckpointData checkpointData) {
328 forceProvider_->writeCheckpointData(checkpointData, DensityFittingModuleInfo::name_);
330 notifiers->checkpointingNotifier_.subscribe(checkpointDataWriting);
332 // reading checkpoint data
333 const auto checkpointDataReading = [this](MDModulesCheckpointReadingDataOnMaster checkpointData) {
334 densityFittingState_.readState(checkpointData.checkpointedData_,
335 DensityFittingModuleInfo::name_);
337 notifiers->checkpointingNotifier_.subscribe(checkpointDataReading);
339 // broadcasting checkpoint data
340 const auto checkpointDataBroadcast = [this](MDModulesCheckpointReadingBroadcast checkpointData) {
341 densityFittingState_.broadcastState(checkpointData.communicator_, checkpointData.isParallelRun_);
343 notifiers->checkpointingNotifier_.subscribe(checkpointDataBroadcast);
347 IMdpOptionProvider* mdpOptionProvider() override { return &densityFittingOptions_; }
349 //! Add this module to the force providers if active
350 void initForceProviders(ForceProviders* forceProviders) override
352 if (densityFittingOptions_.active())
354 const auto& parameters = densityFittingOptions_.buildParameters();
355 densityFittingSimulationParameters_.readReferenceDensityFromFile(
356 densityFittingOptions_.referenceDensityFileName());
357 if (parameters.normalizeDensities_)
359 densityFittingSimulationParameters_.normalizeReferenceDensity();
361 forceProvider_ = std::make_unique<DensityFittingForceProvider>(
363 densityFittingSimulationParameters_.referenceDensity(),
364 densityFittingSimulationParameters_.transformationToDensityLattice(),
365 densityFittingSimulationParameters_.localAtomSet(),
366 densityFittingSimulationParameters_.periodicBoundaryConditionType(),
367 densityFittingSimulationParameters_.simulationTimeStep(),
368 densityFittingState_);
369 forceProviders->addForceProvider(forceProvider_.get());
373 //! This MDModule provides its own output
374 IMDOutputProvider* outputProvider() override { return &densityFittingOutputProvider_; }
376 /*! \brief Set up the local atom sets that are used by this module.
378 * \note When density fitting is set up with MDModulesNotifier in
379 * the constructor, this function is called back.
381 * \param[in] localAtomSetManager the manager to add local atom sets.
383 void constructLocalAtomSet(LocalAtomSetManager* localAtomSetManager)
385 LocalAtomSet atomSet = localAtomSetManager->add(densityFittingOptions_.buildParameters().indices_);
386 densityFittingSimulationParameters_.setLocalAtomSet(atomSet);
389 /*! \brief Request energy output to energy file during simulation.
391 void setEnergyOutputRequest(MDModulesEnergyOutputToDensityFittingRequestChecker* energyOutputRequest)
393 energyOutputRequest->energyOutputToDensityFitting_ = densityFittingOptions_.active();
397 //! The output provider
398 DensityFittingOutputProvider densityFittingOutputProvider_;
399 //! The options provided for density fitting
400 DensityFittingOptions densityFittingOptions_;
401 //! Object that evaluates the forces
402 std::unique_ptr<DensityFittingForceProvider> forceProvider_;
403 /*! \brief Parameters for density fitting that become available at
404 * simulation setup time.
406 DensityFittingSimulationParameterSetup densityFittingSimulationParameters_;
407 //! The internal parameters of density fitting force provider
408 DensityFittingForceProviderState densityFittingState_;
410 GMX_DISALLOW_COPY_AND_ASSIGN(DensityFitting);
415 std::unique_ptr<IMDModule> DensityFittingModuleInfo::create()
417 return std::make_unique<DensityFitting>();
420 const std::string DensityFittingModuleInfo::name_ = "density-guided-simulation";