Remove inputrec from initialize_lambdas signature
[alexxy/gromacs.git] / src / gromacs / modularsimulator / freeenergyperturbationdata.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 /*! \internal \file
36  * \brief Defines the free energy perturbation element for the modular simulator
37  *
38  * \author Pascal Merz <pascal.merz@me.com>
39  * \ingroup module_modularsimulator
40  */
41
42 #include "gmxpre.h"
43
44 #include "freeenergyperturbationdata.h"
45
46 #include "gromacs/domdec/domdec_network.h"
47 #include "gromacs/mdlib/freeenergyparameters.h"
48 #include "gromacs/mdlib/md_support.h"
49 #include "gromacs/mdlib/mdatoms.h"
50 #include "gromacs/mdtypes/checkpointdata.h"
51 #include "gromacs/mdtypes/commrec.h"
52 #include "gromacs/mdtypes/inputrec.h"
53 #include "gromacs/mdtypes/mdatom.h"
54 #include "gromacs/mdtypes/state.h"
55 #include "gromacs/trajectory/trajectoryframe.h"
56
57 #include "modularsimulator.h"
58 #include "simulatoralgorithm.h"
59
60 namespace gmx
61 {
62
63 FreeEnergyPerturbationData::FreeEnergyPerturbationData(FILE* fplog, const t_inputrec& inputrec, MDAtoms* mdAtoms) :
64     element_(std::make_unique<Element>(this, inputrec.fepvals->delta_lambda)),
65     lambda_(),
66     currentFEPState_(0),
67     fplog_(fplog),
68     inputrec_(inputrec),
69     mdAtoms_(mdAtoms)
70 {
71     std::fill(lambda_.begin(), lambda_.end(), 0);
72     // The legacy implementation only filled the lambda vector in state_global, which is only
73     // available on master. We have the lambda vector available everywhere, so we pass a `true`
74     // for isMaster on all ranks. See #3647.
75     initialize_lambdas(fplog_,
76                        inputrec_.efep,
77                        inputrec_.bSimTemp,
78                        *inputrec_.fepvals,
79                        inputrec_.simtempvals->temperatures,
80                        gmx::arrayRefFromArray(inputrec_.opts.ref_t, inputrec_.opts.ngtc),
81                        true,
82                        &currentFEPState_,
83                        lambda_);
84 }
85
86 void FreeEnergyPerturbationData::Element::scheduleTask(Step step,
87                                                        Time gmx_unused            time,
88                                                        const RegisterRunFunction& registerRunFunction)
89 {
90     if (lambdasChange_)
91     {
92         registerRunFunction([this, step]() { freeEnergyPerturbationData_->updateLambdas(step); });
93     }
94 }
95
96 void FreeEnergyPerturbationData::updateLambdas(Step step)
97 {
98     // at beginning of step (if lambdas change...)
99     lambda_ = currentLambdas(step, *(inputrec_.fepvals), currentFEPState_);
100     updateMDAtoms();
101 }
102
103 ArrayRef<real> FreeEnergyPerturbationData::lambdaView()
104 {
105     return lambda_;
106 }
107
108 ArrayRef<const real> FreeEnergyPerturbationData::constLambdaView()
109 {
110     return lambda_;
111 }
112
113 int FreeEnergyPerturbationData::currentFEPState()
114 {
115     return currentFEPState_;
116 }
117
118 void FreeEnergyPerturbationData::updateMDAtoms()
119 {
120     update_mdatoms(mdAtoms_->mdatoms(), lambda_[FreeEnergyPerturbationCouplingType::Mass]);
121 }
122
123 namespace
124 {
125 /*!
126  * \brief Enum describing the contents FreeEnergyPerturbationData::Element writes to modular checkpoint
127  *
128  * When changing the checkpoint content, add a new element just above Count, and adjust the
129  * checkpoint functionality.
130  */
131 enum class CheckpointVersion
132 {
133     Base, //!< First version of modular checkpointing
134     Count //!< Number of entries. Add new versions right above this!
135 };
136 constexpr auto c_currentVersion = CheckpointVersion(int(CheckpointVersion::Count) - 1);
137 } // namespace
138
139 template<CheckpointDataOperation operation>
140 void FreeEnergyPerturbationData::doCheckpointData(CheckpointData<operation>* checkpointData)
141 {
142     checkpointVersion(checkpointData, "FreeEnergyPerturbationData version", c_currentVersion);
143
144     checkpointData->scalar("current FEP state", &currentFEPState_);
145     checkpointData->arrayRef("lambda vector", makeCheckpointArrayRef<operation>(lambda_));
146 }
147
148 void FreeEnergyPerturbationData::Element::saveCheckpointState(std::optional<WriteCheckpointData> checkpointData,
149                                                               const t_commrec*                   cr)
150 {
151     if (MASTER(cr))
152     {
153         freeEnergyPerturbationData_->doCheckpointData<CheckpointDataOperation::Write>(
154                 &checkpointData.value());
155     }
156 }
157
158 void FreeEnergyPerturbationData::Element::restoreCheckpointState(std::optional<ReadCheckpointData> checkpointData,
159                                                                  const t_commrec* cr)
160 {
161     if (MASTER(cr))
162     {
163         freeEnergyPerturbationData_->doCheckpointData<CheckpointDataOperation::Read>(
164                 &checkpointData.value());
165     }
166     if (DOMAINDECOMP(cr))
167     {
168         dd_bcast(cr->dd, sizeof(int), &freeEnergyPerturbationData_->currentFEPState_);
169         dd_bcast(cr->dd,
170                  ssize(freeEnergyPerturbationData_->lambda_) * int(sizeof(real)),
171                  freeEnergyPerturbationData_->lambda_.data());
172     }
173 }
174
175 const std::string& FreeEnergyPerturbationData::Element::clientID()
176 {
177     return FreeEnergyPerturbationData::checkpointID();
178 }
179
180 FreeEnergyPerturbationData::Element::Element(FreeEnergyPerturbationData* freeEnergyPerturbationElement,
181                                              double                      deltaLambda) :
182     freeEnergyPerturbationData_(freeEnergyPerturbationElement),
183     lambdasChange_(deltaLambda != 0)
184 {
185 }
186
187 void FreeEnergyPerturbationData::Element::elementSetup()
188 {
189     freeEnergyPerturbationData_->updateMDAtoms();
190 }
191
192 FreeEnergyPerturbationData::Element* FreeEnergyPerturbationData::element()
193 {
194     return element_.get();
195 }
196
197 ISimulatorElement* FreeEnergyPerturbationData::Element::getElementPointerImpl(
198         LegacySimulatorData gmx_unused*        legacySimulatorData,
199         ModularSimulatorAlgorithmBuilderHelper gmx_unused* builderHelper,
200         StatePropagatorData gmx_unused* statePropagatorData,
201         EnergyData gmx_unused*      energyData,
202         FreeEnergyPerturbationData* freeEnergyPerturbationData,
203         GlobalCommunicationHelper gmx_unused* globalCommunicationHelper)
204 {
205     return freeEnergyPerturbationData->element();
206 }
207
208 void FreeEnergyPerturbationData::readCheckpointToTrxFrame(t_trxframe* trxFrame,
209                                                           std::optional<ReadCheckpointData> readCheckpointData)
210 {
211     if (readCheckpointData)
212     {
213         FreeEnergyPerturbationData freeEnergyPerturbationData(nullptr, t_inputrec(), nullptr);
214         freeEnergyPerturbationData.doCheckpointData(&readCheckpointData.value());
215         trxFrame->lambda = freeEnergyPerturbationData.lambda_[FreeEnergyPerturbationCouplingType::Fep];
216         trxFrame->fep_state = freeEnergyPerturbationData.currentFEPState_;
217     }
218     else
219     {
220         trxFrame->lambda    = 0;
221         trxFrame->fep_state = 0;
222     }
223     trxFrame->bLambda = true;
224 }
225
226 const std::string& FreeEnergyPerturbationData::checkpointID()
227 {
228     static const std::string identifier = "FreeEnergyPerturbationData";
229     return identifier;
230 }
231
232 } // namespace gmx