0f3c8c28b4c55c6e5436eaebb0bde15c06e469b5
[alexxy/gromacs.git] / src / gromacs / applied_forces / awh / awh.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2015,2016,2017,2018,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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  * \brief
38  * Implements the Awh class.
39  *
40  * \author Viveca Lindahl
41  * \author Berk Hess <hess@kth.se>
42  * \author Magnus Lundborg
43  * \ingroup module_awh
44  */
45
46 #include "gmxpre.h"
47
48 #include "awh.h"
49
50 #include <cassert>
51 #include <cmath>
52 #include <cstdio>
53 #include <cstdlib>
54 #include <cstring>
55
56 #include <algorithm>
57
58 #include "gromacs/fileio/enxio.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/math/units.h"
61 #include "gromacs/mdrunutility/multisim.h"
62 #include "gromacs/mdtypes/awh_history.h"
63 #include "gromacs/mdtypes/awh_params.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/forceoutput.h"
66 #include "gromacs/mdtypes/inputrec.h"
67 #include "gromacs/mdtypes/pull_params.h"
68 #include "gromacs/mdtypes/state.h"
69 #include "gromacs/pbcutil/pbc.h"
70 #include "gromacs/pulling/pull.h"
71 #include "gromacs/timing/wallcycle.h"
72 #include "gromacs/trajectory/energyframe.h"
73 #include "gromacs/utility/arrayref.h"
74 #include "gromacs/utility/exceptions.h"
75 #include "gromacs/utility/gmxassert.h"
76 #include "gromacs/utility/pleasecite.h"
77
78 #include "bias.h"
79 #include "biassharing.h"
80 #include "correlationgrid.h"
81 #include "pointstate.h"
82
83 namespace gmx
84 {
85
86 /*! \internal
87  * \brief A bias and its coupling to the system.
88  *
89  * This struct is used to separate the bias machinery in the Bias class,
90  * which should be independent from the reaction coordinate, from the
91  * obtaining of the reaction coordinate values and passing the computed forces.
92  * Currently the AWH method couples to the system by mapping each
93  * AWH bias to a pull coordinate. This can easily be generalized here.
94  */
95 struct BiasCoupledToSystem
96 {
97     /*! \brief Constructor, couple a bias to a set of pull coordinates.
98      *
99      * \param[in] bias            The bias.
100      * \param[in] pullCoordIndex  The pull coordinate indices.
101      */
102     BiasCoupledToSystem(Bias bias, const std::vector<int>& pullCoordIndex);
103
104     Bias                   bias_;           /**< The bias. */
105     const std::vector<int> pullCoordIndex_; /**< The pull coordinates this bias acts on. */
106
107     /* Here AWH can be extended to work on other coordinates than pull. */
108 };
109
110 /*! \brief Checks whether any dimension uses pulling as a coordinate provider.
111  *
112  * \param[in] awhBiasParams The bias params to check.
113  * \returns true if any dimension of the bias is linked to pulling.
114  */
115 static bool anyDimUsesPull(const AwhBiasParams& awhBiasParams)
116 {
117     for (int d = 0; d < awhBiasParams.ndim; d++)
118     {
119         const AwhDimParams& awhDimParams = awhBiasParams.dimParams[d];
120         if (awhDimParams.eCoordProvider == eawhcoordproviderPULL)
121         {
122             return true;
123         }
124     }
125     return false;
126 }
127
128 /*! \brief Checks whether any dimension uses pulling as a coordinate provider.
129  *
130  * \param[in] awhParams The AWH params to check.
131  * \returns true if any dimension of awh is linked to pulling.
132  */
133 static bool anyDimUsesPull(const AwhParams& awhParams)
134 {
135     for (int k = 0; k < awhParams.numBias; k++)
136     {
137         const AwhBiasParams& awhBiasParams = awhParams.awhBiasParams[k];
138         if (anyDimUsesPull(awhBiasParams))
139         {
140             return true;
141         }
142     }
143     return false;
144 }
145
146 /*! \brief Checks whether any dimension uses pulling as a coordinate provider.
147  *
148  * \param[in] biasCoupledToSystem The AWH biases to check.
149  * \returns true if any dimension of the provided biases is linked to pulling.
150  */
151 static bool anyDimUsesPull(const ArrayRef<BiasCoupledToSystem> biasCoupledToSystem)
152 {
153     for (const auto& biasCts : biasCoupledToSystem)
154     {
155         if (!biasCts.pullCoordIndex_.empty())
156         {
157             return true;
158         }
159     }
160     return false;
161 }
162
163 BiasCoupledToSystem::BiasCoupledToSystem(Bias bias, const std::vector<int>& pullCoordIndex) :
164     bias_(std::move(bias)),
165     pullCoordIndex_(pullCoordIndex)
166 {
167     /* We already checked for this in grompp, but check again here. */
168     GMX_RELEASE_ASSERT(
169             static_cast<size_t>(bias_.ndim()) == pullCoordIndex_.size() + bias_.hasFepLambdaDimension() ? 1 : 0,
170             "The bias dimensionality should match the number of pull and lambda coordinates.");
171 }
172
173 Awh::Awh(FILE*                 fplog,
174          const t_inputrec&     inputRecord,
175          const t_commrec*      commRecord,
176          const gmx_multisim_t* multiSimRecord,
177          const AwhParams&      awhParams,
178          const std::string&    biasInitFilename,
179          pull_t*               pull_work,
180          int                   numFepLambdaStates,
181          int                   fepLambdaState) :
182     seed_(awhParams.seed),
183     nstout_(awhParams.nstOut),
184     commRecord_(commRecord),
185     multiSimRecord_(multiSimRecord),
186     pull_(pull_work),
187     potentialOffset_(0),
188     numFepLambdaStates_(numFepLambdaStates),
189     fepLambdaState_(fepLambdaState)
190 {
191     if (anyDimUsesPull(awhParams))
192     {
193         GMX_RELEASE_ASSERT(inputRecord.pull != nullptr, "With AWH we should have pull parameters");
194         GMX_RELEASE_ASSERT(pull_work != nullptr,
195                            "With AWH pull should be initialized before initializing AWH");
196     }
197
198     if (fplog != nullptr)
199     {
200         please_cite(fplog, "Lindahl2014");
201     }
202
203     if (haveBiasSharingWithinSimulation(awhParams))
204     {
205         /* This has likely been checked by grompp, but throw anyhow. */
206         GMX_THROW(
207                 InvalidInputError("Biases within a simulation are shared, currently sharing of "
208                                   "biases is only supported between simulations"));
209     }
210
211     int numSharingSimulations = 1;
212     if (awhParams.shareBiasMultisim && isMultiSim(multiSimRecord_))
213     {
214         numSharingSimulations = multiSimRecord_->numSimulations_;
215     }
216
217     /* Initialize all the biases */
218     const double beta = 1 / (BOLTZ * inputRecord.opts.ref_t[0]);
219     for (int k = 0; k < awhParams.numBias; k++)
220     {
221         const AwhBiasParams& awhBiasParams = awhParams.awhBiasParams[k];
222
223         std::vector<int>       pullCoordIndex;
224         std::vector<DimParams> dimParams;
225         for (int d = 0; d < awhBiasParams.ndim; d++)
226         {
227             const AwhDimParams& awhDimParams = awhBiasParams.dimParams[d];
228             if (awhDimParams.eCoordProvider != eawhcoordproviderPULL
229                 && awhDimParams.eCoordProvider != eawhcoordproviderFREE_ENERGY_LAMBDA)
230             {
231                 GMX_THROW(
232                         InvalidInputError("Currently only the pull code and lambda are supported "
233                                           "as coordinate providers"));
234             }
235             if (awhDimParams.eCoordProvider == eawhcoordproviderPULL)
236             {
237                 const t_pull_coord& pullCoord = inputRecord.pull->coord[awhDimParams.coordIndex];
238                 if (pullCoord.eGeom == epullgDIRPBC)
239                 {
240                     GMX_THROW(InvalidInputError(
241                             "Pull geometry 'direction-periodic' is not supported by AWH"));
242                 }
243                 double conversionFactor = pull_coordinate_is_angletype(&pullCoord) ? DEG2RAD : 1;
244                 pullCoordIndex.push_back(awhDimParams.coordIndex);
245                 dimParams.push_back(DimParams::pullDimParams(conversionFactor,
246                                                              awhDimParams.forceConstant, beta));
247             }
248             else
249             {
250                 dimParams.push_back(DimParams::fepLambdaDimParams(numFepLambdaStates_, beta));
251             }
252         }
253
254         /* Construct the bias and couple it to the system. */
255         Bias::ThisRankWillDoIO thisRankWillDoIO =
256                 (MASTER(commRecord_) ? Bias::ThisRankWillDoIO::Yes : Bias::ThisRankWillDoIO::No);
257         biasCoupledToSystem_.emplace_back(
258                 Bias(k, awhParams, awhParams.awhBiasParams[k], dimParams, beta, inputRecord.delta_t,
259                      numSharingSimulations, biasInitFilename, thisRankWillDoIO),
260                 pullCoordIndex);
261
262         biasCoupledToSystem_.back().bias_.printInitializationToLog(fplog);
263     }
264
265     /* Need to register the AWH coordinates to be allowed to apply forces to the pull coordinates. */
266     registerAwhWithPull(awhParams, pull_);
267
268     if (numSharingSimulations > 1 && MASTER(commRecord_))
269     {
270         std::vector<size_t> pointSize;
271         for (auto const& biasCts : biasCoupledToSystem_)
272         {
273             pointSize.push_back(biasCts.bias_.state().points().size());
274         }
275         /* Ensure that the shared biased are compatible between simulations */
276         biasesAreCompatibleForSharingBetweenSimulations(awhParams, pointSize, multiSimRecord_);
277     }
278 }
279
280 Awh::~Awh() = default;
281
282 bool Awh::isOutputStep(int64_t step) const
283 {
284     return (nstout_ > 0 && step % nstout_ == 0);
285 }
286
287 real Awh::applyBiasForcesAndUpdateBias(PbcType                pbcType,
288                                        const real*            masses,
289                                        ArrayRef<const double> neighborLambdaEnergies,
290                                        ArrayRef<const double> neighborLambdaDhdl,
291                                        const matrix           box,
292                                        gmx::ForceWithVirial*  forceWithVirial,
293                                        double                 t,
294                                        int64_t                step,
295                                        gmx_wallcycle*         wallcycle,
296                                        FILE*                  fplog)
297 {
298     if (anyDimUsesPull(biasCoupledToSystem_))
299     {
300         GMX_ASSERT(forceWithVirial, "Need a valid ForceWithVirial object");
301     }
302
303     wallcycle_start(wallcycle, ewcAWH);
304
305     t_pbc pbc;
306     set_pbc(&pbc, pbcType, box);
307
308     /* During the AWH update the potential can instantaneously jump due to either
309        an bias update or moving the umbrella. The jumps are kept track of and
310        subtracted from the potential in order to get a useful conserved energy quantity. */
311     double awhPotential = potentialOffset_;
312
313     for (auto& biasCts : biasCoupledToSystem_)
314     {
315         /* Update the AWH coordinate values with those of the corresponding
316          * pull coordinates.
317          */
318         awh_dvec coordValue           = { 0, 0, 0, 0 };
319         int      numLambdaDimsCounted = 0;
320         for (int d = 0; d < biasCts.bias_.ndim(); d++)
321         {
322             if (biasCts.bias_.dimParams()[d].isPullDimension())
323             {
324                 coordValue[d] = get_pull_coord_value(
325                         pull_, biasCts.pullCoordIndex_[d - numLambdaDimsCounted], &pbc);
326             }
327             else
328             {
329                 coordValue[d] = fepLambdaState_;
330                 numLambdaDimsCounted += 1;
331             }
332         }
333
334         /* Perform an AWH biasing step: this means, at regular intervals,
335          * sampling observables based on the input pull coordinate value,
336          * setting the bias force and/or updating the AWH bias state.
337          */
338         double biasPotential;
339         double biasPotentialJump;
340         /* Note: In the near future this call will be split in calls
341          *       to supports bias sharing within a single simulation.
342          */
343         gmx::ArrayRef<const double> biasForce = biasCts.bias_.calcForceAndUpdateBias(
344                 coordValue, neighborLambdaEnergies, neighborLambdaDhdl, &biasPotential,
345                 &biasPotentialJump, commRecord_, multiSimRecord_, t, step, seed_, fplog);
346
347         awhPotential += biasPotential;
348
349         /* Keep track of the total potential shift needed to remove the potential jumps. */
350         potentialOffset_ -= biasPotentialJump;
351
352         /* Communicate the bias force to the pull struct.
353          * The bias potential is returned at the end of this function,
354          * so that it can be added externally to the correct energy data block.
355          */
356         numLambdaDimsCounted = 0;
357         for (int d = 0; d < biasCts.bias_.ndim(); d++)
358         {
359             if (biasCts.bias_.dimParams()[d].isPullDimension())
360             {
361                 apply_external_pull_coord_force(pull_, biasCts.pullCoordIndex_[d - numLambdaDimsCounted],
362                                                 biasForce[d], masses, forceWithVirial);
363             }
364             else
365             {
366                 int umbrellaGridpointIndex = biasCts.bias_.state().coordState().umbrellaGridpoint();
367                 fepLambdaState_ = biasCts.bias_.getGridCoordValue(umbrellaGridpointIndex)[d];
368                 numLambdaDimsCounted += 1;
369             }
370         }
371
372         if (isOutputStep(step))
373         {
374             /* We might have skipped updates for part of the grid points.
375              * Ensure all points are updated before writing out their data.
376              */
377             biasCts.bias_.doSkippedUpdatesForAllPoints();
378         }
379     }
380
381     wallcycle_stop(wallcycle, ewcAWH);
382
383     return MASTER(commRecord_) ? static_cast<real>(awhPotential) : 0;
384 }
385
386 std::shared_ptr<AwhHistory> Awh::initHistoryFromState() const
387 {
388     if (MASTER(commRecord_))
389     {
390         std::shared_ptr<AwhHistory> awhHistory(new AwhHistory);
391         awhHistory->bias.clear();
392         awhHistory->bias.resize(biasCoupledToSystem_.size());
393
394         for (size_t k = 0; k < awhHistory->bias.size(); k++)
395         {
396             biasCoupledToSystem_[k].bias_.initHistoryFromState(&awhHistory->bias[k]);
397         }
398
399         return awhHistory;
400     }
401     else
402     {
403         /* Return an empty pointer */
404         return std::shared_ptr<AwhHistory>();
405     }
406 }
407
408 void Awh::restoreStateFromHistory(const AwhHistory* awhHistory)
409 {
410     /* Restore the history to the current state */
411     if (MASTER(commRecord_))
412     {
413         GMX_RELEASE_ASSERT(awhHistory != nullptr,
414                            "The master rank should have a valid awhHistory when restoring the "
415                            "state from history.");
416
417         if (awhHistory->bias.size() != biasCoupledToSystem_.size())
418         {
419             GMX_THROW(InvalidInputError(
420                     "AWH state and history contain different numbers of biases. Likely you "
421                     "provided a checkpoint from a different simulation."));
422         }
423
424         potentialOffset_ = awhHistory->potentialOffset;
425     }
426     if (PAR(commRecord_))
427     {
428         gmx_bcast(sizeof(potentialOffset_), &potentialOffset_, commRecord_->mpi_comm_mygroup);
429     }
430
431     for (size_t k = 0; k < biasCoupledToSystem_.size(); k++)
432     {
433         biasCoupledToSystem_[k].bias_.restoreStateFromHistory(
434                 awhHistory ? &awhHistory->bias[k] : nullptr, commRecord_);
435     }
436 }
437
438 void Awh::updateHistory(AwhHistory* awhHistory) const
439 {
440     if (!MASTER(commRecord_))
441     {
442         return;
443     }
444
445     /* This assert will also catch a non-master rank calling this function. */
446     GMX_RELEASE_ASSERT(awhHistory->bias.size() == biasCoupledToSystem_.size(),
447                        "AWH state and history bias count should match");
448
449     awhHistory->potentialOffset = potentialOffset_;
450
451     for (size_t k = 0; k < awhHistory->bias.size(); k++)
452     {
453         biasCoupledToSystem_[k].bias_.updateHistory(&awhHistory->bias[k]);
454     }
455 }
456
457 const char* Awh::externalPotentialString()
458 {
459     return "AWH";
460 }
461
462 void Awh::registerAwhWithPull(const AwhParams& awhParams, pull_t* pull_work)
463 {
464     GMX_RELEASE_ASSERT(!anyDimUsesPull(awhParams) || pull_work, "Need a valid pull object");
465
466     for (int k = 0; k < awhParams.numBias; k++)
467     {
468         const AwhBiasParams& biasParams = awhParams.awhBiasParams[k];
469
470         for (int d = 0; d < biasParams.ndim; d++)
471         {
472             if (biasParams.dimParams[d].eCoordProvider == eawhcoordproviderPULL)
473             {
474                 register_external_pull_potential(pull_work, biasParams.dimParams[d].coordIndex,
475                                                  Awh::externalPotentialString());
476             }
477         }
478     }
479 }
480
481 /* Fill the AWH data block of an energy frame with data (if there is any). */
482 void Awh::writeToEnergyFrame(int64_t step, t_enxframe* frame) const
483 {
484     GMX_ASSERT(MASTER(commRecord_), "writeToEnergyFrame should only be called on the master rank");
485     GMX_ASSERT(frame != nullptr, "Need a valid energy frame");
486
487     if (!isOutputStep(step))
488     {
489         /* This is not an AWH output step, don't write any AWH data */
490         return;
491     }
492
493     /* Get the total number of energy subblocks that AWH needs */
494     int numSubblocks = 0;
495     for (auto& biasCoupledToSystem : biasCoupledToSystem_)
496     {
497         numSubblocks += biasCoupledToSystem.bias_.numEnergySubblocksToWrite();
498     }
499     GMX_ASSERT(numSubblocks > 0, "We should always have data to write");
500
501     /* Add 1 energy block */
502     add_blocks_enxframe(frame, frame->nblock + 1);
503
504     /* Take the block that was just added and set the number of subblocks. */
505     t_enxblock* awhEnergyBlock = &(frame->block[frame->nblock - 1]);
506     add_subblocks_enxblock(awhEnergyBlock, numSubblocks);
507
508     /* Claim it as an AWH block. */
509     awhEnergyBlock->id = enxAWH;
510
511     /* Transfer AWH data blocks to energy sub blocks */
512     int energySubblockCount = 0;
513     for (auto& biasCoupledToSystem : biasCoupledToSystem_)
514     {
515         energySubblockCount += biasCoupledToSystem.bias_.writeToEnergySubblocks(
516                 &(awhEnergyBlock->sub[energySubblockCount]));
517     }
518 }
519
520 bool Awh::hasFepLambdaDimension() const
521 {
522     return std::any_of(
523             std::begin(biasCoupledToSystem_), std::end(biasCoupledToSystem_),
524             [](const auto& coupledBias) { return coupledBias.bias_.hasFepLambdaDimension(); });
525 }
526
527 bool Awh::needForeignEnergyDifferences(const int64_t step) const
528 {
529     /* If there is no FEP lambda dimension at all in any bias there will be no need for foreign
530      * energy differences */
531     if (!hasFepLambdaDimension())
532     {
533         return false;
534     }
535     if (step == 0)
536     {
537         return true;
538     }
539     /* Check whether the bias(es) that has/have a FEP lambda dimension should sample coordinates
540      * this step. Since the biases may have different sampleCoordStep it is necessary to check
541      * this combination. */
542     for (const auto& biasCts : biasCoupledToSystem_)
543     {
544         if (biasCts.bias_.hasFepLambdaDimension() && biasCts.bias_.isSampleCoordStep(step))
545         {
546             return true;
547         }
548     }
549     return false;
550 }
551
552 std::unique_ptr<Awh> prepareAwhModule(FILE*                 fplog,
553                                       const t_inputrec&     inputRecord,
554                                       t_state*              stateGlobal,
555                                       const t_commrec*      commRecord,
556                                       const gmx_multisim_t* multiSimRecord,
557                                       const bool            startingFromCheckpoint,
558                                       const bool            usingShellParticles,
559                                       const std::string&    biasInitFilename,
560                                       pull_t*               pull_work)
561 {
562     if (!inputRecord.bDoAwh)
563     {
564         return nullptr;
565     }
566     if (usingShellParticles)
567     {
568         GMX_THROW(InvalidInputError("AWH biasing does not support shell particles."));
569     }
570
571     auto awh = std::make_unique<Awh>(
572             fplog, inputRecord, commRecord, multiSimRecord, *inputRecord.awhParams, biasInitFilename,
573             pull_work, inputRecord.fepvals->n_lambda, inputRecord.fepvals->init_fep_state);
574
575     if (startingFromCheckpoint)
576     {
577         // Restore the AWH history read from checkpoint
578         awh->restoreStateFromHistory(MASTER(commRecord) ? stateGlobal->awhHistory.get() : nullptr);
579     }
580     else if (MASTER(commRecord))
581     {
582         // Initialize the AWH history here
583         stateGlobal->awhHistory = awh->initHistoryFromState();
584     }
585     return awh;
586 }
587
588 } // namespace gmx