Move responsibility for checking bondeds in reverse topology
[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,2020,2021, 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/gmxmpi.h"
57 #include "gromacs/utility/range.h"
58 #include "gromacs/utility/real.h"
59
60 //! Max number of zones in domain decomposition
61 #define DD_MAXZONE 8
62 //! Max number of izones in domain decomposition
63 #define DD_MAXIZONE 4
64
65 struct AtomDistribution;
66 struct gmx_domdec_comm_t;
67 struct gmx_domdec_constraints_t;
68 struct gmx_domdec_specat_comm_t;
69 class gmx_ga2la_t;
70 struct gmx_pme_comm_n_box_t;
71 struct t_inputrec;
72 struct gmx_reverse_top_t;
73 struct gmx_mtop_t;
74 struct ReverseTopOptions;
75
76 namespace gmx
77 {
78 template<typename T>
79 class HashedMap;
80 class LocalAtomSetManager;
81 class GpuHaloExchange;
82 } // namespace gmx
83
84 /*! \internal
85  * \brief Pair interaction zone and atom range for an i-zone
86  */
87 struct DDPairInteractionRanges
88 {
89     //! The index of this i-zone in the i-zone list
90     int iZoneIndex = -1;
91     //! The range of j-zones
92     gmx::Range<int> jZoneRange;
93     //! The i-atom range
94     gmx::Range<int> iAtomRange;
95     //! The j-atom range
96     gmx::Range<int> jAtomRange;
97     //! Minimum shifts to consider
98     gmx::IVec shift0 = { 0, 0, 0 };
99     //! Maximum shifts to consider
100     gmx::IVec shift1 = { 0, 0, 0 };
101 };
102
103 typedef struct gmx_domdec_zone_size
104 {
105     /* Zone lower corner in triclinic coordinates         */
106     gmx::RVec x0 = { 0, 0, 0 };
107     /* Zone upper corner in triclinic coordinates         */
108     gmx::RVec x1 = { 0, 0, 0 };
109     /* Zone bounding box lower corner in Cartesian coords */
110     gmx::RVec bb_x0 = { 0, 0, 0 };
111     /* Zone bounding box upper corner in Cartesian coords */
112     gmx::RVec bb_x1 = { 0, 0, 0 };
113 } gmx_domdec_zone_size_t;
114
115 struct gmx_domdec_zones_t
116 {
117     /* The number of zones including the home zone */
118     int n = 0;
119     /* The shift of the zones with respect to the home zone */
120     std::array<ivec, DD_MAXZONE> shift;
121     /* The charge group boundaries for the zones */
122     std::array<int, DD_MAXZONE + 1> cg_range;
123     /* The pair interaction zone and atom ranges per each i-zone */
124     std::vector<DDPairInteractionRanges> iZones;
125     /* Boundaries of the zones */
126     std::array<gmx_domdec_zone_size_t, DD_MAXZONE> size;
127     /* The cg density of the home zone */
128     real dens_zone0 = 0;
129 };
130
131 struct gmx_ddbox_t
132 {
133     int       npbcdim;
134     int       nboundeddim;
135     gmx::RVec box0     = { 0, 0, 0 };
136     gmx::RVec box_size = { 0, 0, 0 };
137     /* Tells if the box is skewed for each of the three cartesian directions */
138     gmx::IVec tric_dir = { 0, 0, 0 };
139     gmx::RVec skew_fac = { 0, 0, 0 };
140     /* Orthogonal vectors for triclinic cells, Cartesian index */
141     rvec v[DIM][DIM];
142     /* Normal vectors for the cells walls */
143     rvec normal[DIM];
144 };
145
146 /*! \internal \brief Provides information about properties of the unit cell */
147 struct UnitCellInfo
148 {
149     //! Constructor
150     UnitCellInfo(const t_inputrec& ir);
151
152     //! We have PBC from dim 0 (X) up to npbcdim
153     int npbcdim;
154     //! The system is bounded from 0 (X) to numBoundedDimensions
155     int numBoundedDimensions;
156     //! Tells whether the box bounding the atoms is dynamic
157     bool ddBoxIsDynamic;
158     //! Screw PBC?
159     bool haveScrewPBC;
160 };
161
162 struct gmx_domdec_t
163 { //NOLINT(clang-analyzer-optin.performance.Padding)
164     //! Constructor, only partial for now
165     gmx_domdec_t(const t_inputrec& ir);
166     ~gmx_domdec_t();
167
168     /* The DD particle-particle nodes only */
169     /* The communication setup within the communicator all
170      * defined in dd->comm in domdec.c
171      */
172     int      nnodes       = 0;
173     MPI_Comm mpi_comm_all = MPI_COMM_NULL;
174     /* The local DD cell index and rank */
175     gmx::IVec ci         = { 0, 0, 0 };
176     int       rank       = 0;
177     gmx::IVec master_ci  = { 0, 0, 0 };
178     int       masterrank = 0;
179     /* Communication with the PME only nodes */
180     int                   pme_nodeid           = 0;
181     gmx_bool              pme_receive_vir_ener = false;
182     gmx_pme_comm_n_box_t* cnb                  = nullptr;
183     int                   nreq_pme             = 0;
184     MPI_Request           req_pme[8];
185
186     /* Properties of the unit cell */
187     UnitCellInfo unitCellInfo;
188
189     /* The communication setup, identical for each cell, cartesian index */
190     //! Todo: refactor nbnxm to not rely on this sometimes being a nullptr so this can be IVec
191     ivec      numCells = { 0, 0, 0 };
192     int       ndim     = 0;
193     gmx::IVec dim      = { 0, 0, 0 }; /* indexed by 0 to ndim */
194
195     /* Forward and backward neighboring cells, indexed by 0 to ndim */
196     int neighbor[DIM][2] = { { 0, 0 }, { 0, 0 }, { 0, 0 } };
197
198     /* Only available on the master node */
199     std::unique_ptr<AtomDistribution> ma;
200
201     /* Global atom number to interaction list */
202     std::unique_ptr<gmx_reverse_top_t> reverse_top;
203
204     /* Whether we have non-self exclusion */
205     bool haveExclusions = false;
206
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;
211
212     /* Constraint stuff */
213     gmx_domdec_constraints_t* constraints     = nullptr;
214     gmx_domdec_specat_comm_t* constraint_comm = nullptr;
215
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;
220
221     /* Index from the local atoms to the global atoms, covers home and received zones */
222     std::vector<int> globalAtomIndices;
223
224     /* Global atom number to local atom number list */
225     gmx_ga2la_t* ga2la = nullptr;
226
227     /* Communication stuff */
228     gmx_domdec_comm_t* comm = nullptr;
229
230     /* The partioning count, to keep track of the state */
231     int64_t ddp_count = 0;
232
233     /* The managed atom sets that are updated in domain decomposition */
234     gmx::LocalAtomSetManager* atomSets = nullptr;
235
236     /* gmx_pme_recv_f buffer */
237     std::vector<gmx::RVec> pmeForceReceiveBuffer;
238
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];
241 };
242
243 //! Are we the master node for domain decomposition
244 static inline bool DDMASTER(const gmx_domdec_t& dd)
245 {
246     return dd.rank == dd.masterrank;
247 };
248
249 //! Are we the master node for domain decomposition, deprecated
250 static inline bool DDMASTER(const gmx_domdec_t* dd)
251 {
252     return dd->rank == dd->masterrank;
253 };
254
255 #endif