Fixes for clang-tidy-9
[alexxy/gromacs.git] / src / gromacs / restraint / restraintmdmodule_impl.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,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
36 #ifndef GROMACS_RESTRAINTMDMODULE_IMPL_H
37 #define GROMACS_RESTRAINTMDMODULE_IMPL_H
38
39 /*! \libinternal \file
40  * \brief Implementation details for RestraintMDModule
41  *
42  * \author M. Eric Irrgang <ericirrgang@gmail.com>
43  *
44  * \ingroup module_restraint
45  */
46
47 #include <iostream>
48 #include <mutex>
49
50 #include "gromacs/domdec/domdec_struct.h"
51 #include "gromacs/domdec/ga2la.h"
52 #include "gromacs/gmxlib/network.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdtypes/commrec.h"
55 #include "gromacs/mdtypes/iforceprovider.h"
56 #include "gromacs/mdtypes/imdmodule.h"
57 #include "gromacs/mdtypes/imdoutputprovider.h"
58 #include "gromacs/mdtypes/imdpoptionprovider.h"
59 #include "gromacs/mdtypes/mdatom.h"
60 #include "gromacs/pbcutil/pbc.h"
61 #include "gromacs/restraint/restraintpotential.h"
62 #include "gromacs/utility/arrayref.h"
63
64 namespace gmx
65 {
66
67 /*! \libinternal
68  * \brief Abstraction for a restraint interaction site.
69  *
70  * A restraint may operate on a single atom or some other entity, such as a selection of atoms.
71  * The Restraint implementation is very independent from how coordinates are provided or what they mean.
72  *
73  * First implementation can only represent single atoms in a global context.
74  *
75  * Ultimately, this should be replaced with a more universal facility for acting
76  * on distributed atom data or simple transformations thereof.
77  *
78  * \inlibraryapi
79  * \ingroup module_restraint
80  */
81 class Site
82 {
83 public:
84     /*! \brief Construct from global atom indices
85      *
86      * \param globalIndex Atom index in the global state (as input to the simulation)
87      */
88     explicit Site(int globalIndex) : index_(globalIndex), r_(0, 0, 0) {}
89
90     /*!
91      * \brief Explicitly define copies.
92      *
93      * Implicit definition is not possible because of the mutex member,
94      * and a copy constructor is necessary to use Site in a std::vector. Really, we should make
95      * a copy point to the same implementation object to reuse its cache.
96
97      */
98     Site(const Site& site) : index_(site.index_), r_(site.r_) {}
99
100     /*! \brief Disallow assignment.
101      *
102      * Assignment doesn't make sense because it implies that a site's meaning is fuzzy.
103      * If the definition of a site is changing, just make a new site.
104      * There's nothing to be gained by reusing one or by creating it uninitialized.
105      */
106     Site& operator=(const Site&) = delete;
107
108     /*!
109      * \brief Get the global atom index of an atomic site.
110      *
111      * \return global index provided at construction.
112      *
113      */
114     int index() const { return index_; }
115
116     /*!
117      * \brief Get the position of this site at time t.
118      *
119      * \param cr Communications record.
120      * \param nx Number of locally available atoms (size of local atom data arrays)
121      * \param x Array of locally available atom coordinates.
122      * \param t the current time.
123      * \return position vector.
124      *
125      * \internal
126      * By providing the current time, we can cache results in order to use them once per timestep.
127      * In the long term, we would prefer to also allow client code to preregister interest in a
128      * position at a given time, or issue "futures".
129      */
130     RVec centerOfMass(const t_commrec& cr, size_t nx, ArrayRef<const RVec> x, double gmx_unused t)
131     {
132         // Center of mass to return for the site. Currently the only form of site
133         // implemented is as a global atomic coordinate.
134         gmx::RVec r = { 0, 0, 0 };
135         if (DOMAINDECOMP(&cr)) // Domain decomposition
136         {
137             // Get global-to-local indexing structure
138             auto crossRef = cr.dd->ga2la;
139             GMX_ASSERT(crossRef, "Domain decomposition must provide global/local cross-reference.");
140             if (const auto localIndex = crossRef->findHome(index_))
141             {
142                 GMX_ASSERT(localIndex,
143                            "Expect not to reach this point if findHome does not find index_.");
144                 GMX_ASSERT(*localIndex < static_cast<decltype(*localIndex)>(nx),
145                            "We assume that the local index cannot be larger than the number of "
146                            "atoms.");
147                 GMX_ASSERT(*localIndex >= 0, "localIndex is a signed type, but is assumed >0.");
148                 // If atom is local, get its location
149                 copy_rvec(x[*localIndex], r);
150             }
151             else
152             {
153                 // Nothing to contribute on this rank. Leave position == [0,0,0].
154             }
155             // AllReduce across the ranks of the simulation to get the center-of-mass
156             // of the site locally available everywhere. For single-atom sites, this
157             // is trivial: exactly one rank should have a non-zero position.
158             // For future multi-atom selections,
159             // we will receive weighted center-of-mass contributions from
160             // each rank and combine to get the global center of mass.
161             // \todo use generalized "pull group" facility when available.
162             std::array<double, 3> buffer{ { r[0], r[1], r[2] } };
163             // This should be an all-reduce sum, which gmx_sumd appears to be.
164             gmx_sumd(3, buffer.data(), &cr);
165             r[0] = static_cast<real>(buffer[0]);
166             r[1] = static_cast<real>(buffer[1]);
167             r[2] = static_cast<real>(buffer[2]);
168
169         } // end domain decomposition branch
170         else
171         {
172             // No DD so all atoms are local.
173             copy_rvec(x[index_], r);
174             (void)nx;
175         }
176         // Update cache and cache status.
177         copy_rvec(r, r_);
178
179         return r_;
180     }
181
182 private:
183     /*!
184      * \brief Global index of the single-atom site.
185      *
186      * \todo This class should be a specialization of a more general Site data source.
187      * \todo use LocalAtomSet
188      */
189     const int index_;
190
191     /*!
192      * \brief Last known value of the center-of-mass.
193      *
194      * Updated with centerOfMass().
195      */
196     RVec r_;
197 };
198
199 /*! \internal
200  * \brief Provide IForceProvider for RestraintMDModuleImpl
201  *
202  * Adapter class from IForceProvider to IRestraintPotential.
203  * Objects of this type are uniquely owned by instances of RestraintMDModuleImpl. The object will
204  * dispatch calls to IForceProvider->calculateForces() to the functor managed by
205  * RestraintMDModuleImpl. \ingroup module_restraint
206  */
207 class RestraintForceProvider final : public gmx::IForceProvider
208 {
209 public:
210     /*!
211      * \brief Can only be constructed when initialized from a restraint.
212      */
213     RestraintForceProvider() = delete;
214
215     ~RestraintForceProvider() = default;
216
217     /*!
218      * \brief RAII construction with an IRestraintPotential
219      *
220      * Note, this object must outlive the pointer that will be provided to ForceProviders.
221      * \param restraint handle to an object providing restraint potential calculation
222      * \param sites List of atomic site indices
223      */
224     explicit RestraintForceProvider(std::shared_ptr<gmx::IRestraintPotential> restraint,
225                                     const std::vector<int>&                   sites);
226
227     /*!
228      * \brief Implement the IForceProvider interface.
229      *
230      * Update the force array with restraint contribution(s) for local atoms.
231      *
232      * RestraintForceProvider is implemented with the assumption that few
233      * restraints apply to many atoms.
234      * That is, the number of restraints affecting a large number of atoms is small,
235      * though there may be several restraints that apply to few atoms each.
236      * Under this assumption, it is considered computationally inexpensive to iterate
237      * over restraints in an outer loop and iterate over atoms within each restraint.
238      * This would be an invalid assumption if, say, several restraints applied
239      * to an entire membrane or the entire solvent group.
240      *
241      * If the assumption causes performance problems, we can look for a good
242      * way to reduce from several restraints in a single pass or a very
243      * lightweight way to determine whether a given restraint applies to a given atom.
244      * There is also the notion in the pulling code of a limited number of
245      * "pull groups" used by the "pull coordinates".
246      * The right optimization will depend on how the code is being used.
247      *
248      * Call the evaluator(s) for the restraints for the configured sites.
249      * Forces are applied to atoms in the first and last site listed.
250      * Intermediate sites are used as reference coordinates when the relevant
251      * vector between sites is on the order of half a box length or otherwise
252      * ambiguous in the case of periodic boundary conditions.
253      */
254     void calculateForces(const ForceProviderInput& forceProviderInput,
255                          ForceProviderOutput*      forceProviderOutput) override;
256
257 private:
258     std::shared_ptr<gmx::IRestraintPotential> restraint_;
259     std::vector<Site>                         sites_;
260 };
261
262 /*! \internal
263  * \brief IMDModule implementation for RestraintMDModule.
264  *
265  * Provides IMDModule interface.
266  *
267  * \ingroup module_restraint
268  */
269 class RestraintMDModuleImpl final
270 {
271 public:
272     RestraintMDModuleImpl() = delete;
273     /*!
274      * \brief Wrap an object implementing IRestraintPotential
275      *
276      * \param restraint handle to restraint to wrap.
277      * \param sites list of sites for framework to process for restraint force calculator.
278      */
279     RestraintMDModuleImpl(std::shared_ptr<gmx::IRestraintPotential> restraint,
280                           const std::vector<int>&                   sites);
281
282     /*!
283      * \brief Allow moves.
284      *
285      * \{
286      */
287     RestraintMDModuleImpl(RestraintMDModuleImpl&&) noexcept = default;
288     RestraintMDModuleImpl& operator=(RestraintMDModuleImpl&&) noexcept = default;
289     /*! \} */
290
291     ~RestraintMDModuleImpl();
292
293     /*!
294      * \brief Implement IMDModule interface.
295      *
296      * \param forceProviders force module manager in the force record that will call this.
297      *
298      * The calling code must ensure that this object stays alive as long as forceProviders needs
299      * the RestraintForceProvider, since forceProviders can't. Typically that is the duration of a do_md() call.
300      */
301     void initForceProviders(ForceProviders* forceProviders);
302
303     //! handle to RestraintForceProvider implementation
304     std::unique_ptr<RestraintForceProvider> forceProvider_;
305 };
306
307 } // end namespace gmx
308
309 #endif // GROMACS_RESTRAINTMDMODULE_IMPL_H