68e07759b453f08087b585a424114815368be330
[alexxy/gromacs.git] / src / gromacs / nbnxm / internal.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,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 Declares internal nbnxm module details
39  *
40  * \author Berk Hess <hess@kth.se>
41  *
42  * \ingroup module_nbnxm
43  */
44
45 #ifndef GMX_NBNXM_INTERNAL_H
46 #define GMX_NBNXM_INTERNAL_H
47
48 #include <memory>
49 #include <vector>
50
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/math/vectypes.h"
53 #include "gromacs/nbnxm/pairlist.h"
54 #include "gromacs/simd/simd.h"
55 #include "gromacs/timing/cyclecounter.h"
56 #include "gromacs/utility/alignedallocator.h"
57 #include "gromacs/utility/arrayref.h"
58 #include "gromacs/utility/real.h"
59
60 struct gmx_domdec_zones_t;
61 struct nbnxn_grid_t;
62
63
64 // TODO Document after refactoring
65 #ifndef DOXYGEN
66
67 /* Strides for x/f with xyz and xyzq coordinate (and charge) storage */
68 #define STRIDE_XYZ         3
69 #define STRIDE_XYZQ        4
70 /* Size of packs of x, y or z with SIMD packed coords/forces */
71 static const int c_packX4 = 4;
72 static const int c_packX8 = 8;
73 /* Strides for a pack of 4 and 8 coordinates/forces */
74 #define STRIDE_P4         (DIM*c_packX4)
75 #define STRIDE_P8         (DIM*c_packX8)
76
77 /* Returns the index in a coordinate array corresponding to atom a */
78 template<int packSize> static inline int atom_to_x_index(int a)
79 {
80     return DIM*(a & ~(packSize - 1)) + (a & (packSize - 1));
81 }
82
83
84 #if GMX_SIMD
85 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
86 #define NBNXN_MEM_ALIGN  (GMX_SIMD_REAL_WIDTH*sizeof(real))
87 #else
88 /* No alignment required, but set it so we can call the same routines */
89 #define NBNXN_MEM_ALIGN  32
90 #endif
91
92 #endif // !DOXYGEN
93
94
95 /*! \brief Convenience declaration for an std::vector with aligned memory */
96 template <class T>
97 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
98
99
100 /* Local cycle count struct for profiling */
101 typedef struct {
102     int          count;
103     gmx_cycles_t c;
104     gmx_cycles_t start;
105 } nbnxn_cycle_t;
106
107 /* Local cycle count enum for profiling */
108 enum {
109     enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
110 };
111
112 /* Thread-local work struct, contains part of nbnxn_grid_t */
113 struct nbnxn_search_work_t
114 {
115     nbnxn_search_work_t();
116
117     ~nbnxn_search_work_t();
118
119     gmx_cache_protect_t       cp0;          /* Buffer to avoid cache polution */
120
121     std::vector<int>          cxy_na;       /* Grid column atom counts temporary buffer */
122
123     std::vector<int>          sortBuffer;   /* Temporary buffer for sorting atoms within a grid column */
124
125     nbnxn_buffer_flags_t      buffer_flags; /* Flags for force buffer access */
126
127     int                       ndistc;       /* Number of distance checks for flop counting */
128
129
130     std::unique_ptr<t_nblist> nbl_fep;      /* Temporary FEP list for load balancing */
131
132     nbnxn_cycle_t             cc[enbsCCnr]; /* Counters for thread-local cycles */
133
134     gmx_cache_protect_t       cp1;          /* Buffer to avoid cache polution */
135 };
136
137 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
138 struct nbnxn_search
139 {
140     /* \brief Constructor
141      *
142      * \param[in] ePBC         The periodic boundary conditions
143      * \param[in] n_dd_cells   The number of domain decomposition cells per dimension, without DD nullptr should be passed
144      * \param[in] zones        The domain decomposition zone setup, without DD nullptr should be passed
145      * \param[in] bFEP         Tells whether non-bonded interactions are perturbed
146      * \param[in] nthread_max  The maximum number of threads used in the search
147      */
148
149     nbnxn_search(int                       ePBC,
150                  const ivec               *n_dd_cells,
151                  const gmx_domdec_zones_t *zones,
152                  gmx_bool                  bFEP,
153                  int                       nthread_max);
154
155     gmx_bool                   bFEP;            /* Do we have perturbed atoms? */
156     int                        ePBC;            /* PBC type enum                              */
157     matrix                     box;             /* The periodic unit-cell                     */
158
159     gmx_bool                   DomDec;          /* Are we doing domain decomposition?         */
160     ivec                       dd_dim;          /* Are we doing DD in x,y,z?                  */
161     const gmx_domdec_zones_t  *zones;           /* The domain decomposition zones        */
162
163     std::vector<nbnxn_grid_t>  grid;            /* Array of grids, size ngrid                 */
164     std::vector<int>           cell;            /* Actual allocated cell array for all grids  */
165     std::vector<int>           a;               /* Atom index for grid, the inverse of cell   */
166
167     int                        natoms_local;    /* The local atoms run from 0 to natoms_local */
168     int                        natoms_nonlocal; /* The non-local atoms run from natoms_local
169                                                  * to natoms_nonlocal */
170
171     gmx_bool             print_cycles;
172     int                  search_count;
173     nbnxn_cycle_t        cc[enbsCCnr];
174
175     /* Thread-local work data */
176     mutable std::vector<nbnxn_search_work_t> work; /* Work array, one entry for each thread */
177 };
178
179
180 /*! \brief Start an nbnxn cycle counter */
181 static inline void nbs_cycle_start(nbnxn_cycle_t *cc)
182 {
183     cc->start = gmx_cycles_read();
184 }
185
186 /*! \brief Stop an nbnxn cycle counter */
187 static inline void nbs_cycle_stop(nbnxn_cycle_t *cc)
188 {
189     cc->c += gmx_cycles_read() - cc->start;
190     cc->count++;
191 }
192
193 #endif