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