Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / domdec / domdec_struct.h
blob4121a84eb5abf767a2bf185d85d54672cbdee2c3
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 <memory>
51 #include <vector>
53 #include "gromacs/math/vectypes.h"
54 #include "gromacs/topology/block.h"
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/gmxmpi.h"
57 #include "gromacs/utility/real.h"
59 //! Max number of zones in domain decomposition
60 #define DD_MAXZONE 8
61 //! Max number of izones in domain decomposition
62 #define DD_MAXIZONE 4
64 struct AtomDistribution;
65 struct gmx_domdec_comm_t;
66 struct gmx_domdec_constraints_t;
67 struct gmx_domdec_specat_comm_t;
68 class gmx_ga2la_t;
69 struct gmx_pme_comm_n_box_t;
70 struct gmx_reverse_top_t;
72 namespace gmx
74 template <typename T> class HashedMap;
75 class LocalAtomSetManager;
78 typedef struct {
79 /* j-zone start */
80 int j0 = 0;
81 /* j-zone end */
82 int j1 = 0;
83 /* i-charge-group end */
84 int cg1 = 0;
85 /* j-charge-group start */
86 int jcg0 = 0;
87 /* j-charge-group end */
88 int jcg1 = 0;
89 /* Minimum shifts to consider */
90 ivec shift0 = { };
91 /* Maximum shifts to consider */
92 ivec shift1 = { };
93 } gmx_domdec_ns_ranges_t;
95 typedef struct {
96 /* Zone lower corner in triclinic coordinates */
97 rvec x0 = { };
98 /* Zone upper corner in triclinic coordinates */
99 rvec x1 = { };
100 /* Zone bounding box lower corner in Cartesian coords */
101 rvec bb_x0 = { };
102 /* Zone bounding box upper corner in Cartesian coords */
103 rvec bb_x1 = { };
104 } gmx_domdec_zone_size_t;
106 struct gmx_domdec_zones_t {
107 /* The number of zones including the home zone */
108 int n = 0;
109 /* The shift of the zones with respect to the home zone */
110 ivec shift[DD_MAXZONE] = { };
111 /* The charge group boundaries for the zones */
112 int cg_range[DD_MAXZONE+1] = { };
113 /* The number of neighbor search zones with i-particles */
114 int nizone = 0;
115 /* The neighbor search charge group ranges for each i-zone */
116 gmx_domdec_ns_ranges_t izone[DD_MAXIZONE];
117 /* Boundaries of the zones */
118 gmx_domdec_zone_size_t size[DD_MAXZONE];
119 /* The cg density of the home zone */
120 real dens_zone0 = 0;
123 struct gmx_ddbox_t {
124 int npbcdim;
125 int nboundeddim;
126 rvec box0;
127 rvec box_size;
128 /* Tells if the box is skewed for each of the three cartesian directions */
129 ivec tric_dir;
130 rvec skew_fac;
131 /* Orthogonal vectors for triclinic cells, Cartesian index */
132 rvec v[DIM][DIM];
133 /* Normal vectors for the cells walls */
134 rvec normal[DIM];
138 struct gmx_domdec_t { //NOLINT(clang-analyzer-optin.performance.Padding)
139 /* The DD particle-particle nodes only */
140 /* The communication setup within the communicator all
141 * defined in dd->comm in domdec.c
143 int nnodes;
144 MPI_Comm mpi_comm_all;
145 /* Use MPI_Sendrecv communication instead of non-blocking calls */
146 gmx_bool bSendRecv2;
147 /* The local DD cell index and rank */
148 ivec ci;
149 int rank;
150 ivec master_ci;
151 int masterrank;
152 /* Communication with the PME only nodes */
153 int pme_nodeid;
154 gmx_bool pme_receive_vir_ener;
155 gmx_pme_comm_n_box_t *cnb = nullptr;
156 int nreq_pme = 0;
157 MPI_Request req_pme[8];
160 /* The communication setup, identical for each cell, cartesian index */
161 ivec nc;
162 int ndim;
163 ivec dim; /* indexed by 0 to ndim */
165 /* TODO: Move the next 4, and more from domdec_internal.h, to a simulation system */
167 /* PBC from dim 0 (X) to npbcdim */
168 int npbcdim;
169 /* The system is bounded from 0 (X) to numBoundedDimensions */
170 int numBoundedDimensions;
171 /* Does the box size change during the simulaton? */
172 bool haveDynamicBox;
174 /* Screw PBC? */
175 gmx_bool bScrewPBC;
177 /* Forward and backward neighboring cells, indexed by 0 to ndim */
178 int neighbor[DIM][2];
180 /* Only available on the master node */
181 std::unique_ptr<AtomDistribution> ma;
183 /* Can atoms connected by constraints be assigned to different domains? */
184 bool splitConstraints;
185 /* Can atoms connected by settles be assigned to different domains? */
186 bool splitSettles;
188 /* Global atom number to interaction list */
189 gmx_reverse_top_t *reverse_top;
190 int nbonded_global;
191 int nbonded_local;
193 /* The number of inter charge-group exclusions */
194 int n_intercg_excl;
196 /* Vsite stuff */
197 gmx::HashedMap<int> *ga2la_vsite = nullptr;
198 gmx_domdec_specat_comm_t *vsite_comm = nullptr;
199 std::vector<int> vsite_requestedGlobalAtomIndices;
201 /* Constraint stuff */
202 gmx_domdec_constraints_t *constraints = nullptr;
203 gmx_domdec_specat_comm_t *constraint_comm = nullptr;
205 /* The number of home atom groups */
206 int ncg_home = 0;
207 /* Global atom group indices for the home and all non-home groups */
208 std::vector<int> globalAtomGroupIndices;
210 /* Index from the local atoms to the global atoms, covers home and received zones */
211 std::vector<int> globalAtomIndices;
213 /* Global atom number to local atom number list */
214 gmx_ga2la_t *ga2la = nullptr;
216 /* Communication stuff */
217 gmx_domdec_comm_t *comm;
219 /* The partioning count, to keep track of the state */
220 int64_t ddp_count;
222 /* The managed atom sets that are updated in domain decomposition */
223 gmx::LocalAtomSetManager * atomSets;
225 /* gmx_pme_recv_f buffer */
226 int pme_recv_f_alloc = 0;
227 rvec *pme_recv_f_buf = nullptr;
230 //! Are we the master node for domain decomposition
231 static inline bool DDMASTER(const gmx_domdec_t &dd)
233 return dd.rank == dd.masterrank;
236 //! Are we the master node for domain decomposition, deprecated
237 static inline bool DDMASTER(const gmx_domdec_t *dd)
239 return dd->rank == dd->masterrank;
242 #endif