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