Use ObservablesReducer for check of DD bonded interaction count.
[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,2021, 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/domdec.h"
47 #include "gromacs/gmxlib/network.h"
48 #include "gromacs/gmxlib/nrnb.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/mdlib/constr.h"
51 #include "gromacs/mdlib/md_support.h"
52 #include "gromacs/mdlib/mdatoms.h"
53 #include "gromacs/mdlib/stat.h"
54 #include "gromacs/mdlib/update.h"
55 #include "gromacs/mdlib/vcm.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/group.h"
58 #include "gromacs/mdtypes/inputrec.h"
59 #include "gromacs/mdtypes/mdatom.h"
60 #include "gromacs/topology/topology.h"
61
62 #include "freeenergyperturbationdata.h"
63 #include "modularsimulator.h"
64 #include "simulatoralgorithm.h"
65
66 namespace gmx
67 {
68 template<ComputeGlobalsAlgorithm algorithm>
69 ComputeGlobalsElement<algorithm>::ComputeGlobalsElement(StatePropagatorData* statePropagatorData,
70                                                         EnergyData*          energyData,
71                                                         FreeEnergyPerturbationData* freeEnergyPerturbationData,
72                                                         SimulationSignals*  signals,
73                                                         int                 nstglobalcomm,
74                                                         FILE*               fplog,
75                                                         const MDLogger&     mdlog,
76                                                         t_commrec*          cr,
77                                                         const t_inputrec*   inputrec,
78                                                         const MDAtoms*      mdAtoms,
79                                                         t_nrnb*             nrnb,
80                                                         gmx_wallcycle*      wcycle,
81                                                         t_forcerec*         fr,
82                                                         const gmx_mtop_t&   global_top,
83                                                         Constraints*        constr,
84                                                         ObservablesReducer* observablesReducer) :
85     energyReductionStep_(-1),
86     virialReductionStep_(-1),
87     vvSchedulingStep_(-1),
88     doStopCM_(inputrec->comm_mode != ComRemovalAlgorithm::No),
89     nstcomm_(inputrec->nstcomm),
90     nstglobalcomm_(nstglobalcomm),
91     lastStep_(inputrec->nsteps + inputrec->init_step),
92     initStep_(inputrec->init_step),
93     nullSignaller_(std::make_unique<SimulationSignaller>(nullptr, nullptr, nullptr, false, false)),
94     statePropagatorData_(statePropagatorData),
95     energyData_(energyData),
96     freeEnergyPerturbationData_(freeEnergyPerturbationData),
97     vcm_(global_top.groups, *inputrec),
98     signals_(signals),
99     fplog_(fplog),
100     mdlog_(mdlog),
101     cr_(cr),
102     inputrec_(inputrec),
103     top_global_(global_top),
104     mdAtoms_(mdAtoms),
105     constr_(constr),
106     nrnb_(nrnb),
107     wcycle_(wcycle),
108     fr_(fr),
109     observablesReducer_(observablesReducer)
110 {
111     reportComRemovalInfo(fplog, vcm_);
112     gstat_ = global_stat_init(inputrec_);
113 }
114
115 template<ComputeGlobalsAlgorithm algorithm>
116 ComputeGlobalsElement<algorithm>::~ComputeGlobalsElement()
117 {
118     global_stat_destroy(gstat_);
119 }
120
121 template<ComputeGlobalsAlgorithm algorithm>
122 void ComputeGlobalsElement<algorithm>::elementSetup()
123 {
124     if (doStopCM_ && !inputrec_->bContinuation)
125     {
126         // To minimize communication, compute_globals computes the COM velocity
127         // and the kinetic energy for the velocities without COM motion removed.
128         // Thus to get the kinetic energy without the COM contribution, we need
129         // to call compute_globals twice.
130
131         compute(-1, CGLO_GSTAT | CGLO_STOPCM, nullSignaller_.get(), false, true);
132
133         auto v = statePropagatorData_->velocitiesView();
134         // At initialization, do not pass x with acceleration-correction mode
135         // to avoid (incorrect) correction of the initial coordinates.
136         auto x = vcm_.mode == ComRemovalAlgorithm::LinearAccelerationCorrection
137                          ? ArrayRefWithPadding<RVec>()
138                          : statePropagatorData_->positionsView();
139         process_and_stopcm_grp(
140                 fplog_, &vcm_, *mdAtoms_->mdatoms(), x.unpaddedArrayRef(), v.unpaddedArrayRef());
141         inc_nrnb(nrnb_, eNR_STOPCM, mdAtoms_->mdatoms()->homenr);
142     }
143
144     unsigned int cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
145                                | (energyData_->hasReadEkinFromCheckpoint() ? 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(energyData_->ekindata()->tcstat[i].ekinh, energyData_->ekindata()->tcstat[i].ekinh_old);
158     }
159 }
160
161 template<ComputeGlobalsAlgorithm algorithm>
162 void ComputeGlobalsElement<algorithm>::scheduleTask(Step                       step,
163                                                     Time gmx_unused            time,
164                                                     const RegisterRunFunction& registerRunFunction)
165 {
166     const bool needComReduction    = doStopCM_ && do_per_step(step, nstcomm_);
167     const bool needGlobalReduction = step == energyReductionStep_ || step == virialReductionStep_
168                                      || needComReduction || do_per_step(step, nstglobalcomm_)
169                                      || (EI_VV(inputrec_->eI) && inputrecNvtTrotter(inputrec_)
170                                          && do_per_step(step - 1, 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    = (needGlobalReduction ? CGLO_GSTAT : 0) | (doEnergy ? CGLO_ENERGY : 0)
193                     | (needComReduction ? CGLO_STOPCM : 0) | CGLO_TEMPERATURE | CGLO_PRESSURE
194                     | CGLO_CONSTRAINT;
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>(
207                 signals_, cr_, nullptr, doInterSimSignal, doIntraSimSignal);
208
209         registerRunFunction([this, step, flags, signaller = std::move(signaller)]() {
210             compute(step, flags, signaller.get(), true);
211         });
212     }
213     else if (algorithm == ComputeGlobalsAlgorithm::VelocityVerlet)
214     {
215         // For VV, we schedule two calls to compute globals per step.
216         if (step != vvSchedulingStep_)
217         {
218             // This is the first scheduling call for this step (positions & velocities at full time
219             // step) Set this as the current scheduling step
220             vvSchedulingStep_ = step;
221
222             // For vv, the state at the beginning of the step is positions at time t, velocities at time t - dt/2
223             // The first velocity propagation (+dt/2) therefore actually corresponds to the previous step.
224             // So we need information from the last step in the first half of the integration
225             if (!needGlobalReduction && !do_per_step(step - 1, nstglobalcomm_))
226             {
227                 return;
228             }
229
230             const bool doTemperature = step != initStep_ || inputrec_->bContinuation;
231             const bool doEnergy      = step == energyReductionStep_;
232
233             int flags = (needGlobalReduction ? CGLO_GSTAT : 0) | (doEnergy ? CGLO_ENERGY : 0)
234                         | (doTemperature ? CGLO_TEMPERATURE : 0) | CGLO_PRESSURE | CGLO_CONSTRAINT
235                         | (needComReduction ? CGLO_STOPCM : 0) | CGLO_SCALEEKIN;
236
237             registerRunFunction(
238                     [this, step, flags]() { compute(step, flags, nullSignaller_.get(), false); });
239         }
240         else
241         {
242             // second call to compute_globals for this step
243             // Reset the scheduling step to avoid confusion if scheduling needs
244             // to be repeated (in case of unexpected simulation termination)
245             vvSchedulingStep_ = -1;
246
247             if (!needGlobalReduction)
248             {
249                 return;
250             }
251             int flags = CGLO_GSTAT | CGLO_CONSTRAINT;
252
253             // Since we're already communicating at this step, we
254             // can propagate intra-simulation signals. Note that
255             // check_nstglobalcomm has the responsibility for
256             // choosing the value of nstglobalcomm which satisfies
257             // the need of the different signallers.
258             const bool doIntraSimSignal = true;
259             // Disable functionality
260             const bool doInterSimSignal = false;
261
262             auto signaller = std::make_shared<SimulationSignaller>(
263                     signals_, cr_, nullptr, doInterSimSignal, doIntraSimSignal);
264
265             registerRunFunction([this, step, flags, signaller = std::move(signaller)]() {
266                 compute(step, flags, signaller.get(), true);
267             });
268         }
269     }
270 }
271
272 template<ComputeGlobalsAlgorithm algorithm>
273 void ComputeGlobalsElement<algorithm>::compute(gmx::Step            step,
274                                                unsigned int         flags,
275                                                SimulationSignaller* signaller,
276                                                bool                 useLastBox,
277                                                bool                 isInit)
278 {
279     auto        x       = statePropagatorData_->positionsView().unpaddedArrayRef();
280     auto        v       = statePropagatorData_->velocitiesView().unpaddedArrayRef();
281     const auto* box     = statePropagatorData_->constBox();
282     const auto* lastbox = useLastBox ? statePropagatorData_->constPreviousBox()
283                                      : statePropagatorData_->constBox();
284
285     compute_globals(gstat_,
286                     cr_,
287                     inputrec_,
288                     fr_,
289                     energyData_->ekindata(),
290                     x,
291                     v,
292                     box,
293                     mdAtoms_->mdatoms(),
294                     nrnb_,
295                     &vcm_,
296                     step != -1 ? wcycle_ : nullptr,
297                     energyData_->enerdata(),
298                     energyData_->forceVirial(step),
299                     energyData_->constraintVirial(step),
300                     energyData_->totalVirial(step),
301                     energyData_->pressure(step),
302                     (((flags & CGLO_ENERGY) != 0) && constr_ != nullptr) ? constr_->rmsdData()
303                                                                          : gmx::ArrayRef<real>{},
304                     signaller,
305                     lastbox,
306                     energyData_->needToSumEkinhOld(),
307                     flags,
308                     step,
309                     observablesReducer_);
310     if (flags & CGLO_STOPCM && !isInit)
311     {
312         process_and_stopcm_grp(fplog_, &vcm_, *mdAtoms_->mdatoms(), x, v);
313         inc_nrnb(nrnb_, eNR_STOPCM, mdAtoms_->mdatoms()->homenr);
314     }
315 }
316
317 template<ComputeGlobalsAlgorithm algorithm>
318 std::optional<SignallerCallback> ComputeGlobalsElement<algorithm>::registerEnergyCallback(EnergySignallerEvent event)
319 {
320     if (event == EnergySignallerEvent::EnergyCalculationStep)
321     {
322         return [this](Step step, Time /*unused*/) { energyReductionStep_ = step; };
323     }
324     if (event == EnergySignallerEvent::VirialCalculationStep)
325     {
326         return [this](Step step, Time /*unused*/) { virialReductionStep_ = step; };
327     }
328     return std::nullopt;
329 }
330
331 template<ComputeGlobalsAlgorithm algorithm>
332 std::optional<SignallerCallback>
333 ComputeGlobalsElement<algorithm>::registerTrajectorySignallerCallback(TrajectoryEvent event)
334 {
335     if (event == TrajectoryEvent::EnergyWritingStep)
336     {
337         return [this](Step step, Time /*unused*/) { energyReductionStep_ = step; };
338     }
339     return std::nullopt;
340 }
341
342 //! Explicit template instantiation
343 //! \{
344 template class ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>;
345 template class ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>;
346 //! \}
347
348 template<>
349 ISimulatorElement* ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>::getElementPointerImpl(
350         LegacySimulatorData*                    legacySimulatorData,
351         ModularSimulatorAlgorithmBuilderHelper* builderHelper,
352         StatePropagatorData*                    statePropagatorData,
353         EnergyData*                             energyData,
354         FreeEnergyPerturbationData*             freeEnergyPerturbationData,
355         GlobalCommunicationHelper*              globalCommunicationHelper,
356         ObservablesReducer*                     observablesReducer)
357 {
358     auto* element = builderHelper->storeElement(
359             std::make_unique<ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>>(
360                     statePropagatorData,
361                     energyData,
362                     freeEnergyPerturbationData,
363                     globalCommunicationHelper->simulationSignals(),
364                     globalCommunicationHelper->nstglobalcomm(),
365                     legacySimulatorData->fplog,
366                     legacySimulatorData->mdlog,
367                     legacySimulatorData->cr,
368                     legacySimulatorData->inputrec,
369                     legacySimulatorData->mdAtoms,
370                     legacySimulatorData->nrnb,
371                     legacySimulatorData->wcycle,
372                     legacySimulatorData->fr,
373                     legacySimulatorData->top_global,
374                     legacySimulatorData->constr,
375                     observablesReducer));
376
377     return element;
378 }
379
380 template<>
381 ISimulatorElement* ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>::getElementPointerImpl(
382         LegacySimulatorData*                    simulator,
383         ModularSimulatorAlgorithmBuilderHelper* builderHelper,
384         StatePropagatorData*                    statePropagatorData,
385         EnergyData*                             energyData,
386         FreeEnergyPerturbationData*             freeEnergyPerturbationData,
387         GlobalCommunicationHelper*              globalCommunicationHelper,
388         ObservablesReducer*                     observablesReducer)
389 {
390     // We allow this element to be added multiple times to the call list, but we only want one
391     // actual element built
392     static const std::string key("vvComputeGlobalsElement");
393
394     const std::optional<std::any> cachedValue = builderHelper->builderData(key);
395
396     if (cachedValue)
397     {
398         return std::any_cast<ISimulatorElement*>(cachedValue.value());
399     }
400     else
401     {
402         ISimulatorElement* vvComputeGlobalsElement = builderHelper->storeElement(
403                 std::make_unique<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>(
404                         statePropagatorData,
405                         energyData,
406                         freeEnergyPerturbationData,
407                         globalCommunicationHelper->simulationSignals(),
408                         globalCommunicationHelper->nstglobalcomm(),
409                         simulator->fplog,
410                         simulator->mdlog,
411                         simulator->cr,
412                         simulator->inputrec,
413                         simulator->mdAtoms,
414                         simulator->nrnb,
415                         simulator->wcycle,
416                         simulator->fr,
417                         simulator->top_global,
418                         simulator->constr,
419                         observablesReducer));
420         builderHelper->storeBuilderData(key, vvComputeGlobalsElement);
421         return vvComputeGlobalsElement;
422     }
423 }
424 } // namespace gmx