55a6a3e9690b253df20089df981867b88cbb3808
[alexxy/gromacs.git] / src / gromacs / domdec / domdec_struct.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /*! \libinternal \file
38  * \brief Declares structures related to domain decomposition.
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42  * \inlibraryapi
43  * \ingroup module_domdec
44  */
45 #ifndef GMX_DOMDEC_DOMDEC_STRUCT_H
46 #define GMX_DOMDEC_DOMDEC_STRUCT_H
47
48 #include <cstddef>
49
50 #include <array>
51 #include <memory>
52 #include <vector>
53
54 #include "gromacs/math/vectypes.h"
55 #include "gromacs/topology/block.h"
56 #include "gromacs/utility/basedefinitions.h"
57 #include "gromacs/utility/gmxmpi.h"
58 #include "gromacs/utility/range.h"
59 #include "gromacs/utility/real.h"
60
61 //! Max number of zones in domain decomposition
62 #define DD_MAXZONE 8
63 //! Max number of izones in domain decomposition
64 #define DD_MAXIZONE 4
65
66 struct AtomDistribution;
67 struct gmx_domdec_comm_t;
68 struct gmx_domdec_constraints_t;
69 struct gmx_domdec_specat_comm_t;
70 class gmx_ga2la_t;
71 struct gmx_pme_comm_n_box_t;
72 struct gmx_reverse_top_t;
73 struct t_inputrec;
74
75 namespace gmx
76 {
77 template<typename T>
78 class HashedMap;
79 class LocalAtomSetManager;
80 class GpuHaloExchange;
81 } // namespace gmx
82
83 /*! \internal
84  * \brief Pair interaction zone and atom range for an i-zone
85  */
86 struct DDPairInteractionRanges
87 {
88     //! The index of this i-zone in the i-zone list
89     int iZoneIndex = -1;
90     //! The range of j-zones
91     gmx::Range<int> jZoneRange;
92     //! The i-atom range
93     gmx::Range<int> iAtomRange;
94     //! The j-atom range
95     gmx::Range<int> jAtomRange;
96     //! Minimum shifts to consider
97     ivec shift0 = {};
98     //! Maximum shifts to consider
99     ivec shift1 = {};
100 };
101
102 typedef struct
103 {
104     /* Zone lower corner in triclinic coordinates         */
105     rvec x0 = {};
106     /* Zone upper corner in triclinic coordinates         */
107     rvec x1 = {};
108     /* Zone bounding box lower corner in Cartesian coords */
109     rvec bb_x0 = {};
110     /* Zone bounding box upper corner in Cartesian coords */
111     rvec bb_x1 = {};
112 } gmx_domdec_zone_size_t;
113
114 struct gmx_domdec_zones_t
115 {
116     /* The number of zones including the home zone */
117     int n = 0;
118     /* The shift of the zones with respect to the home zone */
119     ivec shift[DD_MAXZONE] = {};
120     /* The charge group boundaries for the zones */
121     int cg_range[DD_MAXZONE + 1] = {};
122     /* The pair interaction zone and atom ranges per each i-zone */
123     std::vector<DDPairInteractionRanges> iZones;
124     /* Boundaries of the zones */
125     gmx_domdec_zone_size_t size[DD_MAXZONE];
126     /* The cg density of the home zone */
127     real dens_zone0 = 0;
128 };
129
130 struct gmx_ddbox_t
131 {
132     int  npbcdim;
133     int  nboundeddim;
134     rvec box0;
135     rvec box_size;
136     /* Tells if the box is skewed for each of the three cartesian directions */
137     ivec tric_dir;
138     rvec skew_fac;
139     /* Orthogonal vectors for triclinic cells, Cartesian index */
140     rvec v[DIM][DIM];
141     /* Normal vectors for the cells walls */
142     rvec normal[DIM];
143 };
144
145 /*! \internal \brief Provides information about properties of the unit cell */
146 struct UnitCellInfo
147 {
148     //! Constructor
149     UnitCellInfo(const t_inputrec& ir);
150
151     //! We have PBC from dim 0 (X) up to npbcdim
152     int npbcdim;
153     //! The system is bounded from 0 (X) to numBoundedDimensions
154     int numBoundedDimensions;
155     //! Tells whether the box bounding the atoms is dynamic
156     bool ddBoxIsDynamic;
157     //! Screw PBC?
158     bool haveScrewPBC;
159 };
160
161 struct gmx_domdec_t
162 { //NOLINT(clang-analyzer-optin.performance.Padding)
163     //! Constructor, only partial for now
164     gmx_domdec_t(const t_inputrec& ir);
165
166     /* The DD particle-particle nodes only */
167     /* The communication setup within the communicator all
168      * defined in dd->comm in domdec.c
169      */
170     int      nnodes       = 0;
171     MPI_Comm mpi_comm_all = MPI_COMM_NULL;
172     /* The local DD cell index and rank */
173     ivec ci         = { 0, 0, 0 };
174     int  rank       = 0;
175     ivec master_ci  = { 0, 0, 0 };
176     int  masterrank = 0;
177     /* Communication with the PME only nodes */
178     int                   pme_nodeid           = 0;
179     gmx_bool              pme_receive_vir_ener = false;
180     gmx_pme_comm_n_box_t* cnb                  = nullptr;
181     int                   nreq_pme             = 0;
182     MPI_Request           req_pme[8];
183
184     /* Properties of the unit cell */
185     UnitCellInfo unitCellInfo;
186
187     /* The communication setup, identical for each cell, cartesian index */
188     ivec nc   = { 0, 0, 0 };
189     int  ndim = 0;
190     ivec dim  = { 0, 0, 0 }; /* indexed by 0 to ndim */
191
192     /* Forward and backward neighboring cells, indexed by 0 to ndim */
193     int neighbor[DIM][2] = { { 0, 0 }, { 0, 0 }, { 0, 0 } };
194
195     /* Only available on the master node */
196     std::unique_ptr<AtomDistribution> ma;
197
198     /* Global atom number to interaction list */
199     gmx_reverse_top_t* reverse_top    = nullptr;
200     int                nbonded_global = 0;
201     int                nbonded_local  = 0;
202
203     /* Whether we have non-self exclusion */
204     bool haveExclusions = false;
205
206     /* Vsite stuff */
207     gmx::HashedMap<int>*      ga2la_vsite = nullptr;
208     gmx_domdec_specat_comm_t* vsite_comm  = nullptr;
209     std::vector<int>          vsite_requestedGlobalAtomIndices;
210
211     /* Constraint stuff */
212     gmx_domdec_constraints_t* constraints     = nullptr;
213     gmx_domdec_specat_comm_t* constraint_comm = nullptr;
214
215     /* The number of home atom groups */
216     int ncg_home = 0;
217     /* Global atom group indices for the home and all non-home groups */
218     std::vector<int> globalAtomGroupIndices;
219
220     /* Index from the local atoms to the global atoms, covers home and received zones */
221     std::vector<int> globalAtomIndices;
222
223     /* Global atom number to local atom number list */
224     gmx_ga2la_t* ga2la = nullptr;
225
226     /* Communication stuff */
227     gmx_domdec_comm_t* comm = nullptr;
228
229     /* The partioning count, to keep track of the state */
230     int64_t ddp_count = 0;
231
232     /* The managed atom sets that are updated in domain decomposition */
233     gmx::LocalAtomSetManager* atomSets = nullptr;
234
235     /* gmx_pme_recv_f buffer */
236     std::vector<gmx::RVec> pmeForceReceiveBuffer;
237
238     /* GPU halo exchange object */
239     std::unique_ptr<gmx::GpuHaloExchange> gpuHaloExchange;
240 };
241
242 //! Are we the master node for domain decomposition
243 static inline bool DDMASTER(const gmx_domdec_t& dd)
244 {
245     return dd.rank == dd.masterrank;
246 };
247
248 //! Are we the master node for domain decomposition, deprecated
249 static inline bool DDMASTER(const gmx_domdec_t* dd)
250 {
251     return dd->rank == dd->masterrank;
252 };
253
254 #endif