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