Update instructions in containers.rst
[gromacs.git] / src / gromacs / domdec / domdec_struct.h
blobdb67680fb1465e67889d9f6ce471a3dba34668af
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,2020, 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>
78 class HashedMap;
79 class LocalAtomSetManager;
80 class GpuHaloExchange;
81 } // namespace gmx
83 /*! \internal
84 * \brief Pair interaction zone and atom range for an i-zone
86 struct DDPairInteractionRanges
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 gmx::IVec shift0 = { 0, 0, 0 };
98 //! Maximum shifts to consider
99 gmx::IVec shift1 = { 0, 0, 0 };
102 typedef struct gmx_domdec_zone_size
104 /* Zone lower corner in triclinic coordinates */
105 gmx::RVec x0 = { 0, 0, 0 };
106 /* Zone upper corner in triclinic coordinates */
107 gmx::RVec x1 = { 0, 0, 0 };
108 /* Zone bounding box lower corner in Cartesian coords */
109 gmx::RVec bb_x0 = { 0, 0, 0 };
110 /* Zone bounding box upper corner in Cartesian coords */
111 gmx::RVec bb_x1 = { 0, 0, 0 };
112 } gmx_domdec_zone_size_t;
114 struct gmx_domdec_zones_t
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;
130 struct gmx_ddbox_t
132 int npbcdim;
133 int nboundeddim;
134 gmx::RVec box0 = { 0, 0, 0 };
135 gmx::RVec box_size = { 0, 0, 0 };
136 /* Tells if the box is skewed for each of the three cartesian directions */
137 gmx::IVec tric_dir = { 0, 0, 0 };
138 gmx::RVec skew_fac = { 0, 0, 0 };
139 /* Orthogonal vectors for triclinic cells, Cartesian index */
140 rvec v[DIM][DIM];
141 /* Normal vectors for the cells walls */
142 rvec normal[DIM];
145 /*! \internal \brief Provides information about properties of the unit cell */
146 struct UnitCellInfo
148 //! Constructor
149 UnitCellInfo(const t_inputrec& ir);
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;
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);
166 /* The DD particle-particle nodes only */
167 /* The communication setup within the communicator all
168 * defined in dd->comm in domdec.c
170 int nnodes = 0;
171 MPI_Comm mpi_comm_all = MPI_COMM_NULL;
172 /* The local DD cell index and rank */
173 gmx::IVec ci = { 0, 0, 0 };
174 int rank = 0;
175 gmx::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];
184 /* Properties of the unit cell */
185 UnitCellInfo unitCellInfo;
187 /* The communication setup, identical for each cell, cartesian index */
188 //! Todo: refactor nbnxm to not rely on this sometimes being a nullptr so this can be IVec
189 ivec numCells = { 0, 0, 0 };
190 int ndim = 0;
191 gmx::IVec dim = { 0, 0, 0 }; /* indexed by 0 to ndim */
193 /* Forward and backward neighboring cells, indexed by 0 to ndim */
194 int neighbor[DIM][2] = { { 0, 0 }, { 0, 0 }, { 0, 0 } };
196 /* Only available on the master node */
197 std::unique_ptr<AtomDistribution> ma;
199 /* Global atom number to interaction list */
200 gmx_reverse_top_t* reverse_top = nullptr;
201 int nbonded_global = 0;
202 int nbonded_local = 0;
204 /* Whether we have non-self exclusion */
205 bool haveExclusions = false;
207 /* Vsite stuff */
208 gmx::HashedMap<int>* ga2la_vsite = nullptr;
209 gmx_domdec_specat_comm_t* vsite_comm = nullptr;
210 std::vector<int> vsite_requestedGlobalAtomIndices;
212 /* Constraint stuff */
213 gmx_domdec_constraints_t* constraints = nullptr;
214 gmx_domdec_specat_comm_t* constraint_comm = nullptr;
216 /* The number of home atom groups */
217 int ncg_home = 0;
218 /* Global atom group indices for the home and all non-home groups */
219 std::vector<int> globalAtomGroupIndices;
221 /* Index from the local atoms to the global atoms, covers home and received zones */
222 std::vector<int> globalAtomIndices;
224 /* Global atom number to local atom number list */
225 gmx_ga2la_t* ga2la = nullptr;
227 /* Communication stuff */
228 gmx_domdec_comm_t* comm = nullptr;
230 /* The partioning count, to keep track of the state */
231 int64_t ddp_count = 0;
233 /* The managed atom sets that are updated in domain decomposition */
234 gmx::LocalAtomSetManager* atomSets = nullptr;
236 /* gmx_pme_recv_f buffer */
237 std::vector<gmx::RVec> pmeForceReceiveBuffer;
239 /* GPU halo exchange objects: this structure supports a vector of pulses for each dimension */
240 std::vector<std::unique_ptr<gmx::GpuHaloExchange>> gpuHaloExchange[DIM];
243 //! Are we the master node for domain decomposition
244 static inline bool DDMASTER(const gmx_domdec_t& dd)
246 return dd.rank == dd.masterrank;
249 //! Are we the master node for domain decomposition, deprecated
250 static inline bool DDMASTER(const gmx_domdec_t* dd)
252 return dd->rank == dd->masterrank;
255 #endif