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