Update energy and trajectory frame types
[gromacs.git] / src / gromacs / awh / awh.cpp
blobe03bea12084354c01e9f46733073ef07d38800d4
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,2018, 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.
36 /*! \internal \file
37 * \brief
38 * Implements the Awh class.
40 * \author Viveca Lindahl
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_awh
45 #include "gmxpre.h"
47 #include "awh.h"
49 #include <assert.h>
51 #include <cmath>
52 #include <cstdio>
53 #include <cstdlib>
54 #include <cstring>
56 #include <algorithm>
58 #include "gromacs/fileio/enxio.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/math/units.h"
61 #include "gromacs/mdtypes/awh-history.h"
62 #include "gromacs/mdtypes/awh-params.h"
63 #include "gromacs/mdtypes/commrec.h"
64 #include "gromacs/mdtypes/forceoutput.h"
65 #include "gromacs/mdtypes/inputrec.h"
66 #include "gromacs/mdtypes/pull-params.h"
67 #include "gromacs/mdtypes/state.h"
68 #include "gromacs/pbcutil/pbc.h"
69 #include "gromacs/pulling/pull.h"
70 #include "gromacs/timing/wallcycle.h"
71 #include "gromacs/trajectory/energyframe.h"
72 #include "gromacs/utility/exceptions.h"
73 #include "gromacs/utility/gmxassert.h"
74 #include "gromacs/utility/pleasecite.h"
76 #include "bias.h"
77 #include "biassharing.h"
78 #include "correlationgrid.h"
79 #include "pointstate.h"
81 namespace gmx
84 /*! \internal
85 * \brief A bias and its coupling to the system.
87 * This struct is used to separate the bias machinery in the Bias class,
88 * which should be independent from the reaction coordinate, from the
89 * obtaining of the reaction coordinate values and passing the computed forces.
90 * Currently the AWH method couples to the system by mapping each
91 * AWH bias to a pull coordinate. This can easily be generalized here.
93 struct BiasCoupledToSystem
95 /*! \brief Constructor, couple a bias to a set of pull coordinates.
97 * \param[in] bias The bias.
98 * \param[in] pullCoordIndex The pull coordinate indices.
100 BiasCoupledToSystem(Bias bias,
101 const std::vector<int> &pullCoordIndex);
103 Bias bias; /**< The bias. */
104 const std::vector<int> pullCoordIndex; /**< The pull coordinates this bias acts on. */
106 /* Here AWH can be extended to work on other coordinates than pull. */
109 BiasCoupledToSystem::BiasCoupledToSystem(Bias bias,
110 const std::vector<int> &pullCoordIndex) :
111 bias(std::move(bias)),
112 pullCoordIndex(pullCoordIndex)
114 /* We already checked for this in grompp, but check again here. */
115 GMX_RELEASE_ASSERT(static_cast<size_t>(bias.ndim()) == pullCoordIndex.size(), "The bias dimensionality should match the number of pull coordinates.");
118 Awh::Awh(FILE *fplog,
119 const t_inputrec &inputRecord,
120 const t_commrec *commRecord,
121 const gmx_multisim_t *multiSimRecord,
122 const AwhParams &awhParams,
123 const std::string &biasInitFilename,
124 pull_t *pull_work) :
125 seed_(awhParams.seed),
126 nstout_(awhParams.nstOut),
127 commRecord_(commRecord),
128 multiSimRecord_(multiSimRecord),
129 pull_(pull_work),
130 potentialOffset_(0)
132 /* We already checked for this in grompp, but check again here. */
133 GMX_RELEASE_ASSERT(inputRecord.pull != nullptr, "With AWH we should have pull parameters");
134 GMX_RELEASE_ASSERT(pull_work != nullptr, "With AWH pull should be initialized before initializing AWH");
136 if (fplog != nullptr)
138 please_cite(fplog, "Lindahl2014");
141 if (haveBiasSharingWithinSimulation(awhParams))
143 /* This has likely been checked by grompp, but throw anyhow. */
144 GMX_THROW(InvalidInputError("Biases within a simulation are shared, currently sharing of biases is only supported between simulations"));
147 int numSharingSimulations = 1;
148 if (awhParams.shareBiasMultisim && isMultiSim(multiSimRecord_))
150 numSharingSimulations = multiSimRecord_->nsim;
153 /* Initialize all the biases */
154 const double beta = 1/(BOLTZ*inputRecord.opts.ref_t[0]);
155 for (int k = 0; k < awhParams.numBias; k++)
157 const AwhBiasParams &awhBiasParams = awhParams.awhBiasParams[k];
159 std::vector<int> pullCoordIndex;
160 std::vector<DimParams> dimParams;
161 for (int d = 0; d < awhBiasParams.ndim; d++)
163 const AwhDimParams &awhDimParams = awhBiasParams.dimParams[d];
164 GMX_RELEASE_ASSERT(awhDimParams.eCoordProvider == eawhcoordproviderPULL, "Currently only the pull code is supported as coordinate provider");
165 const t_pull_coord &pullCoord = inputRecord.pull->coord[awhDimParams.coordIndex];
166 double conversionFactor = pull_coordinate_is_angletype(&pullCoord) ? DEG2RAD : 1;
167 dimParams.push_back(DimParams(conversionFactor, awhDimParams.forceConstant, beta));
169 pullCoordIndex.push_back(awhDimParams.coordIndex);
172 /* Construct the bias and couple it to the system. */
173 Bias::ThisRankWillDoIO thisRankWillDoIO = (MASTER(commRecord_) ? Bias::ThisRankWillDoIO::Yes : Bias::ThisRankWillDoIO::No);
174 biasCoupledToSystem_.emplace_back(Bias(k, awhParams, awhParams.awhBiasParams[k], dimParams, beta, inputRecord.delta_t, numSharingSimulations, biasInitFilename, thisRankWillDoIO),
175 pullCoordIndex);
177 biasCoupledToSystem_.back().bias.printInitializationToLog(fplog);
180 /* Need to register the AWH coordinates to be allowed to apply forces to the pull coordinates. */
181 registerAwhWithPull(awhParams, pull_);
183 if (numSharingSimulations > 1 && MASTER(commRecord_))
185 std::vector<size_t> pointSize;
186 for (auto const &biasCts : biasCoupledToSystem_)
188 pointSize.push_back(biasCts.bias.state().points().size());
190 /* Ensure that the shared biased are compatible between simulations */
191 biasesAreCompatibleForSharingBetweenSimulations(awhParams, pointSize, multiSimRecord_);
195 Awh::~Awh() = default;
197 bool Awh::isOutputStep(gmx_int64_t step) const
199 return (nstout_ > 0 && step % nstout_ == 0);
202 real Awh::applyBiasForcesAndUpdateBias(int ePBC,
203 const t_mdatoms &mdatoms,
204 const matrix box,
205 gmx::ForceWithVirial *forceWithVirial,
206 double t,
207 gmx_int64_t step,
208 gmx_wallcycle *wallcycle,
209 FILE *fplog)
211 GMX_ASSERT(forceWithVirial, "Need a valid ForceWithVirial object");
213 wallcycle_start(wallcycle, ewcAWH);
215 t_pbc pbc;
216 set_pbc(&pbc, ePBC, box);
218 /* During the AWH update the potential can instantaneously jump due to either
219 an bias update or moving the umbrella. The jumps are kept track of and
220 subtracted from the potential in order to get a useful conserved energy quantity. */
221 double awhPotential = potentialOffset_;
223 for (auto &biasCts : biasCoupledToSystem_)
225 /* Update the AWH coordinate values with those of the corresponding
226 * pull coordinates.
228 awh_dvec coordValue = { 0, 0, 0, 0 };
229 for (int d = 0; d < biasCts.bias.ndim(); d++)
231 coordValue[d] = get_pull_coord_value(pull_, biasCts.pullCoordIndex[d], &pbc);
234 /* Perform an AWH biasing step: this means, at regular intervals,
235 * sampling observables based on the input pull coordinate value,
236 * setting the bias force and/or updating the AWH bias state.
238 double biasPotential;
239 double biasPotentialJump;
240 /* Note: In the near future this call will be split in calls
241 * to supports bias sharing within a single simulation.
243 gmx::ArrayRef<const double> biasForce =
244 biasCts.bias.calcForceAndUpdateBias(coordValue,
245 &biasPotential, &biasPotentialJump,
246 commRecord_,
247 multiSimRecord_,
248 t, step, seed_, fplog);
250 awhPotential += biasPotential;
252 /* Keep track of the total potential shift needed to remove the potential jumps. */
253 potentialOffset_ -= biasPotentialJump;
255 /* Communicate the bias force to the pull struct.
256 * The bias potential is returned at the end of this function,
257 * so that it can be added externally to the correct energy data block.
259 for (int d = 0; d < biasCts.bias.ndim(); d++)
261 apply_external_pull_coord_force(pull_, biasCts.pullCoordIndex[d],
262 biasForce[d], &mdatoms,
263 forceWithVirial);
266 if (isOutputStep(step))
268 /* We might have skipped updates for part of the grid points.
269 * Ensure all points are updated before writing out their data.
271 biasCts.bias.doSkippedUpdatesForAllPoints();
275 wallcycle_stop(wallcycle, ewcAWH);
277 return MASTER(commRecord_) ? static_cast<real>(awhPotential) : 0;
280 std::shared_ptr<AwhHistory> Awh::initHistoryFromState() const
282 if (MASTER(commRecord_))
284 std::shared_ptr<AwhHistory> awhHistory(new AwhHistory);
285 awhHistory->bias.clear();
286 awhHistory->bias.resize(biasCoupledToSystem_.size());
288 for (size_t k = 0; k < awhHistory->bias.size(); k++)
290 biasCoupledToSystem_[k].bias.initHistoryFromState(&awhHistory->bias[k]);
293 return awhHistory;
295 else
297 /* Return an empty pointer */
298 return std::shared_ptr<AwhHistory>();
302 void Awh::restoreStateFromHistory(const AwhHistory *awhHistory)
304 /* Restore the history to the current state */
305 if (MASTER(commRecord_))
307 GMX_RELEASE_ASSERT(awhHistory != nullptr, "The master rank should have a valid awhHistory when restoring the state from history.");
309 if (awhHistory->bias.size() != biasCoupledToSystem_.size())
311 GMX_THROW(InvalidInputError("AWH state and history contain different numbers of biases. Likely you provided a checkpoint from a different simulation."));
314 potentialOffset_ = awhHistory->potentialOffset;
316 if (PAR(commRecord_))
318 gmx_bcast(sizeof(potentialOffset_), &potentialOffset_, commRecord_);
321 for (size_t k = 0; k < biasCoupledToSystem_.size(); k++)
323 biasCoupledToSystem_[k].bias.restoreStateFromHistory(awhHistory ? &awhHistory->bias[k] : nullptr, commRecord_);
327 void Awh::updateHistory(AwhHistory *awhHistory) const
329 if (!MASTER(commRecord_))
331 return;
334 /* This assert will also catch a non-master rank calling this function. */
335 GMX_RELEASE_ASSERT(awhHistory->bias.size() == biasCoupledToSystem_.size(), "AWH state and history bias count should match");
337 awhHistory->potentialOffset = potentialOffset_;
339 for (size_t k = 0; k < awhHistory->bias.size(); k++)
341 biasCoupledToSystem_[k].bias.updateHistory(&awhHistory->bias[k]);
345 const char * Awh::externalPotentialString()
347 return "AWH";
350 void Awh::registerAwhWithPull(const AwhParams &awhParams,
351 pull_t *pull_work)
353 GMX_RELEASE_ASSERT(pull_work, "Need a valid pull object");
355 for (int k = 0; k < awhParams.numBias; k++)
357 const AwhBiasParams &biasParams = awhParams.awhBiasParams[k];
359 for (int d = 0; d < biasParams.ndim; d++)
361 register_external_pull_potential(pull_work, biasParams.dimParams[d].coordIndex, Awh::externalPotentialString());
366 /* Fill the AWH data block of an energy frame with data (if there is any). */
367 void Awh::writeToEnergyFrame(gmx_int64_t step,
368 t_enxframe *frame) const
370 GMX_ASSERT(MASTER(commRecord_), "writeToEnergyFrame should only be called on the master rank");
371 GMX_ASSERT(frame != nullptr, "Need a valid energy frame");
373 if (!isOutputStep(step))
375 /* This is not an AWH output step, don't write any AWH data */
376 return;
379 /* Get the total number of energy subblocks that AWH needs */
380 int numSubblocks = 0;
381 for (auto &biasCoupledToSystem : biasCoupledToSystem_)
383 numSubblocks += biasCoupledToSystem.bias.numEnergySubblocksToWrite();
385 GMX_ASSERT(numSubblocks > 0, "We should always have data to write");
387 /* Add 1 energy block */
388 add_blocks_enxframe(frame, frame->nblock + 1);
390 /* Take the block that was just added and set the number of subblocks. */
391 t_enxblock *awhEnergyBlock = &(frame->block[frame->nblock - 1]);
392 add_subblocks_enxblock(awhEnergyBlock, numSubblocks);
394 /* Claim it as an AWH block. */
395 awhEnergyBlock->id = enxAWH;
397 /* Transfer AWH data blocks to energy sub blocks */
398 int energySubblockCount = 0;
399 for (auto &biasCoupledToSystem : biasCoupledToSystem_)
401 energySubblockCount += biasCoupledToSystem.bias.writeToEnergySubblocks(&(awhEnergyBlock->sub[energySubblockCount]));
405 std::unique_ptr<Awh>
406 prepareAwhModule(FILE *fplog,
407 const t_inputrec &inputRecord,
408 t_state *stateGlobal,
409 const t_commrec *commRecord,
410 const gmx_multisim_t *multiSimRecord,
411 const bool startingFromCheckpoint,
412 const bool usingShellParticles,
413 const std::string &biasInitFilename,
414 pull_t *pull_work)
416 if (!inputRecord.bDoAwh)
418 return nullptr;
420 if (usingShellParticles)
422 GMX_THROW(InvalidInputError("AWH biasing does not support shell particles."));
425 auto awh = compat::make_unique<Awh>(fplog, inputRecord, commRecord, multiSimRecord, *inputRecord.awhParams,
426 biasInitFilename, pull_work);
428 if (startingFromCheckpoint)
430 // Restore the AWH history read from checkpoint
431 awh->restoreStateFromHistory(MASTER(commRecord) ? stateGlobal->awhHistory.get() : nullptr);
433 else if (MASTER(commRecord))
435 // Initialize the AWH history here
436 stateGlobal->awhHistory = awh->initHistoryFromState();
438 return awh;
441 } // namespace gmx