734ba4eda49981a0b06b87c7083bd0bdcca1ec45
[alexxy/gromacs.git] / src / gromacs / modularsimulator / forceelement.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 force 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 "forceelement.h"
45
46 #include "gromacs/domdec/mdsetup.h"
47 #include "gromacs/math/vectypes.h"
48 #include "gromacs/mdlib/constr.h"
49 #include "gromacs/mdlib/force.h"
50 #include "gromacs/mdlib/force_flags.h"
51 #include "gromacs/mdlib/mdatoms.h"
52 #include "gromacs/mdrun/shellfc.h"
53 #include "gromacs/mdtypes/inputrec.h"
54 #include "gromacs/mdtypes/mdatom.h"
55 #include "gromacs/pbcutil/pbc.h"
56
57 #include "energyelement.h"
58 #include "freeenergyperturbationelement.h"
59 #include "statepropagatordata.h"
60
61 struct gmx_edsam;
62 struct gmx_enfrot;
63 struct gmx_multisim_t;
64 class history_t;
65
66 namespace gmx
67 {
68 ForceElement::ForceElement(StatePropagatorData*           statePropagatorData,
69                            EnergyElement*                 energyElement,
70                            FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
71                            bool                           isVerbose,
72                            bool                           isDynamicBox,
73                            FILE*                          fplog,
74                            const t_commrec*               cr,
75                            const t_inputrec*              inputrec,
76                            const MDAtoms*                 mdAtoms,
77                            t_nrnb*                        nrnb,
78                            t_forcerec*                    fr,
79                            t_fcdata*                      fcd,
80                            gmx_wallcycle*                 wcycle,
81                            MdrunScheduleWorkload*         runScheduleWork,
82                            gmx_vsite_t*                   vsite,
83                            ImdSession*                    imdSession,
84                            pull_t*                        pull_work,
85                            Constraints*                   constr,
86                            const gmx_mtop_t*              globalTopology,
87                            gmx_enfrot*                    enforcedRotation) :
88     shellfc_(init_shell_flexcon(fplog,
89                                 globalTopology,
90                                 constr ? constr->numFlexibleConstraints() : 0,
91                                 inputrec->nstcalcenergy,
92                                 DOMAINDECOMP(cr))),
93     doShellFC_(shellfc_ != nullptr),
94     nextNSStep_(-1),
95     nextEnergyCalculationStep_(-1),
96     nextVirialCalculationStep_(-1),
97     nextFreeEnergyCalculationStep_(-1),
98     statePropagatorData_(statePropagatorData),
99     energyElement_(energyElement),
100     freeEnergyPerturbationElement_(freeEnergyPerturbationElement),
101     localTopology_(nullptr),
102     isDynamicBox_(isDynamicBox),
103     isVerbose_(isVerbose),
104     nShellRelaxationSteps_(0),
105     ddBalanceRegionHandler_(cr),
106     lambda_(),
107     fplog_(fplog),
108     cr_(cr),
109     inputrec_(inputrec),
110     mdAtoms_(mdAtoms),
111     nrnb_(nrnb),
112     wcycle_(wcycle),
113     fr_(fr),
114     vsite_(vsite),
115     imdSession_(imdSession),
116     pull_work_(pull_work),
117     fcd_(fcd),
118     runScheduleWork_(runScheduleWork),
119     constr_(constr),
120     enforcedRotation_(enforcedRotation)
121 {
122     lambda_.fill(0);
123
124     if (doShellFC_ && !DOMAINDECOMP(cr))
125     {
126         // This was done in mdAlgorithmsSetupAtomData(), but shellfc
127         // won't be available outside this element.
128         make_local_shells(cr, mdAtoms->mdatoms(), shellfc_);
129     }
130 }
131
132 void ForceElement::scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction)
133 {
134     unsigned int flags =
135             (GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | (isDynamicBox_ ? GMX_FORCE_DYNAMICBOX : 0)
136              | (nextVirialCalculationStep_ == step ? GMX_FORCE_VIRIAL : 0)
137              | (nextEnergyCalculationStep_ == step ? GMX_FORCE_ENERGY : 0)
138              | (nextFreeEnergyCalculationStep_ == step ? GMX_FORCE_DHDL : 0)
139              | (!doShellFC_ && nextNSStep_ == step ? GMX_FORCE_NS : 0));
140
141     (*registerRunFunction)(std::make_unique<SimulatorRunFunction>([this, step, time, flags]() {
142         if (doShellFC_)
143         {
144             run<true>(step, time, flags);
145         }
146         else
147         {
148             run<false>(step, time, flags);
149         }
150     }));
151 }
152
153 void ForceElement::elementSetup()
154 {
155     GMX_ASSERT(localTopology_, "Setup called before local topology was set.");
156 }
157
158 template<bool doShellFC>
159 void ForceElement::run(Step step, Time time, unsigned int flags)
160 {
161     // Disabled functionality
162     gmx_multisim_t* ms = nullptr;
163
164
165     if (!DOMAINDECOMP(cr_) && (flags & GMX_FORCE_NS) && inputrecDynamicBox(inputrec_))
166     {
167         // TODO: Correcting the box is done in DomDecHelper (if using DD) or here (non-DD simulations).
168         //       Think about unifying this responsibility, could this be done in one place?
169         auto box = statePropagatorData_->box();
170         correct_box(fplog_, step, box);
171     }
172
173     /* The coordinates (x) are shifted (to get whole molecules)
174      * in do_force.
175      * This is parallelized as well, and does communication too.
176      * Check comments in sim_util.c
177      */
178     auto       x      = statePropagatorData_->positionsView();
179     auto       forces = statePropagatorData_->forcesView();
180     auto       box    = statePropagatorData_->constBox();
181     history_t* hist   = nullptr; // disabled
182
183     tensor force_vir = { { 0 } };
184     // TODO: Make lambda const (needs some adjustments in lower force routines)
185     ArrayRef<real> lambda =
186             freeEnergyPerturbationElement_ ? freeEnergyPerturbationElement_->lambdaView() : lambda_;
187
188     if (doShellFC)
189     {
190         auto v = statePropagatorData_->velocitiesView();
191
192         relax_shell_flexcon(
193                 fplog_, cr_, ms, isVerbose_, enforcedRotation_, step, inputrec_, imdSession_,
194                 pull_work_, step == nextNSStep_, static_cast<int>(flags), localTopology_, constr_,
195                 energyElement_->enerdata(), fcd_, statePropagatorData_->localNumAtoms(), x, v, box,
196                 lambda, hist, forces, force_vir, mdAtoms_->mdatoms(), nrnb_, wcycle_, shellfc_, fr_,
197                 runScheduleWork_, time, energyElement_->muTot(), vsite_, ddBalanceRegionHandler_);
198         nShellRelaxationSteps_++;
199     }
200     else
201     {
202         // Disabled functionality
203         Awh*       awh = nullptr;
204         gmx_edsam* ed  = nullptr;
205
206         do_force(fplog_, cr_, ms, inputrec_, awh, enforcedRotation_, imdSession_, pull_work_, step,
207                  nrnb_, wcycle_, localTopology_, box, x, hist, forces, force_vir, mdAtoms_->mdatoms(),
208                  energyElement_->enerdata(), fcd_, lambda, fr_, runScheduleWork_, vsite_,
209                  energyElement_->muTot(), time, ed, static_cast<int>(flags), ddBalanceRegionHandler_);
210     }
211     energyElement_->addToForceVirial(force_vir, step);
212 }
213
214 void ForceElement::elementTeardown()
215 {
216     if (doShellFC_)
217     {
218         done_shellfc(fplog_, shellfc_, nShellRelaxationSteps_);
219     }
220 }
221
222 void ForceElement::setTopology(const gmx_localtop_t* top)
223 {
224     localTopology_ = top;
225 }
226
227 SignallerCallbackPtr ForceElement::registerNSCallback()
228 {
229     return std::make_unique<SignallerCallback>(
230             [this](Step step, Time gmx_unused time) { this->nextNSStep_ = step; });
231 }
232
233 SignallerCallbackPtr ForceElement::registerEnergyCallback(EnergySignallerEvent event)
234 {
235     if (event == EnergySignallerEvent::EnergyCalculationStep)
236     {
237         return std::make_unique<SignallerCallback>(
238                 [this](Step step, Time /*unused*/) { nextEnergyCalculationStep_ = step; });
239     }
240     if (event == EnergySignallerEvent::VirialCalculationStep)
241     {
242         return std::make_unique<SignallerCallback>(
243                 [this](Step step, Time /*unused*/) { nextVirialCalculationStep_ = step; });
244     }
245     if (event == EnergySignallerEvent::FreeEnergyCalculationStep)
246     {
247         return std::make_unique<SignallerCallback>(
248                 [this](Step step, Time /*unused*/) { nextFreeEnergyCalculationStep_ = step; });
249     }
250     return nullptr;
251 }
252 } // namespace gmx