Merge branch 'origin/release-2020' into master
[alexxy/gromacs.git] / src / gromacs / nbnxm / grid.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 /*! \internal \file
37  *
38  * \brief
39  * Declares the Grid class.
40  *
41  * This class provides functionality for setting up and accessing atoms
42  * on a grid for one domain decomposition zone. This grid is used for
43  * generating cluster pair lists for computing non-bonded pair interactions.
44  * The grid consists of a regular array of columns along dimensions x and y.
45  * Along z the number of cells and their boundaries vary between the columns.
46  * Each cell can hold one or more clusters of atoms, depending on the grid
47  * geometry, which is set by the pair-list type.
48  *
49  * \author Berk Hess <hess@kth.se>
50  * \ingroup module_nbnxm
51  */
52
53 #ifndef GMX_NBNXM_GRID_H
54 #define GMX_NBNXM_GRID_H
55
56 #include <memory>
57 #include <vector>
58
59 #include "gromacs/gpu_utils/hostallocator.h"
60 #include "gromacs/math/vectypes.h"
61 #include "gromacs/utility/alignedallocator.h"
62 #include "gromacs/utility/arrayref.h"
63 #include "gromacs/utility/range.h"
64
65 struct gmx_domdec_zones_t;
66 struct nbnxn_atomdata_t;
67 struct nbnxn_search;
68 enum class PairlistType;
69
70 namespace gmx
71 {
72 class UpdateGroupsCog;
73 } // namespace gmx
74
75 namespace Nbnxm
76 {
77
78 struct GridSetData;
79 struct GridWork;
80
81 /*! \internal
82  * \brief Bounding box for a nbnxm atom cluster
83  *
84  * \note Should be aligned in memory to enable 4-wide SIMD operations.
85  */
86 struct BoundingBox
87 {
88     /*! \internal
89      * \brief Corner for the bounding box, padded with one element to enable 4-wide SIMD operations
90      */
91     struct Corner
92     {
93         //! Returns a corner with the minimum coordinates along each dimension
94         static Corner min(const Corner& c1, const Corner& c2)
95         {
96             Corner cMin;
97
98             cMin.x = std::min(c1.x, c2.x);
99             cMin.y = std::min(c1.y, c2.y);
100             cMin.z = std::min(c1.z, c2.z);
101             /* This value of the padding is irrelevant, as long as it
102              * is initialized. We use min to allow auto-vectorization.
103              */
104             cMin.padding = std::min(c1.padding, c2.padding);
105
106             return cMin;
107         }
108
109         //! Returns a corner with the maximum coordinates along each dimension
110         static Corner max(const Corner& c1, const Corner& c2)
111         {
112             Corner cMax;
113
114             cMax.x       = std::max(c1.x, c2.x);
115             cMax.y       = std::max(c1.y, c2.y);
116             cMax.z       = std::max(c1.z, c2.z);
117             cMax.padding = std::max(c1.padding, c2.padding);
118
119             return cMax;
120         }
121
122         //! Returns a pointer for SIMD loading of a Corner object
123         const float* ptr() const { return &x; }
124
125         //! Returns a pointer for SIMD storing of a Corner object
126         float* ptr() { return &x; }
127
128         //! x coordinate
129         float x;
130         //! y coordinate
131         float y;
132         //! z coordinate
133         float z;
134         //! padding, unused, but should be set to avoid operations on unitialized data
135         float padding;
136     };
137
138     //! lower, along x and y and z, corner
139     Corner lower;
140     //! upper, along x and y and z, corner
141     Corner upper;
142 };
143
144 /*! \internal
145  * \brief Bounding box for one dimension of a grid cell
146  */
147 struct BoundingBox1D
148 {
149     //! lower bound
150     float lower;
151     //! upper bound
152     float upper;
153 };
154
155 } // namespace Nbnxm
156
157 namespace Nbnxm
158 {
159
160 /*! \internal
161  * \brief A pair-search grid object for one domain decomposition zone
162  *
163  * This is a rectangular 3D grid covering a potentially non-rectangular
164  * volume which is either the whole unit cell or the local zone or part
165  * of a non-local zone when using domain decomposition. The grid cells
166  * are even spaced along x/y and irregular along z. Each cell is sub-divided
167  * into atom clusters. With a CPU geometry, each cell contains 1 or 2 clusters.
168  * With a GPU geometry, each cell contains up to 8 clusters. The geometry is
169  * set by the pairlist type which is the only argument of the constructor.
170  *
171  * When multiple grids are used, i.e. with domain decomposition, we want
172  * to avoid the overhead of multiple coordinate arrays or extra indexing.
173  * Therefore each grid stores a cell offset, so a contiguous cell index
174  * can be used to index atom arrays. All methods returning atom indices
175  * return indices which index into a full atom array.
176  *
177  * Note that when atom groups, instead of individual atoms, are assigned
178  * to grid cells, individual atoms can be geometrically outside the cell
179  * and grid that they have been assigned to (as determined by the center
180  * or geometry of the atom group they belong to).
181  */
182 class Grid
183 {
184 public:
185     /*! \internal
186      * \brief The cluster and cell geometry of a grid
187      */
188     struct Geometry
189     {
190         //! Constructs the cluster/cell geometry given the type of pairlist
191         Geometry(PairlistType pairlistType);
192
193         //! Is this grid simple (CPU) or hierarchical (GPU)
194         bool isSimple;
195         //! Number of atoms per cluster
196         int numAtomsICluster;
197         //! Number of atoms for list j-clusters
198         int numAtomsJCluster;
199         //! Number of atoms per cell
200         int numAtomsPerCell;
201         //! 2log of na_c
202         int numAtomsICluster2Log;
203     };
204
205     //! The physical dimensions of a grid \internal
206     struct Dimensions
207     {
208         //! The lower corner of the (local) grid
209         rvec lowerCorner;
210         //! The upper corner of the (local) grid
211         rvec upperCorner;
212         //! The physical grid size: upperCorner - lowerCorner
213         rvec gridSize;
214         //! An estimate for the atom number density of the region targeted by the grid
215         real atomDensity;
216         //! The maximum distance an atom can be outside of a cell and outside of the grid
217         real maxAtomGroupRadius;
218         //! Size of cell along dimension x and y
219         real cellSize[DIM - 1];
220         //! 1/size of a cell along dimensions x and y
221         real invCellSize[DIM - 1];
222         //! The number of grid cells along dimensions x and y
223         int numCells[DIM - 1];
224     };
225
226     //! Constructs a grid given the type of pairlist
227     Grid(PairlistType pairlistType, const bool& haveFep);
228
229     //! Returns the geometry of the grid cells
230     const Geometry& geometry() const { return geometry_; }
231
232     //! Returns the dimensions of the grid
233     const Dimensions& dimensions() const { return dimensions_; }
234
235     //! Returns the total number of grid columns
236     int numColumns() const { return dimensions_.numCells[XX] * dimensions_.numCells[YY]; }
237
238     //! Returns the total number of grid cells
239     int numCells() const { return numCellsTotal_; }
240
241     //! Returns the cell offset of (the first cell of) this grid in the list of cells combined over all grids
242     int cellOffset() const { return cellOffset_; }
243
244     //! Returns the maximum number of grid cells in a column
245     int numCellsColumnMax() const { return numCellsColumnMax_; }
246
247     //! Returns the start of the source atom range mapped to this grid
248     int srcAtomBegin() const { return srcAtomBegin_; }
249
250     //! Returns the end of the source atom range mapped to this grid
251     int srcAtomEnd() const { return srcAtomEnd_; }
252
253     //! Returns the first cell index in the grid, starting at 0 in this grid
254     int firstCellInColumn(int columnIndex) const { return cxy_ind_[columnIndex]; }
255
256     //! Returns the number of cells in the column
257     int numCellsInColumn(int columnIndex) const
258     {
259         return cxy_ind_[columnIndex + 1LL] - cxy_ind_[columnIndex];
260     }
261
262     //! Returns the index of the first atom in the column
263     int firstAtomInColumn(int columnIndex) const
264     {
265         return (cellOffset_ + cxy_ind_[columnIndex]) * geometry_.numAtomsPerCell;
266     }
267
268     //! Returns the number of real atoms in the column
269     int numAtomsInColumn(int columnIndex) const { return cxy_na_[columnIndex]; }
270
271     /*! \brief Returns a view of the number of non-filler, atoms for each grid column
272      *
273      * \todo Needs a useful name. */
274     gmx::ArrayRef<const int> cxy_na() const { return cxy_na_; }
275     /*! \brief Returns a view of the grid-local cell index for each grid column
276      *
277      * \todo Needs a useful name. */
278     gmx::ArrayRef<const int> cxy_ind() const { return cxy_ind_; }
279
280     //! Returns the number of real atoms in the column
281     int numAtomsPerCell() const { return geometry_.numAtomsPerCell; }
282
283     //! Returns the number of atoms in the column including padding
284     int paddedNumAtomsInColumn(int columnIndex) const
285     {
286         return numCellsInColumn(columnIndex) * geometry_.numAtomsPerCell;
287     }
288
289     //! Returns the end of the atom index range on the grid, including padding
290     int atomIndexEnd() const { return (cellOffset_ + numCellsTotal_) * geometry_.numAtomsPerCell; }
291
292     //! Returns whether any atom in the cluster is perturbed
293     bool clusterIsPerturbed(int clusterIndex) const { return fep_[clusterIndex] != 0U; }
294
295     //! Returns whether the given atom in the cluster is perturbed
296     bool atomIsPerturbed(int clusterIndex, int atomIndexInCluster) const
297     {
298         return (fep_[clusterIndex] & (1 << atomIndexInCluster)) != 0U;
299     }
300
301     //! Returns the free-energy perturbation bits for the cluster
302     unsigned int fepBits(int clusterIndex) const { return fep_[clusterIndex]; }
303
304     //! Returns the i-bounding boxes for all clusters on the grid
305     gmx::ArrayRef<const BoundingBox> iBoundingBoxes() const { return bb_; }
306
307     //! Returns the j-bounding boxes for all clusters on the grid
308     gmx::ArrayRef<const BoundingBox> jBoundingBoxes() const { return bbj_; }
309
310     //! Returns the packed bounding boxes for all clusters on the grid, empty with a CPU list
311     gmx::ArrayRef<const float> packedBoundingBoxes() const { return pbb_; }
312
313     //! Returns the bounding boxes along z for all cells on the grid
314     gmx::ArrayRef<const BoundingBox1D> zBoundingBoxes() const { return bbcz_; }
315
316     //! Returns the flags for all clusters on the grid
317     gmx::ArrayRef<const int> clusterFlags() const { return flags_; }
318
319     //! Returns the number of clusters for all cells on the grid, empty with a CPU geometry
320     gmx::ArrayRef<const int> numClustersPerCell() const { return numClusters_; }
321
322     //! Returns the cluster index for an atom
323     int atomToCluster(int atomIndex) const { return (atomIndex >> geometry_.numAtomsICluster2Log); }
324
325     //! Returns the total number of clusters on the grid
326     int numClusters() const
327     {
328         if (geometry_.isSimple)
329         {
330             return numCellsTotal_;
331         }
332         else
333         {
334             return numClustersTotal_;
335         }
336     }
337
338     //! Sets the grid dimensions
339     void setDimensions(int                ddZone,
340                        int                numAtoms,
341                        gmx::RVec          lowerCorner,
342                        gmx::RVec          upperCorner,
343                        real               atomDensity,
344                        real               maxAtomGroupRadius,
345                        bool               haveFep,
346                        gmx::PinningPolicy pinningPolicy);
347
348     //! Sets the cell indices using indices in \p gridSetData and \p gridWork
349     void setCellIndices(int                            ddZone,
350                         int                            cellOffset,
351                         GridSetData*                   gridSetData,
352                         gmx::ArrayRef<GridWork>        gridWork,
353                         gmx::Range<int>                atomRange,
354                         const int*                     atinfo,
355                         gmx::ArrayRef<const gmx::RVec> x,
356                         int                            numAtomsMoved,
357                         nbnxn_atomdata_t*              nbat);
358
359     //! Determine in which grid columns atoms should go, store cells and atom counts in \p cell and \p cxy_na
360     static void calcColumnIndices(const Grid::Dimensions&        gridDims,
361                                   const gmx::UpdateGroupsCog*    updateGroupsCog,
362                                   gmx::Range<int>                atomRange,
363                                   gmx::ArrayRef<const gmx::RVec> x,
364                                   int                            dd_zone,
365                                   const int*                     move,
366                                   int                            thread,
367                                   int                            nthread,
368                                   gmx::ArrayRef<int>             cell,
369                                   gmx::ArrayRef<int>             cxy_na);
370
371 private:
372     /*! \brief Fill a pair search cell with atoms
373      *
374      * Potentially sorts atoms and sets the interaction flags.
375      */
376     void fillCell(GridSetData*                   gridSetData,
377                   nbnxn_atomdata_t*              nbat,
378                   int                            atomStart,
379                   int                            atomEnd,
380                   const int*                     atinfo,
381                   gmx::ArrayRef<const gmx::RVec> x,
382                   BoundingBox gmx_unused* bb_work_aligned);
383
384     //! Spatially sort the atoms within the given column range, for CPU geometry
385     void sortColumnsCpuGeometry(GridSetData*                   gridSetData,
386                                 int                            dd_zone,
387                                 const int*                     atinfo,
388                                 gmx::ArrayRef<const gmx::RVec> x,
389                                 nbnxn_atomdata_t*              nbat,
390                                 gmx::Range<int>                columnRange,
391                                 gmx::ArrayRef<int>             sort_work);
392
393     //! Spatially sort the atoms within the given column range, for GPU geometry
394     void sortColumnsGpuGeometry(GridSetData*                   gridSetData,
395                                 int                            dd_zone,
396                                 const int*                     atinfo,
397                                 gmx::ArrayRef<const gmx::RVec> x,
398                                 nbnxn_atomdata_t*              nbat,
399                                 gmx::Range<int>                columnRange,
400                                 gmx::ArrayRef<int>             sort_work);
401
402     /* Data members */
403     //! The geometry of the grid clusters and cells
404     Geometry geometry_;
405     //! The physical dimensions of the grid
406     Dimensions dimensions_;
407
408     //! The total number of cells in this grid
409     int numCellsTotal_;
410     //! Index in nbs->cell corresponding to cell 0
411     int cellOffset_;
412     //! The maximum number of cells in a column
413     int numCellsColumnMax_;
414
415     //! The start of the source atom range mapped to this grid
416     int srcAtomBegin_;
417     //! The end of the source atom range mapped to this grid
418     int srcAtomEnd_;
419
420     /* Grid data */
421     /*! \brief The number of, non-filler, atoms for each grid column.
422      *
423      * \todo Needs a useful name. */
424     gmx::HostVector<int> cxy_na_;
425     /*! \brief The grid-local cell index for each grid column
426      *
427      * \todo Needs a useful name. */
428     gmx::HostVector<int> cxy_ind_;
429
430     //! The number of cluster for each cell
431     std::vector<int> numClusters_;
432
433     /* Bounding boxes */
434     //! Bounding boxes in z for the cells
435     std::vector<BoundingBox1D> bbcz_;
436     //! 3D bounding boxes for the sub cells
437     std::vector<BoundingBox, gmx::AlignedAllocator<BoundingBox>> bb_;
438     //! 3D j-bounding boxes for the case where the i- and j-cluster sizes are different
439     std::vector<BoundingBox, gmx::AlignedAllocator<BoundingBox>> bbjStorage_;
440     //! 3D j-bounding boxes
441     gmx::ArrayRef<BoundingBox> bbj_;
442     //! 3D bounding boxes in packed xxxx format per cell
443     std::vector<float, gmx::AlignedAllocator<float>> pbb_;
444
445     //! Tells whether we have perturbed interactions, authorative source is in GridSet (never modified)
446     const bool& haveFep_;
447
448     /* Bit-flag information */
449     //! Flags for properties of clusters in each cell
450     std::vector<int> flags_;
451     //! Signal bits for atoms in each cell that tell whether an atom is perturbed
452     std::vector<unsigned int> fep_;
453
454     /* Statistics */
455     //! Total number of clusters, used for printing
456     int numClustersTotal_;
457 };
458
459 } // namespace Nbnxm
460
461 #endif