Add nbnxm GridSet class
[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 An object holding a set of search grids for the local + non-local DD zones
75  */
76 class GridSet
77 {
78     public:
79         //! Constructs a grid set for 1 or multiple DD zones, when numDDCells!=nullptr
80         GridSet(const ivec   *numDDCells,
81                 PairlistType  pairlistType,
82                 bool          haveFep,
83                 int           numThreads);
84
85         //! Puts the atoms in \p ddZone on the grid and copies the coordinates to \p nbat
86         void putOnGrid(const matrix                    box,
87                        int                             ddZone,
88                        const rvec                      lowerCorner,
89                        const rvec                      upperCorner,
90                        const gmx::UpdateGroupsCog     *updateGroupsCog,
91                        int                             atomStart,
92                        int                             atomEnd,
93                        real                            atomDensity,
94                        const int                      *atinfo,
95                        gmx::ArrayRef<const gmx::RVec>  x,
96                        int                             numAtomsMoved,
97                        const int                      *move,
98                        nbnxn_atomdata_t               *nbat);
99
100         //! Returns the number of cells along x and y for the local grid
101         void getLocalNumCells(int *numCellsX,
102                               int *numCellsY) const
103         {
104             *numCellsX = grids_[0].dimensions().numCells[XX];
105             *numCellsY = grids_[0].dimensions().numCells[YY];
106         }
107
108         //! Returns the total number of atoms in the grid set, including padding
109         int numGridAtomsTotal() const
110         {
111             return grids_.back().atomIndexEnd();
112         }
113
114         //! Returns the number of local real atoms, i.e. without padded atoms
115         int numRealAtomsLocal() const
116         {
117             return numRealAtomsLocal_;
118         }
119
120         //! Returns the number of total real atoms, i.e. without padded atoms
121         int numRealAtomsTotal() const
122         {
123             return numRealAtomsTotal_;
124         }
125
126         //! Returns the atom order on the grid for the local atoms
127         gmx::ArrayRef<const int> getLocalAtomorder() const
128         {
129             /* Return the atom order for the home cell (index 0) */
130             const int numIndices = grids_[0].atomIndexEnd() - grids_[0].firstAtomInColumn(0);
131
132             return gmx::constArrayRefFromArray(atomIndices_.data(), numIndices);
133         }
134
135         //! Sets the order of the local atoms to the order grid atom ordering
136         void setLocalAtomOrder();
137
138         //! Returns the list of grids
139         gmx::ArrayRef<const Grid> grids() const
140         {
141             return grids_;
142         }
143
144         //! Returns the grid atom indices covering all grids
145         gmx::ArrayRef<const int> cells() const
146         {
147             return cells_;
148         }
149
150         //! Returns the grid atom indices covering all grids
151         gmx::ArrayRef<const int> atomIndices() const
152         {
153             return atomIndices_;
154         }
155
156         //! Returns whether we have perturbed non-bonded interactions
157         bool haveFep() const
158         {
159             return haveFep_;
160         }
161
162         //! Returns the unit cell in \p box
163         void getBox(matrix box) const
164         {
165             copy_mat(box_, box);
166         }
167
168     private:
169         //! Returns collection of the data that covers all grids
170         const GridSetData getGridSetData()
171         {
172             GridSetData gridSetData = { cells_, atomIndices_, haveFep_ };
173
174             return gridSetData;
175         }
176
177         /* Data members */
178         //! The search grids
179         std::vector<Grid>     grids_;
180         //! The actual cell indices for all atoms, covering all grids
181         std::vector<int>      cells_;
182         //! The actual array of atom indices, covering all grids
183         std::vector<int>      atomIndices_;
184         //! Tells whether we have perturbed non-bonded interactions
185         bool                  haveFep_;
186         //! The periodic unit-cell
187         matrix                box_;
188         //! The number of local real atoms, i.e. without padded atoms, local atoms: 0 to numAtomsLocal_
189         int                   numRealAtomsLocal_;
190         //! The total number of real atoms, i.e. without padded atoms
191         int                   numRealAtomsTotal_;
192         //! Working data for constructing a single grid, one entry per thread
193         std::vector<GridWork> gridWork_;
194 };
195
196 } // namespace Nbnxm
197
198 #endif