c60a9f5f684294f8c97ac7f5abd25b9180e6d60c
[alexxy/gromacs.git] / src / gromacs / modularsimulator / computeglobalselement.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 global reduction 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 "computeglobalselement.h"
45
46 #include "gromacs/domdec/partition.h"
47 #include "gromacs/gmxlib/nrnb.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdlib/md_support.h"
50 #include "gromacs/mdlib/mdatoms.h"
51 #include "gromacs/mdlib/stat.h"
52 #include "gromacs/mdtypes/group.h"
53 #include "gromacs/mdtypes/inputrec.h"
54 #include "gromacs/mdtypes/md_enums.h"
55 #include "gromacs/topology/topology.h"
56
57 #include "freeenergyperturbationelement.h"
58
59 namespace gmx
60 {
61 template<ComputeGlobalsAlgorithm algorithm>
62 ComputeGlobalsElement<algorithm>::ComputeGlobalsElement(StatePropagatorData* statePropagatorData,
63                                                         EnergyElement*       energyElement,
64                                                         FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
65                                                         SimulationSignals* signals,
66                                                         int                nstglobalcomm,
67                                                         FILE*              fplog,
68                                                         const MDLogger&    mdlog,
69                                                         t_commrec*         cr,
70                                                         const t_inputrec*  inputrec,
71                                                         const MDAtoms*     mdAtoms,
72                                                         t_nrnb*            nrnb,
73                                                         gmx_wallcycle*     wcycle,
74                                                         t_forcerec*        fr,
75                                                         const gmx_mtop_t*  global_top,
76                                                         Constraints*       constr,
77                                                         bool               hasReadEkinState) :
78     energyReductionStep_(-1),
79     virialReductionStep_(-1),
80     vvSchedulingStep_(-1),
81     doStopCM_(inputrec->comm_mode != ecmNO),
82     nstcomm_(inputrec->nstcomm),
83     nstglobalcomm_(nstglobalcomm),
84     lastStep_(inputrec->nsteps + inputrec->init_step),
85     initStep_(inputrec->init_step),
86     nullSignaller_(std::make_unique<SimulationSignaller>(nullptr, nullptr, nullptr, false, false)),
87     hasReadEkinState_(hasReadEkinState),
88     totalNumberOfBondedInteractions_(0),
89     shouldCheckNumberOfBondedInteractions_(false),
90     statePropagatorData_(statePropagatorData),
91     energyElement_(energyElement),
92     localTopology_(nullptr),
93     freeEnergyPerturbationElement_(freeEnergyPerturbationElement),
94     vcm_(global_top->groups, *inputrec),
95     signals_(signals),
96     fplog_(fplog),
97     mdlog_(mdlog),
98     cr_(cr),
99     inputrec_(inputrec),
100     top_global_(global_top),
101     mdAtoms_(mdAtoms),
102     constr_(constr),
103     nrnb_(nrnb),
104     wcycle_(wcycle),
105     fr_(fr)
106 {
107     reportComRemovalInfo(fplog, vcm_);
108     gstat_ = global_stat_init(inputrec_);
109 }
110
111 template<ComputeGlobalsAlgorithm algorithm>
112 ComputeGlobalsElement<algorithm>::~ComputeGlobalsElement()
113 {
114     global_stat_destroy(gstat_);
115 }
116
117 template<ComputeGlobalsAlgorithm algorithm>
118 void ComputeGlobalsElement<algorithm>::elementSetup()
119 {
120     GMX_ASSERT(localTopology_, "Setup called before local topology was set.");
121
122     if (doStopCM_ && !inputrec_->bContinuation)
123     {
124         // To minimize communication, compute_globals computes the COM velocity
125         // and the kinetic energy for the velocities without COM motion removed.
126         // Thus to get the kinetic energy without the COM contribution, we need
127         // to call compute_globals twice.
128
129         compute(-1, CGLO_GSTAT | CGLO_STOPCM, nullSignaller_.get(), false, true);
130
131         auto v = as_rvec_array(statePropagatorData_->velocitiesView().paddedArrayRef().data());
132         // At initialization, do not pass x with acceleration-correction mode
133         // to avoid (incorrect) correction of the initial coordinates.
134         rvec* xPtr = nullptr;
135         if (vcm_.mode != ecmLINEAR_ACCELERATION_CORRECTION)
136         {
137             xPtr = as_rvec_array(statePropagatorData_->positionsView().paddedArrayRef().data());
138         }
139         process_and_stopcm_grp(fplog_, &vcm_, *mdAtoms_->mdatoms(), xPtr, v);
140         inc_nrnb(nrnb_, eNR_STOPCM, mdAtoms_->mdatoms()->homenr);
141     }
142
143     unsigned int cglo_flags =
144             (CGLO_TEMPERATURE | CGLO_GSTAT
145              | (shouldCheckNumberOfBondedInteractions_ ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
146              | (hasReadEkinState_ ? CGLO_READEKIN : 0));
147
148     if (algorithm == ComputeGlobalsAlgorithm::VelocityVerlet)
149     {
150         cglo_flags |= CGLO_PRESSURE | CGLO_CONSTRAINT;
151     }
152
153     compute(-1, cglo_flags, nullSignaller_.get(), false, true);
154
155     // Calculate the initial half step temperature, and save the ekinh_old
156     for (int i = 0; (i < inputrec_->opts.ngtc); i++)
157     {
158         copy_mat(energyElement_->ekindata()->tcstat[i].ekinh,
159                  energyElement_->ekindata()->tcstat[i].ekinh_old);
160     }
161 }
162
163 template<ComputeGlobalsAlgorithm algorithm>
164 void ComputeGlobalsElement<algorithm>::scheduleTask(Step step,
165                                                     Time gmx_unused               time,
166                                                     const RegisterRunFunctionPtr& registerRunFunction)
167 {
168     const bool needComReduction    = doStopCM_ && do_per_step(step, nstcomm_);
169     const bool needGlobalReduction = step == energyReductionStep_ || step == virialReductionStep_
170                                      || needComReduction || do_per_step(step, nstglobalcomm_);
171
172     // TODO: CGLO_GSTAT is only used for needToSumEkinhOld_, i.e. to signal that we do or do not
173     //       sum the previous kinetic energy. We should simplify / clarify this.
174
175     if (algorithm == ComputeGlobalsAlgorithm::LeapFrog)
176     {
177         // With Leap-Frog we can skip compute_globals at
178         // non-communication steps, but we need to calculate
179         // the kinetic energy one step before communication.
180
181         // With leap-frog we also need to compute the half-step
182         // kinetic energy at the step before we need to write
183         // the full-step kinetic energy
184         const bool needEkinAtNextStep = (do_per_step(step + 1, nstglobalcomm_) || step + 1 == lastStep_);
185
186         if (!needGlobalReduction && !needEkinAtNextStep)
187         {
188             return;
189         }
190
191         const bool doEnergy = step == energyReductionStep_;
192         int        flags =
193                 (needGlobalReduction ? CGLO_GSTAT : 0) | (doEnergy ? CGLO_ENERGY : 0)
194                 | (needComReduction ? CGLO_STOPCM : 0) | CGLO_TEMPERATURE | CGLO_PRESSURE | CGLO_CONSTRAINT
195                 | (shouldCheckNumberOfBondedInteractions_ ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0);
196
197         // Since we're already communicating at this step, we
198         // can propagate intra-simulation signals. Note that
199         // check_nstglobalcomm has the responsibility for
200         // choosing the value of nstglobalcomm which satisfies
201         // the need of the different signallers.
202         const bool doIntraSimSignal = true;
203         // Disable functionality
204         const bool doInterSimSignal = false;
205
206         // Make signaller to signal stop / reset / checkpointing signals
207         auto signaller = std::make_shared<SimulationSignaller>(signals_, cr_, nullptr,
208                                                                doInterSimSignal, doIntraSimSignal);
209
210         (*registerRunFunction)(std::make_unique<SimulatorRunFunction>(
211                 [this, step, flags, signaller = std::move(signaller)]() {
212                     compute(step, flags, signaller.get(), true);
213                 }));
214     }
215     else if (algorithm == ComputeGlobalsAlgorithm::VelocityVerlet)
216     {
217         // For VV, we schedule two calls to compute globals per step.
218         if (step != vvSchedulingStep_)
219         {
220             // This is the first scheduling call for this step (positions & velocities at full time
221             // step) Set this as the current scheduling step
222             vvSchedulingStep_ = step;
223
224             // For vv, the state at the beginning of the step is positions at time t, velocities at time t - dt/2
225             // The first velocity propagation (+dt/2) therefore actually corresponds to the previous step.
226             // So we need information from the last step in the first half of the integration
227             if (!needGlobalReduction && !do_per_step(step - 1, nstglobalcomm_))
228             {
229                 return;
230             }
231
232             const bool doTemperature = step != initStep_ || inputrec_->bContinuation;
233             const bool doEnergy      = step == energyReductionStep_;
234
235             int flags = (needGlobalReduction ? CGLO_GSTAT : 0) | (doEnergy ? CGLO_ENERGY : 0)
236                         | (doTemperature ? CGLO_TEMPERATURE : 0) | CGLO_PRESSURE | CGLO_CONSTRAINT
237                         | (needComReduction ? CGLO_STOPCM : 0)
238                         | (shouldCheckNumberOfBondedInteractions_ ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
239                                                                   : 0)
240                         | CGLO_SCALEEKIN;
241
242             (*registerRunFunction)(std::make_unique<SimulatorRunFunction>(
243                     [this, step, flags]() { compute(step, flags, nullSignaller_.get(), false); }));
244         }
245         else
246         {
247             // second call to compute_globals for this step
248             // Reset the scheduling step to avoid confusion if scheduling needs
249             // to be repeated (in case of unexpected simulation termination)
250             vvSchedulingStep_ = -1;
251
252             if (!needGlobalReduction)
253             {
254                 return;
255             }
256             int flags = CGLO_GSTAT | CGLO_CONSTRAINT
257                         | (shouldCheckNumberOfBondedInteractions_ ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
258                                                                   : 0);
259
260             // Since we're already communicating at this step, we
261             // can propagate intra-simulation signals. Note that
262             // check_nstglobalcomm has the responsibility for
263             // choosing the value of nstglobalcomm which satisfies
264             // the need of the different signallers.
265             const bool doIntraSimSignal = true;
266             // Disable functionality
267             const bool doInterSimSignal = false;
268
269             auto signaller = std::make_shared<SimulationSignaller>(
270                     signals_, cr_, nullptr, doInterSimSignal, doIntraSimSignal);
271
272             (*registerRunFunction)(std::make_unique<SimulatorRunFunction>(
273                     [this, step, flags, signaller = std::move(signaller)]() {
274                         compute(step, flags, signaller.get(), true);
275                     }));
276         }
277     }
278 }
279
280 template<ComputeGlobalsAlgorithm algorithm>
281 void ComputeGlobalsElement<algorithm>::compute(gmx::Step            step,
282                                                unsigned int         flags,
283                                                SimulationSignaller* signaller,
284                                                bool                 useLastBox,
285                                                bool                 isInit)
286 {
287     auto x       = as_rvec_array(statePropagatorData_->positionsView().paddedArrayRef().data());
288     auto v       = as_rvec_array(statePropagatorData_->velocitiesView().paddedArrayRef().data());
289     auto box     = statePropagatorData_->constBox();
290     auto lastbox = useLastBox ? statePropagatorData_->constPreviousBox()
291                               : statePropagatorData_->constBox();
292
293     const real vdwLambda = freeEnergyPerturbationElement_
294                                    ? freeEnergyPerturbationElement_->constLambdaView()[efptVDW]
295                                    : 0;
296
297     compute_globals(gstat_, cr_, inputrec_, fr_, energyElement_->ekindata(), x, v, box, vdwLambda,
298                     mdAtoms_->mdatoms(), nrnb_, &vcm_, step != -1 ? wcycle_ : nullptr,
299                     energyElement_->enerdata(), energyElement_->forceVirial(step),
300                     energyElement_->constraintVirial(step), energyElement_->totalVirial(step),
301                     energyElement_->pressure(step), constr_, signaller, lastbox,
302                     &totalNumberOfBondedInteractions_, energyElement_->needToSumEkinhOld(), flags);
303     checkNumberOfBondedInteractions(mdlog_, cr_, totalNumberOfBondedInteractions_, top_global_,
304                                     localTopology_, x, box, &shouldCheckNumberOfBondedInteractions_);
305     if (flags & CGLO_STOPCM && !isInit)
306     {
307         process_and_stopcm_grp(fplog_, &vcm_, *mdAtoms_->mdatoms(), x, v);
308         inc_nrnb(nrnb_, eNR_STOPCM, mdAtoms_->mdatoms()->homenr);
309     }
310 }
311
312 template<ComputeGlobalsAlgorithm algorithm>
313 CheckBondedInteractionsCallbackPtr ComputeGlobalsElement<algorithm>::getCheckNumberOfBondedInteractionsCallback()
314 {
315     return std::make_unique<CheckBondedInteractionsCallback>(
316             [this]() { needToCheckNumberOfBondedInteractions(); });
317 }
318
319 template<ComputeGlobalsAlgorithm algorithm>
320 void ComputeGlobalsElement<algorithm>::needToCheckNumberOfBondedInteractions()
321 {
322     shouldCheckNumberOfBondedInteractions_ = true;
323 }
324
325 template<ComputeGlobalsAlgorithm algorithm>
326 void ComputeGlobalsElement<algorithm>::setTopology(const gmx_localtop_t* top)
327 {
328     localTopology_ = top;
329 }
330
331 template<ComputeGlobalsAlgorithm algorithm>
332 SignallerCallbackPtr ComputeGlobalsElement<algorithm>::registerEnergyCallback(EnergySignallerEvent event)
333 {
334     if (event == EnergySignallerEvent::EnergyCalculationStep)
335     {
336         return std::make_unique<SignallerCallback>(
337                 [this](Step step, Time /*unused*/) { energyReductionStep_ = step; });
338     }
339     if (event == EnergySignallerEvent::VirialCalculationStep)
340     {
341         return std::make_unique<SignallerCallback>(
342                 [this](Step step, Time /*unused*/) { virialReductionStep_ = step; });
343     }
344     return nullptr;
345 }
346
347 template<ComputeGlobalsAlgorithm algorithm>
348 SignallerCallbackPtr ComputeGlobalsElement<algorithm>::registerTrajectorySignallerCallback(TrajectoryEvent event)
349 {
350     if (event == TrajectoryEvent::EnergyWritingStep)
351     {
352         return std::make_unique<SignallerCallback>(
353                 [this](Step step, Time /*unused*/) { energyReductionStep_ = step; });
354     }
355     return nullptr;
356 }
357
358 //! Explicit template instantiation
359 //! @{
360 template class ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>;
361 template class ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>;
362 //! @}
363 } // namespace gmx