Remove charge groups from domdec and localtop
[alexxy/gromacs.git] / src / gromacs / nbnxm / gridset.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019, 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 /*! \internal \file
37  *
38  * \brief
39  * Declares the GridSet class.
40  *
41  * This class holds the grids for the local and non-local domain decomposition
42  * zones, as well as the cell and atom data that covers all grids.
43  *
44  * \author Berk Hess <hess@kth.se>
45  * \ingroup module_nbnxm
46  */
47
48 #ifndef GMX_NBNXM_GRIDSET_H
49 #define GMX_NBNXM_GRIDSET_H
50
51 #include <memory>
52 #include <vector>
53
54 #include "gromacs/math/vec.h"
55 #include "gromacs/math/vectypes.h"
56 #include "gromacs/utility/alignedallocator.h"
57 #include "gromacs/utility/arrayref.h"
58
59 #include "grid.h"
60
61
62 struct nbnxn_atomdata_t;
63 enum class PairlistType;
64
65 namespace gmx
66 {
67 class UpdateGroupsCog;
68 }
69
70 namespace Nbnxm
71 {
72
73 /*! \internal
74  * \brief Holds a set of search grids for the local + non-local DD zones
75  */
76 class GridSet
77 {
78     public:
79         /*! \internal
80          * \brief Description of the domain setup: PBC and the connections between domains
81          */
82         struct DomainSetup
83         {
84             //! Constructor, without DD \p numDDCells and \p ddZones should be nullptr
85             DomainSetup(int                       ePBC,
86                         const ivec               *numDDCells,
87                         const gmx_domdec_zones_t *ddZones);
88
89             //! The type of PBC
90             int                       ePBC;
91             //! Are there multiple domains?
92             bool                      haveMultipleDomains;
93             //! Are there multiple domains along each dimension?
94             std::array<bool, DIM>     haveMultipleDomainsPerDim;
95             //! The domain decomposition zone setup
96             const gmx_domdec_zones_t *zones;
97         };
98
99         //! Constructs a grid set for 1 or multiple DD zones, when numDDCells!=nullptr
100         GridSet(int                       ePBC,
101                 const ivec               *numDDCells,
102                 const gmx_domdec_zones_t *ddZones,
103                 PairlistType              pairlistType,
104                 bool                      haveFep,
105                 int                       numThreads);
106
107         //! Puts the atoms in \p ddZone on the grid and copies the coordinates to \p nbat
108         void putOnGrid(const matrix                    box,
109                        int                             ddZone,
110                        const rvec                      lowerCorner,
111                        const rvec                      upperCorner,
112                        const gmx::UpdateGroupsCog     *updateGroupsCog,
113                        int                             atomStart,
114                        int                             atomEnd,
115                        real                            atomDensity,
116                        gmx::ArrayRef<const int>        atomInfo,
117                        gmx::ArrayRef<const gmx::RVec>  x,
118                        int                             numAtomsMoved,
119                        const int                      *move,
120                        nbnxn_atomdata_t               *nbat);
121
122         //! Returns the domain setup
123         const DomainSetup domainSetup() const
124         {
125             return domainSetup_;
126         }
127
128         //! Returns the total number of atoms in the grid set, including padding
129         int numGridAtomsTotal() const
130         {
131             return grids_.back().atomIndexEnd();
132         }
133
134         //! Returns the number of local real atoms, i.e. without padded atoms
135         int numRealAtomsLocal() const
136         {
137             return numRealAtomsLocal_;
138         }
139
140         //! Returns the number of total real atoms, i.e. without padded atoms
141         int numRealAtomsTotal() const
142         {
143             return numRealAtomsTotal_;
144         }
145
146         //! Returns the atom order on the grid for the local atoms
147         gmx::ArrayRef<const int> getLocalAtomorder() const
148         {
149             /* Return the atom order for the home cell (index 0) */
150             const int numIndices = grids_[0].atomIndexEnd() - grids_[0].firstAtomInColumn(0);
151
152             return gmx::constArrayRefFromArray(atomIndices_.data(), numIndices);
153         }
154
155         //! Sets the order of the local atoms to the order grid atom ordering
156         void setLocalAtomOrder();
157
158         //! Returns the list of grids
159         gmx::ArrayRef<const Grid> grids() const
160         {
161             return grids_;
162         }
163
164         //! Returns the grid atom indices covering all grids
165         gmx::ArrayRef<const int> cells() const
166         {
167             return cells_;
168         }
169
170         //! Returns the grid atom indices covering all grids
171         gmx::ArrayRef<const int> atomIndices() const
172         {
173             return atomIndices_;
174         }
175
176         //! Returns whether we have perturbed non-bonded interactions
177         bool haveFep() const
178         {
179             return haveFep_;
180         }
181
182         //! Returns the unit cell in \p box
183         void getBox(matrix box) const
184         {
185             copy_mat(box_, box);
186         }
187
188     private:
189         //! Returns collection of the data that covers all grids
190         const GridSetData getGridSetData()
191         {
192             GridSetData gridSetData = { cells_, atomIndices_, haveFep_ };
193
194             return gridSetData;
195         }
196
197         /* Data members */
198         //! The domain setup
199         DomainSetup           domainSetup_;
200         //! The search grids
201         std::vector<Grid>     grids_;
202         //! The actual cell indices for all atoms, covering all grids
203         std::vector<int>      cells_;
204         //! The actual array of atom indices, covering all grids
205         std::vector<int>      atomIndices_;
206         //! Tells whether we have perturbed non-bonded interactions
207         bool                  haveFep_;
208         //! The periodic unit-cell
209         matrix                box_;
210         //! The number of local real atoms, i.e. without padded atoms, local atoms: 0 to numAtomsLocal_
211         int                   numRealAtomsLocal_;
212         //! The total number of real atoms, i.e. without padded atoms
213         int                   numRealAtomsTotal_;
214         //! Working data for constructing a single grid, one entry per thread
215         std::vector<GridWork> gridWork_;
216 };
217
218 } // namespace Nbnxm
219
220 #endif