3980a37c6852246909cc8111a995d02070838fec
[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 #include "gromacs/utility/range.h"
59
60 #include "grid.h"
61 #include "gridsetdata.h"
62
63
64 struct nbnxn_atomdata_t;
65 enum class PairlistType;
66
67 namespace gmx
68 {
69 class UpdateGroupsCog;
70 }
71
72 namespace Nbnxm
73 {
74
75 /*! \internal
76  * \brief Holds a set of search grids for the local + non-local DD zones
77  *
78  * The are three different possible setups:
79  * - a single grid, this is the standard case without domain decomposition
80  * - one grid for each domain decomposition zone
81  * - with test particle insertion there are two grids, one for the system
82  *   to insert in and one for the molecule that is inserted
83  */
84 class GridSet
85 {
86 public:
87     /*! \internal
88      * \brief Description of the domain setup: PBC and the connections between domains
89      */
90     struct DomainSetup
91     {
92         //! Constructor, without DD \p numDDCells and \p ddZones should be nullptr
93         DomainSetup(int                       ePBC,
94                     bool                      doTestParticleInsertion,
95                     const ivec*               numDDCells,
96                     const gmx_domdec_zones_t* ddZones);
97
98         //! The type of PBC
99         int ePBC;
100         //! Tells whether we are doing test-particle insertion
101         bool doTestParticleInsertion;
102         //! Are there multiple domains?
103         bool haveMultipleDomains;
104         //! Are there multiple domains along each dimension?
105         std::array<bool, DIM> haveMultipleDomainsPerDim;
106         //! The domain decomposition zone setup
107         const gmx_domdec_zones_t* zones;
108     };
109
110     //! Constructs a grid set for 1 or multiple DD zones, when numDDCells!=nullptr
111     GridSet(int                       ePBC,
112             bool                      doTestParticleInsertion,
113             const ivec*               numDDCells,
114             const gmx_domdec_zones_t* ddZones,
115             PairlistType              pairlistType,
116             bool                      haveFep,
117             int                       numThreads,
118             gmx::PinningPolicy        pinningPolicy);
119
120     //! Puts the atoms on the grid with index \p gridIndex and copies the coordinates to \p nbat
121     void putOnGrid(const matrix                   box,
122                    int                            gridIndex,
123                    const rvec                     lowerCorner,
124                    const rvec                     upperCorner,
125                    const gmx::UpdateGroupsCog*    updateGroupsCog,
126                    gmx::Range<int>                atomRange,
127                    real                           atomDensity,
128                    gmx::ArrayRef<const int>       atomInfo,
129                    gmx::ArrayRef<const gmx::RVec> x,
130                    int                            numAtomsMoved,
131                    const int*                     move,
132                    nbnxn_atomdata_t*              nbat);
133
134     //! Returns the domain setup
135     DomainSetup domainSetup() const { return domainSetup_; }
136
137     //! Returns the total number of atoms in the grid set, including padding
138     int numGridAtomsTotal() const { return grids_.back().atomIndexEnd(); }
139
140     //! Returns the number of local real atoms, i.e. without padded atoms
141     int numRealAtomsLocal() const { return numRealAtomsLocal_; }
142
143     //! Returns the number of total real atoms, i.e. without padded atoms
144     int numRealAtomsTotal() const { return numRealAtomsTotal_; }
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 { return grids_; }
160
161     //! Returns the grid atom indices covering all grids
162     gmx::ArrayRef<const int> cells() const { return gridSetData_.cells; }
163
164     //! Returns the grid atom indices covering all grids
165     gmx::ArrayRef<const int> atomIndices() const { return gridSetData_.atomIndices; }
166
167     //! Returns whether we have perturbed non-bonded interactions
168     bool haveFep() const { return haveFep_; }
169
170     //! Returns the unit cell in \p box
171     void getBox(matrix box) const { copy_mat(box_, box); }
172
173     //! Returns the maximum number of columns across all grids
174     int numColumnsMax() const { return numColumnsMax_; }
175
176     //! Sets the maximum number of columns across all grids
177     void setNumColumnsMax(int numColumnsMax) { numColumnsMax_ = numColumnsMax; }
178
179 private:
180     /* Data members */
181     //! The domain setup
182     DomainSetup domainSetup_;
183     //! The search grids
184     std::vector<Grid> grids_;
185     //! The cell and atom index data which runs over all grids
186     GridSetData gridSetData_;
187     //! Tells whether we have perturbed non-bonded interactions
188     bool haveFep_;
189     //! The periodic unit-cell
190     matrix box_;
191     //! The number of local real atoms, i.e. without padded atoms, local atoms: 0 to numAtomsLocal_
192     int numRealAtomsLocal_;
193     //! The total number of real atoms, i.e. without padded atoms
194     int numRealAtomsTotal_;
195     //! Working data for constructing a single grid, one entry per thread
196     std::vector<GridWork> gridWork_;
197     //! Maximum number of columns across all grids
198     int numColumnsMax_;
199 };
200
201 } // namespace Nbnxm
202
203 #endif