83db804cc1b26f11e67c1c46a669d5bd2a342f4e
[alexxy/gromacs.git] / src / gromacs / mdtypes / observablesreducer.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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
37  * Defines gmx::ObservablesReducer and its builder
38  *
39  * These are defined in the same translation unit so that the Impl
40  * object of the ObservablesReducer can be built by the builder.
41  *
42  * \inlibraryapi
43  * \ingroup module_mdtypes
44  */
45 #include "gmxpre.h"
46
47 #include "observablesreducer.h"
48
49 #include <algorithm>
50 #include <numeric>
51
52 #include "gromacs/utility/arrayref.h"
53 #include "gromacs/utility/gmxassert.h"
54
55 namespace gmx
56 {
57
58 //! Impl class for ObservablesReducer
59 class ObservablesReducer::Impl
60 {
61 public:
62     //! Constructor
63     Impl(std::vector<double>&&                                            communicationBuffer,
64          std::vector<ObservablesReducerBuilder::CallbackAfterReduction>&& registeredCallbacks) :
65         communicationBuffer_(std::move(communicationBuffer)),
66         registeredCallbacks_(std::move(registeredCallbacks))
67     {
68     }
69
70     /*! \brief May be called by any subscribed module, via the callback
71      * supplied by ObservablesReducerBuilder::build().
72      *
73      * If this method has been called on this rank since the last call
74      * to reduceComplete() with a \c requirement to communicate soon,
75      * then it will make the communication buffer available via
76      * communicationBuffer() to \c compute_globals(), so it can copy
77      * it to the buffer it uses for MPI communication.
78      *
79      * It is the subscribers' responsibility to coordinate so that
80      * all subscribers on all ranks agree on the need to
81      * communicate, e.g. by orchestating communication based on
82      * the current step number or a previous message.
83      *
84      * Does not check that the callback corresponds to a module that
85      * subscribed to the builder().
86      *
87      * Returns the status of the ObservablesReducer about whether
88      * reduction has already been called on this step.
89      */
90     ObservablesReducerStatus requireReduction(int callbackIndex, ReductionRequirement requirement)
91     {
92         if (requirement == ReductionRequirement::Soon)
93         {
94             reduceSoon_ = true;
95         }
96         callbacksAfterReduction_.push_back(callbackIndex);
97         return status_;
98     }
99
100     /*! \brief Storage for communication buffer
101      *
102      * Must never be resized, because that would potentially
103      * invalidate views of it held by the subscribers. */
104     std::vector<double> communicationBuffer_;
105     //! Registered callbacks that might be used after communication
106     std::vector<ObservablesReducerBuilder::CallbackAfterReduction> registeredCallbacks_;
107     //! Indices into registeredCallbacks_ of callbacks to use after the next reduction
108     std::vector<int> callbacksAfterReduction_;
109     /*! \brief Whether the reduction will occur soon because a
110      * module required it
111      *
112      * "Soon" means this step or next, depending when during the step
113      * it was required, as there is only one point during a normal
114      * simulator step where observables reduction might occur. */
115     bool reduceSoon_ = false;
116     //! Whether reduction has taken place this step
117     ObservablesReducerStatus status_ = ObservablesReducerStatus::ReadyToReduce;
118 };
119
120 ObservablesReducer::ObservablesReducer(std::unique_ptr<Impl> impl) : impl_(std::move(impl)) {}
121
122 ObservablesReducer::ObservablesReducer(ObservablesReducer&&) noexcept = default;
123
124 ObservablesReducer::~ObservablesReducer() = default;
125
126 ObservablesReducer& ObservablesReducer::operator=(ObservablesReducer&& other) noexcept
127 {
128     impl_ = std::move(other.impl_);
129     return *this;
130 }
131
132 ArrayRef<double> ObservablesReducer::communicationBuffer()
133 {
134     if (!impl_->reduceSoon_)
135     {
136         // Nothing to reduce
137         return {};
138     }
139     return impl_->communicationBuffer_;
140 }
141
142 void ObservablesReducer::reductionComplete(Step step)
143 {
144     impl_->status_ = ObservablesReducerStatus::AlreadyReducedThisStep;
145     for (auto& callbackIndex : impl_->callbacksAfterReduction_)
146     {
147         impl_->registeredCallbacks_[callbackIndex](step);
148     }
149     // Prepare for the next reduction
150     std::fill(impl_->communicationBuffer_.begin(), impl_->communicationBuffer_.end(), 0.0);
151     impl_->callbacksAfterReduction_.clear();
152     impl_->reduceSoon_ = false;
153 }
154
155 void ObservablesReducer::stepComplete()
156 {
157     impl_->status_ = ObservablesReducerStatus::ReadyToReduce;
158 }
159
160 //! Impl class for ObservablesReducerBuilder
161 class ObservablesReducerBuilder::Impl
162 {
163 public:
164     //! Data required to set up a subscription
165     struct Subscription
166     {
167         //! Number of doubles required to reduce
168         int sizeRequired;
169         //! The callback to notify of the view and future callback to require reduction
170         ObservablesReducerBuilder::CallbackFromBuilder callbackFromBuilder;
171         //! The callback later used by ObservablesReducer after reduction events
172         ObservablesReducerBuilder::CallbackAfterReduction callbackAfterReduction;
173     };
174
175     //! Contains all subscriptions received
176     std::vector<Subscription> subscriptions_;
177     //! Whether build() has already been called on the owner object
178     bool buildHasBeenCalled_ = false;
179 };
180
181 ObservablesReducerBuilder::ObservablesReducerBuilder() : impl_(std::make_unique<Impl>()) {}
182
183 ObservablesReducerBuilder::ObservablesReducerBuilder(ObservablesReducerBuilder&&) noexcept = default;
184
185 ObservablesReducerBuilder& ObservablesReducerBuilder::operator=(ObservablesReducerBuilder&& other) noexcept
186 {
187     impl_ = std::move(other.impl_);
188     return *this;
189 }
190
191 ObservablesReducerBuilder::~ObservablesReducerBuilder() = default;
192
193 void ObservablesReducerBuilder::addSubscriber(const int                sizeRequired,
194                                               CallbackFromBuilder&&    callbackFromBuilder,
195                                               CallbackAfterReduction&& callbackAfterReduction)
196 {
197     GMX_RELEASE_ASSERT(!impl_->buildHasBeenCalled_,
198                        "Cannot add subscribers to a builder once build() has been called");
199     impl_->subscriptions_.emplace_back(Impl::Subscription{
200             sizeRequired, std::move(callbackFromBuilder), std::move(callbackAfterReduction) });
201 }
202
203 ObservablesReducer ObservablesReducerBuilder::build()
204 {
205     GMX_RELEASE_ASSERT(!impl_->buildHasBeenCalled_,
206                        "Cannot build ObservablesReducer again from the same builder");
207
208     // Prepare the communication buffer
209     const int           totalSizeRequired = std::accumulate(impl_->subscriptions_.begin(),
210                                                   impl_->subscriptions_.end(),
211                                                   0.0,
212                                                   [](int subtotal, const auto& subscription) {
213                                                       return subtotal + subscription.sizeRequired;
214                                                   });
215     std::vector<double> communicationBuffer(totalSizeRequired);
216     // Set up a view of the communication buffer that we can use after
217     // ownership has been transferred to the impl object.
218     const ArrayRef<double> bufferView(communicationBuffer);
219
220     std::vector<ObservablesReducerBuilder::CallbackAfterReduction> registeredCallbacksAfterReduction;
221     registeredCallbacksAfterReduction.reserve(impl_->subscriptions_.size());
222     for (const Impl::Subscription& subscription : impl_->subscriptions_)
223     {
224         registeredCallbacksAfterReduction.emplace_back(subscription.callbackAfterReduction);
225     }
226
227     // Make the impl object so we can set up callbacks to it. This is
228     // safe because the impl object is allocated on the heap, so we
229     // have a stable value to share with the subscribers.
230     auto implPtr = std::make_unique<ObservablesReducer::Impl>(
231             std::move(communicationBuffer), std::move(registeredCallbacksAfterReduction));
232     auto* impl = implPtr.get();
233
234     // Then use the impl object to make the real one. Note that there
235     // is no need for an ObservablesReducer to keep track of all the
236     // modules that might notify it, so it doesn't do that.
237     ObservablesReducer observablesReducer(std::move(implPtr));
238
239     // Now let the subscribers know how to require reduction in
240     // future, and which memory they should use for input and output.
241     size_t start                         = 0;
242     int    indexToCallbackAfterReduction = 0;
243     for (const Impl::Subscription& subscription : impl_->subscriptions_)
244     {
245         // Construct the callback that will hereafter be owned by the
246         // subscriber.
247         CallbackToRequireReduction callbackToRequireReduction =
248                 [impl, indexToCallbackAfterReduction](ReductionRequirement requirement) {
249                     return impl->requireReduction(indexToCallbackAfterReduction, requirement);
250                 };
251         subscription.callbackFromBuilder(std::move(callbackToRequireReduction),
252                                          bufferView.subArray(start, subscription.sizeRequired));
253         start += subscription.sizeRequired;
254         ++indexToCallbackAfterReduction;
255     }
256
257     impl_->buildHasBeenCalled_ = true;
258
259     return observablesReducer;
260 }
261
262 } // namespace gmx