Use gmx::Range in domdec
[gromacs.git] / src / gromacs / domdec / domdec_struct.h
blobc2be7bc66ecbd1df384ac5cae63978f40beb06e1
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 /*! \libinternal \file
38 * \brief Declares structures related to domain decomposition.
40 * \author Berk Hess <hess@kth.se>
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \inlibraryapi
43 * \ingroup module_domdec
45 #ifndef GMX_DOMDEC_DOMDEC_STRUCT_H
46 #define GMX_DOMDEC_DOMDEC_STRUCT_H
48 #include <cstddef>
50 #include <array>
51 #include <memory>
52 #include <vector>
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"
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
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;
75 namespace gmx
77 template <typename T> class HashedMap;
78 class LocalAtomSetManager;
79 class GpuHaloExchange;
82 /*! \internal
83 * \brief Pair interaction zone and atom range for an i-zone
85 struct DDPairInteractionRanges
87 //! The index of this i-zone in the i-zone list
88 int iZoneIndex = -1;
89 //! The range of j-zones
90 gmx::Range<int> jZoneRange;
91 //! The i-atom range
92 gmx::Range<int> iAtomRange;
93 //! The j-atom range
94 gmx::Range<int> jAtomRange;
95 //! Minimum shifts to consider
96 ivec shift0 = { };
97 //! Maximum shifts to consider
98 ivec shift1 = { };
101 typedef struct {
102 /* Zone lower corner in triclinic coordinates */
103 rvec x0 = { };
104 /* Zone upper corner in triclinic coordinates */
105 rvec x1 = { };
106 /* Zone bounding box lower corner in Cartesian coords */
107 rvec bb_x0 = { };
108 /* Zone bounding box upper corner in Cartesian coords */
109 rvec bb_x1 = { };
110 } gmx_domdec_zone_size_t;
112 struct gmx_domdec_zones_t {
113 /* The number of zones including the home zone */
114 int n = 0;
115 /* The shift of the zones with respect to the home zone */
116 ivec shift[DD_MAXZONE] = { };
117 /* The charge group boundaries for the zones */
118 int cg_range[DD_MAXZONE+1] = { };
119 /* The pair interaction zone and atom ranges per each i-zone */
120 std::vector<DDPairInteractionRanges> iZones;
121 /* Boundaries of the zones */
122 gmx_domdec_zone_size_t size[DD_MAXZONE];
123 /* The cg density of the home zone */
124 real dens_zone0 = 0;
127 struct gmx_ddbox_t {
128 int npbcdim;
129 int nboundeddim;
130 rvec box0;
131 rvec box_size;
132 /* Tells if the box is skewed for each of the three cartesian directions */
133 ivec tric_dir;
134 rvec skew_fac;
135 /* Orthogonal vectors for triclinic cells, Cartesian index */
136 rvec v[DIM][DIM];
137 /* Normal vectors for the cells walls */
138 rvec normal[DIM];
141 /*! \internal \brief Provides information about properties of the unit cell */
142 struct UnitCellInfo
144 //! Constructor
145 UnitCellInfo(const t_inputrec &ir);
147 //! We have PBC from dim 0 (X) up to npbcdim
148 int npbcdim;
149 //! The system is bounded from 0 (X) to numBoundedDimensions
150 int numBoundedDimensions;
151 //! Tells whether the box bounding the atoms is dynamic
152 bool ddBoxIsDynamic;
153 //! Screw PBC?
154 bool haveScrewPBC;
157 struct gmx_domdec_t { //NOLINT(clang-analyzer-optin.performance.Padding)
158 //! Constructor, only partial for now
159 gmx_domdec_t(const t_inputrec &ir);
161 /* The DD particle-particle nodes only */
162 /* The communication setup within the communicator all
163 * defined in dd->comm in domdec.c
165 int nnodes = 0;
166 MPI_Comm mpi_comm_all = MPI_COMM_NULL;
167 /* The local DD cell index and rank */
168 ivec ci = { 0, 0, 0 };
169 int rank = 0;
170 ivec master_ci = { 0, 0, 0 };
171 int masterrank = 0;
172 /* Communication with the PME only nodes */
173 int pme_nodeid = 0;
174 gmx_bool pme_receive_vir_ener = false;
175 gmx_pme_comm_n_box_t *cnb = nullptr;
176 int nreq_pme = 0;
177 MPI_Request req_pme[8];
179 /* Properties of the unit cell */
180 UnitCellInfo unitCellInfo;
182 /* The communication setup, identical for each cell, cartesian index */
183 ivec nc = { 0, 0, 0 };
184 int ndim = 0;
185 ivec dim = { 0, 0, 0 }; /* indexed by 0 to ndim */
187 /* Forward and backward neighboring cells, indexed by 0 to ndim */
188 int neighbor[DIM][2] = { { 0, 0 }, { 0, 0 }, { 0, 0 } };
190 /* Only available on the master node */
191 std::unique_ptr<AtomDistribution> ma;
193 /* Global atom number to interaction list */
194 gmx_reverse_top_t *reverse_top = nullptr;
195 int nbonded_global = 0;
196 int nbonded_local = 0;
198 /* Whether we have non-self exclusion */
199 bool haveExclusions = false;
201 /* Vsite stuff */
202 gmx::HashedMap<int> *ga2la_vsite = nullptr;
203 gmx_domdec_specat_comm_t *vsite_comm = nullptr;
204 std::vector<int> vsite_requestedGlobalAtomIndices;
206 /* Constraint stuff */
207 gmx_domdec_constraints_t *constraints = nullptr;
208 gmx_domdec_specat_comm_t *constraint_comm = nullptr;
210 /* The number of home atom groups */
211 int ncg_home = 0;
212 /* Global atom group indices for the home and all non-home groups */
213 std::vector<int> globalAtomGroupIndices;
215 /* Index from the local atoms to the global atoms, covers home and received zones */
216 std::vector<int> globalAtomIndices;
218 /* Global atom number to local atom number list */
219 gmx_ga2la_t *ga2la = nullptr;
221 /* Communication stuff */
222 gmx_domdec_comm_t *comm = nullptr;
224 /* The partioning count, to keep track of the state */
225 int64_t ddp_count = 0;
227 /* The managed atom sets that are updated in domain decomposition */
228 gmx::LocalAtomSetManager * atomSets = nullptr;
230 /* gmx_pme_recv_f buffer */
231 std::vector<gmx::RVec> pmeForceReceiveBuffer;
233 /* GPU halo exchange object */
234 std::unique_ptr<gmx::GpuHaloExchange> gpuHaloExchange;
237 //! Are we the master node for domain decomposition
238 static inline bool DDMASTER(const gmx_domdec_t &dd)
240 return dd.rank == dd.masterrank;
243 //! Are we the master node for domain decomposition, deprecated
244 static inline bool DDMASTER(const gmx_domdec_t *dd)
246 return dd->rank == dd->masterrank;
249 #endif