2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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/localatomsetmanager.h"
50 #include "gromacs/fileio/checkpoint.h"
51 #include "gromacs/fileio/mrcdensitymap.h"
52 #include "gromacs/math/multidimarray.h"
53 #include "gromacs/mdlib/broadcaststructs.h"
54 #include "gromacs/mdtypes/imdmodule.h"
55 #include "gromacs/utility/classhelpers.h"
56 #include "gromacs/utility/exceptions.h"
57 #include "gromacs/utility/keyvaluetreebuilder.h"
58 #include "gromacs/utility/mdmodulenotification.h"
60 #include "densityfittingforceprovider.h"
61 #include "densityfittingoptions.h"
62 #include "densityfittingoutputprovider.h"
67 class IMdpOptionProvider
;
68 class DensityFittingForceProvider
;
74 * \brief Collect density fitting parameters only available during simulation setup.
76 * \todo Implement builder pattern that will not use unqiue_ptr to check if
77 * parameters have been set or not.
79 * To build the density fitting force provider during simulation setup,
80 * the DensityFitting class needs access to parameters that become available
81 * only during simulation setup.
83 * This class collects these parameters via MdModuleNotifications in the
84 * simulation setup phase and provides a check if all necessary parameters have
87 class DensityFittingSimulationParameterSetup
90 DensityFittingSimulationParameterSetup() = default;
92 /*! \brief Set the local atom set for the density fitting.
93 * \param[in] localAtomSet of atoms to be fitted
95 void setLocalAtomSet(const LocalAtomSet
& localAtomSet
)
97 localAtomSet_
= std::make_unique
<LocalAtomSet
>(localAtomSet
);
100 /*! \brief Return local atom set for density fitting.
101 * \throws InternalError if local atom set is not set
102 * \returns local atom set for density fitting
104 const LocalAtomSet
& localAtomSet() const
106 if (localAtomSet_
== nullptr)
109 InternalError("Local atom set is not set for density "
110 "guided simulation."));
112 return *localAtomSet_
;
115 /*! \brief Return transformation into density lattice.
116 * \throws InternalError if transformation into density lattice is not set
117 * \returns transformation into density lattice
119 const TranslateAndScale
& transformationToDensityLattice() const
121 if (transformationToDensityLattice_
== nullptr)
123 GMX_THROW(InternalError(
124 "Transformation to reference density not set for density guided simulation."));
126 return *transformationToDensityLattice_
;
128 /*! \brief Return reference density
129 * \throws InternalError if reference density is not set
130 * \returns the reference density
132 basic_mdspan
<const float, dynamicExtents3D
> referenceDensity() const
134 if (referenceDensity_
== nullptr)
136 GMX_THROW(InternalError("Reference density not set for density guided simulation."));
138 return referenceDensity_
->asConstView();
141 /*! \brief Reads the reference density from file.
143 * Reads and check file, then set and communicate the internal
144 * parameters related to the reference density with the file data.
146 * \throws FileIOError if reading from file was not successful
148 void readReferenceDensityFromFile(const std::string
& referenceDensityFileName
)
150 MrcDensityMapOfFloatFromFileReader
reader(referenceDensityFileName
);
151 referenceDensity_
= std::make_unique
<MultiDimArray
<std::vector
<float>, dynamicExtents3D
>>(
152 reader
.densityDataCopy());
153 transformationToDensityLattice_
=
154 std::make_unique
<TranslateAndScale
>(reader
.transformationToDensityLattice());
157 //! Normalize the reference density so that the sum over all voxels is unity
158 void normalizeReferenceDensity()
160 if (referenceDensity_
== nullptr)
162 GMX_THROW(InternalError("Need to set reference density before normalizing it."));
165 const real sumOfDensityData
= std::accumulate(begin(referenceDensity_
->asView()),
166 end(referenceDensity_
->asView()), 0.);
167 for (float& referenceDensityVoxel
: referenceDensity_
->asView())
169 referenceDensityVoxel
/= sumOfDensityData
;
172 /*! \brief Set the periodic boundary condition via MdModuleNotifier.
174 * The pbc type is wrapped in PeriodicBoundaryConditionType to
175 * allow the MdModuleNotifier to statically distinguish the callback
176 * function type from other 'int' function callbacks.
178 * \param[in] pbcType enumerates the periodic boundary condition.
180 void setPeriodicBoundaryConditionType(const PbcType
& pbcType
)
182 pbcType_
= std::make_unique
<PbcType
>(pbcType
);
185 //! Get the periodic boundary conditions
186 PbcType
periodicBoundaryConditionType()
188 if (pbcType_
== nullptr)
190 GMX_THROW(InternalError(
191 "Periodic boundary condition enum not set for density guided simulation."));
196 //! Set the simulation time step
197 void setSimulationTimeStep(double timeStep
) { simulationTimeStep_
= timeStep
; }
199 //! Return the simulation time step
200 double simulationTimeStep() { return simulationTimeStep_
; }
203 //! The reference density to fit to
204 std::unique_ptr
<MultiDimArray
<std::vector
<float>, dynamicExtents3D
>> referenceDensity_
;
205 //! The coordinate transformation into the reference density
206 std::unique_ptr
<TranslateAndScale
> transformationToDensityLattice_
;
207 //! The local atom set to act on
208 std::unique_ptr
<LocalAtomSet
> localAtomSet_
;
209 //! The type of periodic boundary conditions in the simulation
210 std::unique_ptr
<PbcType
> pbcType_
;
211 //! The simulation time step
212 double simulationTimeStep_
= 1;
214 GMX_DISALLOW_COPY_AND_ASSIGN(DensityFittingSimulationParameterSetup
);
218 * \brief Density fitting
220 * Class that implements the density fitting forces and potential
221 * \note the virial calculation is not yet implemented
223 class DensityFitting final
: public IMDModule
227 //! Construct the density fitting module.
228 explicit DensityFitting() = default;
230 /*! \brief Request to be notified during pre-processing.
232 * \param[in] notifier allows the module to subscribe to notifications from MdModules.
234 * The density fitting code subscribes to these notifications:
235 * - setting atom group indices in the densityFittingOptions_ from an
236 * index group string by taking a parmeter const IndexGroupsAndNames &
237 * - storing its internal parameters in a tpr file by writing to a
238 * key-value-tree during pre-processing by a function taking a
239 * KeyValueTreeObjectBuilder as parameter
241 void subscribeToPreProcessingNotifications(MdModulesNotifier
* notifier
) override
243 // Callbacks for several kinds of MdModuleNotification are created
244 // and subscribed, and will be dispatched correctly at run time
245 // based on the type of the parameter required by the lambda.
247 if (!densityFittingOptions_
.active())
252 // Setting the atom group indices from index group string
253 const auto setFitGroupIndicesFunction
= [this](const IndexGroupsAndNames
& indexGroupsAndNames
) {
254 densityFittingOptions_
.setFitGroupIndices(indexGroupsAndNames
);
256 notifier
->preProcessingNotifications_
.subscribe(setFitGroupIndicesFunction
);
258 // Writing internal parameters during pre-processing
259 const auto writeInternalParametersFunction
= [this](KeyValueTreeObjectBuilder treeBuilder
) {
260 densityFittingOptions_
.writeInternalParametersToKvt(treeBuilder
);
262 notifier
->preProcessingNotifications_
.subscribe(writeInternalParametersFunction
);
264 // Checking for consistency with all .mdp options
265 const auto checkEnergyCaluclationFrequencyFunction
=
266 [this](EnergyCalculationFrequencyErrors
* energyCalculationFrequencyErrors
) {
267 densityFittingOptions_
.checkEnergyCaluclationFrequency(energyCalculationFrequencyErrors
);
269 notifier
->preProcessingNotifications_
.subscribe(checkEnergyCaluclationFrequencyFunction
);
272 /*! \brief Request to be notified.
273 * The density fitting code subscribes to these notifications:
274 * - reading its internal parameters from a key-value-tree during
275 * simulation setup by taking a const KeyValueTreeObject & parameter
276 * - constructing local atom sets in the simulation parameter setup
277 * by taking a LocalAtomSetManager * as parameter
278 * - the type of periodic boundary conditions that are used
279 * by taking a PeriodicBoundaryConditionType as parameter
280 * - the writing of checkpoint data
281 * by taking a MdModulesWriteCheckpointData as parameter
282 * - the reading of checkpoint data
283 * by taking a MdModulesCheckpointReadingDataOnMaster as parameter
284 * - the broadcasting of checkpoint data
285 * by taking MdModulesCheckpointReadingBroadcast as parameter
287 void subscribeToSimulationSetupNotifications(MdModulesNotifier
* notifier
) override
289 if (!densityFittingOptions_
.active())
294 // Reading internal parameters during simulation setup
295 const auto readInternalParametersFunction
= [this](const KeyValueTreeObject
& tree
) {
296 densityFittingOptions_
.readInternalParametersFromKvt(tree
);
298 notifier
->simulationSetupNotifications_
.subscribe(readInternalParametersFunction
);
300 // constructing local atom sets during simulation setup
301 const auto setLocalAtomSetFunction
= [this](LocalAtomSetManager
* localAtomSetManager
) {
302 this->constructLocalAtomSet(localAtomSetManager
);
304 notifier
->simulationSetupNotifications_
.subscribe(setLocalAtomSetFunction
);
306 // constructing local atom sets during simulation setup
307 const auto setPeriodicBoundaryContionsFunction
= [this](const PbcType
& pbc
) {
308 this->densityFittingSimulationParameters_
.setPeriodicBoundaryConditionType(pbc
);
310 notifier
->simulationSetupNotifications_
.subscribe(setPeriodicBoundaryContionsFunction
);
312 // setting the simulation time step
313 const auto setSimulationTimeStepFunction
= [this](const SimulationTimeStep
& simulationTimeStep
) {
314 this->densityFittingSimulationParameters_
.setSimulationTimeStep(simulationTimeStep
.delta_t
);
316 notifier
->simulationSetupNotifications_
.subscribe(setSimulationTimeStepFunction
);
318 // adding output to energy file
319 const auto requestEnergyOutput
=
320 [this](MdModulesEnergyOutputToDensityFittingRequestChecker
* energyOutputRequest
) {
321 this->setEnergyOutputRequest(energyOutputRequest
);
323 notifier
->simulationSetupNotifications_
.subscribe(requestEnergyOutput
);
325 // writing checkpoint data
326 const auto checkpointDataWriting
= [this](MdModulesWriteCheckpointData checkpointData
) {
327 this->writeCheckpointData(checkpointData
);
329 notifier
->checkpointingNotifications_
.subscribe(checkpointDataWriting
);
331 // reading checkpoint data
332 const auto checkpointDataReading
= [this](MdModulesCheckpointReadingDataOnMaster checkpointData
) {
333 this->readCheckpointDataOnMaster(checkpointData
);
335 notifier
->checkpointingNotifications_
.subscribe(checkpointDataReading
);
337 // broadcasting checkpoint data
338 const auto checkpointDataBroadcast
= [this](MdModulesCheckpointReadingBroadcast checkpointData
) {
339 this->broadcastCheckpointData(checkpointData
);
341 notifier
->checkpointingNotifications_
.subscribe(checkpointDataBroadcast
);
345 IMdpOptionProvider
* mdpOptionProvider() override
{ return &densityFittingOptions_
; }
347 //! Add this module to the force providers if active
348 void initForceProviders(ForceProviders
* forceProviders
) override
350 if (densityFittingOptions_
.active())
352 const auto& parameters
= densityFittingOptions_
.buildParameters();
353 densityFittingSimulationParameters_
.readReferenceDensityFromFile(
354 densityFittingOptions_
.referenceDensityFileName());
355 if (parameters
.normalizeDensities_
)
357 densityFittingSimulationParameters_
.normalizeReferenceDensity();
359 forceProvider_
= std::make_unique
<DensityFittingForceProvider
>(
360 parameters
, densityFittingSimulationParameters_
.referenceDensity(),
361 densityFittingSimulationParameters_
.transformationToDensityLattice(),
362 densityFittingSimulationParameters_
.localAtomSet(),
363 densityFittingSimulationParameters_
.periodicBoundaryConditionType(),
364 densityFittingSimulationParameters_
.simulationTimeStep(), densityFittingState_
);
365 forceProviders
->addForceProvider(forceProvider_
.get());
369 //! This MDModule provides its own output
370 IMDOutputProvider
* outputProvider() override
{ return &densityFittingOutputProvider_
; }
372 /*! \brief Set up the local atom sets that are used by this module.
374 * \note When density fitting is set up with MdModuleNotification in
375 * the constructor, this function is called back.
377 * \param[in] localAtomSetManager the manager to add local atom sets.
379 void constructLocalAtomSet(LocalAtomSetManager
* localAtomSetManager
)
381 LocalAtomSet atomSet
= localAtomSetManager
->add(densityFittingOptions_
.buildParameters().indices_
);
382 densityFittingSimulationParameters_
.setLocalAtomSet(atomSet
);
385 /*! \brief Request energy output to energy file during simulation.
387 void setEnergyOutputRequest(MdModulesEnergyOutputToDensityFittingRequestChecker
* energyOutputRequest
)
389 energyOutputRequest
->energyOutputToDensityFitting_
= densityFittingOptions_
.active();
392 /*! \brief Write internal density fitting data to checkpoint file.
393 * \param[in] checkpointWriting enables writing to the Key-Value-Tree
394 * that is used for storing the checkpoint
397 * \note The provided state to checkpoint has to change if checkpointing
398 * is moved before the force provider call in the MD-loop.
400 void writeCheckpointData(MdModulesWriteCheckpointData checkpointWriting
)
402 const DensityFittingForceProviderState
& state
= forceProvider_
->stateToCheckpoint();
403 checkpointWriting
.builder_
.addValue
<std::int64_t>(
404 DensityFittingModuleInfo::name_
+ "-stepsSinceLastCalculation",
405 state
.stepsSinceLastCalculation_
);
406 checkpointWriting
.builder_
.addValue
<real
>(
407 DensityFittingModuleInfo::name_
+ "-adaptiveForceConstantScale",
408 state
.adaptiveForceConstantScale_
);
409 KeyValueTreeObjectBuilder exponentialMovingAverageKvtEntry
= checkpointWriting
.builder_
.addObject(
410 DensityFittingModuleInfo::name_
+ "-exponentialMovingAverageState");
411 exponentialMovingAverageStateAsKeyValueTree(exponentialMovingAverageKvtEntry
,
412 state
.exponentialMovingAverageState_
);
415 /*! \brief Read the internal parameters from the checkpoint file on master
416 * \param[in] checkpointReading holding the checkpoint information
418 void readCheckpointDataOnMaster(MdModulesCheckpointReadingDataOnMaster checkpointReading
)
420 if (checkpointReading
.checkpointedData_
.keyExists(DensityFittingModuleInfo::name_
421 + "-stepsSinceLastCalculation"))
423 densityFittingState_
.stepsSinceLastCalculation_
=
425 .checkpointedData_
[DensityFittingModuleInfo::name_
426 + "-stepsSinceLastCalculation"]
427 .cast
<std::int64_t>();
429 if (checkpointReading
.checkpointedData_
.keyExists(DensityFittingModuleInfo::name_
430 + "-adaptiveForceConstantScale"))
432 densityFittingState_
.adaptiveForceConstantScale_
=
434 .checkpointedData_
[DensityFittingModuleInfo::name_
435 + "-adaptiveForceConstantScale"]
438 if (checkpointReading
.checkpointedData_
.keyExists(DensityFittingModuleInfo::name_
439 + "-exponentialMovingAverageState"))
441 densityFittingState_
.exponentialMovingAverageState_
= exponentialMovingAverageStateFromKeyValueTree(
443 .checkpointedData_
[DensityFittingModuleInfo::name_
+ "-exponentialMovingAverageState"]
448 /*! \brief Broadcast the internal parameters.
449 * \param[in] checkpointBroadcast containing the communication record to
450 * broadcast the checkpoint information
452 void broadcastCheckpointData(MdModulesCheckpointReadingBroadcast checkpointBroadcast
)
454 if (checkpointBroadcast
.isParallelRun_
)
456 block_bc(checkpointBroadcast
.communicator_
, densityFittingState_
.stepsSinceLastCalculation_
);
457 block_bc(checkpointBroadcast
.communicator_
, densityFittingState_
.adaptiveForceConstantScale_
);
458 block_bc(checkpointBroadcast
.communicator_
, densityFittingState_
.exponentialMovingAverageState_
);
463 //! The output provider
464 DensityFittingOutputProvider densityFittingOutputProvider_
;
465 //! The options provided for density fitting
466 DensityFittingOptions densityFittingOptions_
;
467 //! Object that evaluates the forces
468 std::unique_ptr
<DensityFittingForceProvider
> forceProvider_
;
469 /*! \brief Parameters for density fitting that become available at
470 * simulation setup time.
472 DensityFittingSimulationParameterSetup densityFittingSimulationParameters_
;
473 //! The internal parameters of density fitting force provider
474 DensityFittingForceProviderState densityFittingState_
;
476 GMX_DISALLOW_COPY_AND_ASSIGN(DensityFitting
);
481 std::unique_ptr
<IMDModule
> DensityFittingModuleInfo::create()
483 return std::make_unique
<DensityFitting
>();
486 const std::string
DensityFittingModuleInfo::name_
= "density-guided-simulation";