Move dispersion correction call to do_force()
[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 = (CGLO_TEMPERATURE | CGLO_GSTAT | (hasReadEkinState_ ? CGLO_READEKIN : 0));
143
144     if (algorithm == ComputeGlobalsAlgorithm::VelocityVerlet)
145     {
146         cglo_flags |= CGLO_PRESSURE | CGLO_CONSTRAINT;
147     }
148
149     compute(-1, cglo_flags, nullSignaller_.get(), false, true);
150
151     // Calculate the initial half step temperature, and save the ekinh_old
152     for (int i = 0; (i < inputrec_->opts.ngtc); i++)
153     {
154         copy_mat(energyElement_->ekindata()->tcstat[i].ekinh,
155                  energyElement_->ekindata()->tcstat[i].ekinh_old);
156     }
157 }
158
159 template<ComputeGlobalsAlgorithm algorithm>
160 void ComputeGlobalsElement<algorithm>::scheduleTask(Step step,
161                                                     Time gmx_unused               time,
162                                                     const RegisterRunFunctionPtr& registerRunFunction)
163 {
164     const bool needComReduction    = doStopCM_ && do_per_step(step, nstcomm_);
165     const bool needGlobalReduction = step == energyReductionStep_ || step == virialReductionStep_
166                                      || needComReduction || do_per_step(step, nstglobalcomm_);
167
168     // TODO: CGLO_GSTAT is only used for needToSumEkinhOld_, i.e. to signal that we do or do not
169     //       sum the previous kinetic energy. We should simplify / clarify this.
170
171     if (algorithm == ComputeGlobalsAlgorithm::LeapFrog)
172     {
173         // With Leap-Frog we can skip compute_globals at
174         // non-communication steps, but we need to calculate
175         // the kinetic energy one step before communication.
176
177         // With leap-frog we also need to compute the half-step
178         // kinetic energy at the step before we need to write
179         // the full-step kinetic energy
180         const bool needEkinAtNextStep = (do_per_step(step + 1, nstglobalcomm_) || step + 1 == lastStep_);
181
182         if (!needGlobalReduction && !needEkinAtNextStep)
183         {
184             return;
185         }
186
187         const bool doEnergy = step == energyReductionStep_;
188         int        flags    = (needGlobalReduction ? CGLO_GSTAT : 0) | (doEnergy ? CGLO_ENERGY : 0)
189                     | (needComReduction ? CGLO_STOPCM : 0) | CGLO_TEMPERATURE | CGLO_PRESSURE
190                     | CGLO_CONSTRAINT;
191
192         // Since we're already communicating at this step, we
193         // can propagate intra-simulation signals. Note that
194         // check_nstglobalcomm has the responsibility for
195         // choosing the value of nstglobalcomm which satisfies
196         // the need of the different signallers.
197         const bool doIntraSimSignal = true;
198         // Disable functionality
199         const bool doInterSimSignal = false;
200
201         // Make signaller to signal stop / reset / checkpointing signals
202         auto signaller = std::make_shared<SimulationSignaller>(signals_, cr_, nullptr,
203                                                                doInterSimSignal, doIntraSimSignal);
204
205         (*registerRunFunction)(std::make_unique<SimulatorRunFunction>(
206                 [this, step, flags, signaller = std::move(signaller)]() {
207                     compute(step, flags, signaller.get(), true);
208                 }));
209     }
210     else if (algorithm == ComputeGlobalsAlgorithm::VelocityVerlet)
211     {
212         // For VV, we schedule two calls to compute globals per step.
213         if (step != vvSchedulingStep_)
214         {
215             // This is the first scheduling call for this step (positions & velocities at full time
216             // step) Set this as the current scheduling step
217             vvSchedulingStep_ = step;
218
219             // For vv, the state at the beginning of the step is positions at time t, velocities at time t - dt/2
220             // The first velocity propagation (+dt/2) therefore actually corresponds to the previous step.
221             // So we need information from the last step in the first half of the integration
222             if (!needGlobalReduction && !do_per_step(step - 1, nstglobalcomm_))
223             {
224                 return;
225             }
226
227             const bool doTemperature = step != initStep_ || inputrec_->bContinuation;
228             const bool doEnergy      = step == energyReductionStep_;
229
230             int flags = (needGlobalReduction ? CGLO_GSTAT : 0) | (doEnergy ? CGLO_ENERGY : 0)
231                         | (doTemperature ? CGLO_TEMPERATURE : 0) | CGLO_PRESSURE | CGLO_CONSTRAINT
232                         | (needComReduction ? CGLO_STOPCM : 0) | CGLO_SCALEEKIN;
233
234             (*registerRunFunction)(std::make_unique<SimulatorRunFunction>(
235                     [this, step, flags]() { compute(step, flags, nullSignaller_.get(), false); }));
236         }
237         else
238         {
239             // second call to compute_globals for this step
240             // Reset the scheduling step to avoid confusion if scheduling needs
241             // to be repeated (in case of unexpected simulation termination)
242             vvSchedulingStep_ = -1;
243
244             if (!needGlobalReduction)
245             {
246                 return;
247             }
248             int flags = CGLO_GSTAT | CGLO_CONSTRAINT;
249
250             // Since we're already communicating at this step, we
251             // can propagate intra-simulation signals. Note that
252             // check_nstglobalcomm has the responsibility for
253             // choosing the value of nstglobalcomm which satisfies
254             // the need of the different signallers.
255             const bool doIntraSimSignal = true;
256             // Disable functionality
257             const bool doInterSimSignal = false;
258
259             auto signaller = std::make_shared<SimulationSignaller>(
260                     signals_, cr_, nullptr, doInterSimSignal, doIntraSimSignal);
261
262             (*registerRunFunction)(std::make_unique<SimulatorRunFunction>(
263                     [this, step, flags, signaller = std::move(signaller)]() {
264                         compute(step, flags, signaller.get(), true);
265                     }));
266         }
267     }
268 }
269
270 template<ComputeGlobalsAlgorithm algorithm>
271 void ComputeGlobalsElement<algorithm>::compute(gmx::Step            step,
272                                                unsigned int         flags,
273                                                SimulationSignaller* signaller,
274                                                bool                 useLastBox,
275                                                bool                 isInit)
276 {
277     auto x       = statePropagatorData_->positionsView().unpaddedArrayRef();
278     auto v       = statePropagatorData_->velocitiesView().unpaddedArrayRef();
279     auto box     = statePropagatorData_->constBox();
280     auto lastbox = useLastBox ? statePropagatorData_->constPreviousBox()
281                               : statePropagatorData_->constBox();
282
283     compute_globals(
284             gstat_, cr_, inputrec_, fr_, energyElement_->ekindata(), x, v, box, mdAtoms_->mdatoms(),
285             nrnb_, &vcm_, step != -1 ? wcycle_ : nullptr, energyElement_->enerdata(),
286             energyElement_->forceVirial(step), energyElement_->constraintVirial(step),
287             energyElement_->totalVirial(step), energyElement_->pressure(step), constr_, signaller,
288             lastbox, &totalNumberOfBondedInteractions_, energyElement_->needToSumEkinhOld(),
289             flags | (shouldCheckNumberOfBondedInteractions_ ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
290     checkNumberOfBondedInteractions(mdlog_, cr_, totalNumberOfBondedInteractions_, top_global_,
291                                     localTopology_, x, box, &shouldCheckNumberOfBondedInteractions_);
292     if (flags & CGLO_STOPCM && !isInit)
293     {
294         process_and_stopcm_grp(fplog_, &vcm_, *mdAtoms_->mdatoms(), x, v);
295         inc_nrnb(nrnb_, eNR_STOPCM, mdAtoms_->mdatoms()->homenr);
296     }
297 }
298
299 template<ComputeGlobalsAlgorithm algorithm>
300 CheckBondedInteractionsCallbackPtr ComputeGlobalsElement<algorithm>::getCheckNumberOfBondedInteractionsCallback()
301 {
302     return std::make_unique<CheckBondedInteractionsCallback>(
303             [this]() { needToCheckNumberOfBondedInteractions(); });
304 }
305
306 template<ComputeGlobalsAlgorithm algorithm>
307 void ComputeGlobalsElement<algorithm>::needToCheckNumberOfBondedInteractions()
308 {
309     shouldCheckNumberOfBondedInteractions_ = true;
310 }
311
312 template<ComputeGlobalsAlgorithm algorithm>
313 void ComputeGlobalsElement<algorithm>::setTopology(const gmx_localtop_t* top)
314 {
315     localTopology_ = top;
316 }
317
318 template<ComputeGlobalsAlgorithm algorithm>
319 SignallerCallbackPtr ComputeGlobalsElement<algorithm>::registerEnergyCallback(EnergySignallerEvent event)
320 {
321     if (event == EnergySignallerEvent::EnergyCalculationStep)
322     {
323         return std::make_unique<SignallerCallback>(
324                 [this](Step step, Time /*unused*/) { energyReductionStep_ = step; });
325     }
326     if (event == EnergySignallerEvent::VirialCalculationStep)
327     {
328         return std::make_unique<SignallerCallback>(
329                 [this](Step step, Time /*unused*/) { virialReductionStep_ = step; });
330     }
331     return nullptr;
332 }
333
334 template<ComputeGlobalsAlgorithm algorithm>
335 SignallerCallbackPtr ComputeGlobalsElement<algorithm>::registerTrajectorySignallerCallback(TrajectoryEvent event)
336 {
337     if (event == TrajectoryEvent::EnergyWritingStep)
338     {
339         return std::make_unique<SignallerCallback>(
340                 [this](Step step, Time /*unused*/) { energyReductionStep_ = step; });
341     }
342     return nullptr;
343 }
344
345 //! Explicit template instantiation
346 //! @{
347 template class ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>;
348 template class ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>;
349 //! @}
350 } // namespace gmx