Add nbnxm GridSet class
[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/atomdata.h"
54 #include "gromacs/nbnxm/pairlist.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 #include "gridset.h"
61
62 struct gmx_domdec_zones_t;
63
64
65 /*! \brief Convenience declaration for an std::vector with aligned memory */
66 template <class T>
67 using AlignedVector = std::vector < T, gmx::AlignedAllocator < T>>;
68
69
70 /* Local cycle count struct for profiling */
71 typedef struct {
72     int          count;
73     gmx_cycles_t c;
74     gmx_cycles_t start;
75 } nbnxn_cycle_t;
76
77 /* Local cycle count enum for profiling */
78 enum {
79     enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
80 };
81
82 /* Thread-local work struct, contains working data for Grid */
83 struct nbnxn_search_work_t
84 {
85     nbnxn_search_work_t();
86
87     ~nbnxn_search_work_t();
88
89     gmx_cache_protect_t       cp0;          /* Buffer to avoid cache polution */
90
91     std::vector<int>          sortBuffer;   /* Temporary buffer for sorting atoms within a grid column */
92
93     nbnxn_buffer_flags_t      buffer_flags; /* Flags for force buffer access */
94
95     int                       ndistc;       /* Number of distance checks for flop counting */
96
97
98     std::unique_ptr<t_nblist> nbl_fep;      /* Temporary FEP list for load balancing */
99
100     nbnxn_cycle_t             cc[enbsCCnr]; /* Counters for thread-local cycles */
101
102     gmx_cache_protect_t       cp1;          /* Buffer to avoid cache polution */
103 };
104
105 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
106 struct nbnxn_search
107 {
108     /* \brief Constructor
109      *
110      * \param[in] ePBC            The periodic boundary conditions
111      * \param[in] n_dd_cells      The number of domain decomposition cells per dimension, without DD nullptr should be passed
112      * \param[in] zones           The domain decomposition zone setup, without DD nullptr should be passed
113      * \param[in] haveFep         Tells whether non-bonded interactions are perturbed
114      * \param[in] maxNumThreads   The maximum number of threads used in the search
115      */
116     nbnxn_search(int                       ePBC,
117                  const ivec               *n_dd_cells,
118                  const gmx_domdec_zones_t *zones,
119                  PairlistType              pairlistType,
120                  bool                      haveFep,
121                  int                       maxNumthreads);
122
123     const Nbnxm::GridSet &gridSet() const
124     {
125         return gridSet_;
126     }
127
128     int                        ePBC;            /* PBC type enum                              */
129     gmx_bool                   DomDec;          /* Are we doing domain decomposition?         */
130     ivec                       dd_dim;          /* Are we doing DD in x,y,z?                  */
131     const gmx_domdec_zones_t  *zones;           /* The domain decomposition zones        */
132
133     //! The local and non-local grids
134     Nbnxm::GridSet       gridSet_;
135
136     gmx_bool             print_cycles;
137     int                  search_count;
138     nbnxn_cycle_t        cc[enbsCCnr];
139
140     /* Thread-local work data */
141     mutable std::vector<nbnxn_search_work_t> work; /* Work array, one entry for each thread */
142 };
143
144
145 /*! \brief Start an nbnxn cycle counter */
146 static inline void nbs_cycle_start(nbnxn_cycle_t *cc)
147 {
148     cc->start = gmx_cycles_read();
149 }
150
151 /*! \brief Stop an nbnxn cycle counter */
152 static inline void nbs_cycle_stop(nbnxn_cycle_t *cc)
153 {
154     cc->c += gmx_cycles_read() - cc->start;
155     cc->count++;
156 }
157
158 #endif