e8f1d56bba05e0301d50f24551ad1e59070afe4a
[alexxy/gromacs.git] / src / gromacs / modularsimulator / freeenergyperturbationelement.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019, 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 "freeenergyperturbationelement.h"
45
46 #include "gromacs/mdlib/md_support.h"
47 #include "gromacs/mdlib/mdatoms.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/mdtypes/state.h"
50
51 namespace gmx
52 {
53
54 FreeEnergyPerturbationElement::FreeEnergyPerturbationElement(FILE*             fplog,
55                                                              const t_inputrec* inputrec,
56                                                              MDAtoms*          mdAtoms) :
57     lambda_(),
58     lambda0_(),
59     currentFEPState_(0),
60     lambdasChange_(inputrec->fepvals->delta_lambda != 0),
61     fplog_(fplog),
62     inputrec_(inputrec),
63     mdAtoms_(mdAtoms)
64 {
65     lambda_.fill(0);
66     lambda0_.fill(0);
67     initialize_lambdas(fplog_, *inputrec_, true, &currentFEPState_, lambda_, lambda0_.data());
68     update_mdatoms(mdAtoms_->mdatoms(), lambda_[efptMASS]);
69 }
70
71 void FreeEnergyPerturbationElement::scheduleTask(Step step,
72                                                  Time gmx_unused               time,
73                                                  const RegisterRunFunctionPtr& registerRunFunction)
74 {
75     if (lambdasChange_)
76     {
77         (*registerRunFunction)(
78                 std::make_unique<SimulatorRunFunction>([this, step]() { updateLambdas(step); }));
79     }
80 }
81
82 void FreeEnergyPerturbationElement::updateLambdas(Step step)
83 {
84     // at beginning of step (if lambdas change...)
85     setCurrentLambdasLocal(step, inputrec_->fepvals, lambda0_.data(), lambda_, currentFEPState_);
86     update_mdatoms(mdAtoms_->mdatoms(), lambda_[efptMASS]);
87 }
88
89 ArrayRef<real> FreeEnergyPerturbationElement::lambdaView()
90 {
91     return lambda_;
92 }
93
94 ArrayRef<const real> FreeEnergyPerturbationElement::constLambdaView()
95 {
96     return lambda_;
97 }
98
99 int FreeEnergyPerturbationElement::currentFEPState()
100 {
101     return currentFEPState_;
102 }
103
104 void FreeEnergyPerturbationElement::writeCheckpoint(t_state* localState, t_state gmx_unused* globalState)
105 {
106     localState->fep_state = currentFEPState_;
107     localState->lambda    = lambda_;
108     localState->flags |= (1U << estLAMBDA) | (1U << estFEPSTATE);
109 }
110
111 } // namespace gmx