Split lines with many copyright years
[gromacs.git] / src / gromacs / domdec / domdec.h
blob6608dc7d4e678cd36e0fd93192fa81efc45b90c9
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2005 - 2014, The GROMACS development team.
5 * Copyright (c) 2015,2016,2017,2018,2019 by the GROMACS development team.
6 * Copyright (c) 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 /*! \defgroup module_domdec Spatial domain decomposition (for parallelization over MPI)
38 * \ingroup group_mdrun
40 * \brief Manages the decomposition of the simulation volume over MPI
41 * ranks to try to distribute work evenly with minimal communication
42 * overheads.
44 * \todo Get domdec stuff out of mdtypes/commrec.h
46 * \author Berk Hess <hess@kth.se>
50 /*! \libinternal \file
52 * \brief This file declares functions for mdrun to call to manage the
53 * details of its domain decomposition.
55 * \author Berk Hess <hess@kth.se>
56 * \inlibraryapi
57 * \ingroup module_domdec
60 #ifndef GMX_DOMDEC_DOMDEC_H
61 #define GMX_DOMDEC_DOMDEC_H
63 #include <vector>
65 #include "gromacs/math/vectypes.h"
66 #include "gromacs/utility/arrayref.h"
67 #include "gromacs/utility/basedefinitions.h"
68 #include "gromacs/utility/real.h"
70 struct cginfo_mb_t;
71 struct gmx_domdec_t;
72 struct gmx_ddbox_t;
73 struct gmx_domdec_zones_t;
74 struct gmx_localtop_t;
75 struct gmx_mtop_t;
76 struct gmx_vsite_t;
77 struct t_block;
78 struct t_blocka;
79 struct t_commrec;
80 struct t_forcerec;
81 struct t_inputrec;
82 struct t_mdatoms;
83 struct t_nrnb;
84 struct gmx_wallcycle;
85 class t_state;
87 namespace gmx
89 class ForceWithShiftForces;
90 class MDLogger;
91 class RangePartitioning;
92 } // namespace gmx
94 /*! \brief Returns the global topology atom number belonging to local atom index i.
96 * This function is intended for writing ASCII output
97 * and returns atom numbers starting at 1.
98 * When dd=NULL returns i+1.
100 int ddglatnr(const gmx_domdec_t* dd, int i);
102 /*! \brief Returns a list of update group partitioning for each molecule type or empty when update groups are not used */
103 gmx::ArrayRef<const gmx::RangePartitioning> getUpdateGroupingPerMoleculetype(const gmx_domdec_t& dd);
105 /*! \brief Store the global cg indices of the home cgs in state,
107 * This means it can be reset, even after a new DD partitioning.
109 void dd_store_state(struct gmx_domdec_t* dd, t_state* state);
111 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
112 struct gmx_domdec_zones_t* domdec_zones(struct gmx_domdec_t* dd);
114 /*! \brief Returns the range for atoms in zones*/
115 int dd_numAtomsZones(const gmx_domdec_t& dd);
117 /*! \brief Returns the number of home atoms */
118 int dd_numHomeAtoms(const gmx_domdec_t& dd);
120 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
121 int dd_natoms_mdatoms(const gmx_domdec_t* dd);
123 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
124 int dd_natoms_vsite(const gmx_domdec_t* dd);
126 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
127 void dd_get_constraint_range(const gmx_domdec_t* dd, int* at_start, int* at_end);
129 /*! \libinternal \brief Struct for passing around the number of PME domains */
130 struct NumPmeDomains
132 int x; //!< The number of PME domains along dimension x
133 int y; //!< The number of PME domains along dimension y
136 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
137 NumPmeDomains getNumPmeDomains(const gmx_domdec_t* dd);
139 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
140 std::vector<int> get_pme_ddranks(const t_commrec* cr, int pmenodeid);
142 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
143 int dd_pme_maxshift_x(const gmx_domdec_t* dd);
145 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
146 int dd_pme_maxshift_y(const gmx_domdec_t* dd);
148 /*! \brief Return whether constraints, not including settles, cross domain boundaries */
149 bool ddHaveSplitConstraints(const gmx_domdec_t& dd);
151 /*! \brief Return whether update groups are used */
152 bool ddUsesUpdateGroups(const gmx_domdec_t& dd);
154 /*! \brief Return whether the DD has a single dimension with a single pulse
156 * The GPU halo exchange code requires a 1D single-pulse DD, and its
157 * setup code can use the returned value to understand what it should
158 * do. */
159 bool is1DAnd1PulseDD(const gmx_domdec_t& dd);
161 /*! \brief Initialize data structures for bonded interactions */
162 void dd_init_bondeds(FILE* fplog,
163 gmx_domdec_t* dd,
164 const gmx_mtop_t& mtop,
165 const gmx_vsite_t* vsite,
166 const t_inputrec* ir,
167 gmx_bool bBCheck,
168 gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
170 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
171 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t& dd);
173 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
174 gmx_bool dd_bonded_molpbc(const gmx_domdec_t* dd, int ePBC);
176 /*! \brief Change the DD non-bonded communication cut-off.
178 * This could fail when trying to increase the cut-off,
179 * then FALSE will be returned and the cut-off is not modified.
181 * \param[in] cr Communication recrod
182 * \param[in] box Box matrix, used for computing the dimensions of the system
183 * \param[in] x Position vector, used for computing the dimensions of the system
184 * \param[in] cutoffRequested The requested atom to atom cut-off distance, usually the pair-list cutoff distance
186 gmx_bool change_dd_cutoff(t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested);
188 /*! \brief Set up communication for averaging GPU wait times over domains
190 * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
191 * are meaningless, as it depends on the order in which tasks on the same
192 * GPU finish. Therefore there wait times need to be averaged over the ranks
193 * sharing the same GPU. This function sets up the communication for that.
195 void dd_setup_dlb_resource_sharing(const t_commrec* cr, int gpu_id);
197 /*! \brief Cycle counter indices used internally in the domain decomposition */
198 enum
200 ddCyclStep,
201 ddCyclPPduringPME,
202 ddCyclF,
203 ddCyclWaitGPU,
204 ddCyclPME,
205 ddCyclNr
208 /*! \brief Add the wallcycle count to the DD counter */
209 void dd_cycles_add(const gmx_domdec_t* dd, float cycles, int ddCycl);
211 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
212 void dd_move_x(struct gmx_domdec_t* dd, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_wallcycle* wcycle);
214 /*! \brief Sum the forces over the neighboring cells.
216 * When fshift!=NULL the shift forces are updated to obtain
217 * the correct virial from the single sum including f.
219 void dd_move_f(struct gmx_domdec_t* dd, gmx::ForceWithShiftForces* forceWithShiftForces, gmx_wallcycle* wcycle);
221 /*! \brief Communicate a real for each atom to the neighboring cells. */
222 void dd_atom_spread_real(struct gmx_domdec_t* dd, real v[]);
224 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
225 void dd_atom_sum_real(struct gmx_domdec_t* dd, real v[]);
227 /*! \brief Reset all the statistics and counters for total run counting */
228 void reset_dd_statistics_counters(struct gmx_domdec_t* dd);
230 /* In domdec_con.c */
232 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
233 void dd_move_f_vsites(struct gmx_domdec_t* dd, rvec* f, rvec* fshift);
235 /*! \brief Clears the forces for virtual sites */
236 void dd_clear_f_vsites(struct gmx_domdec_t* dd, rvec* f);
238 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
239 void dd_move_x_constraints(struct gmx_domdec_t* dd, const matrix box, rvec* x0, rvec* x1, gmx_bool bX1IsCoord);
241 /*! \brief Communicates the coordinates involved in virtual sites */
242 void dd_move_x_vsites(struct gmx_domdec_t* dd, const matrix box, rvec* x);
244 /*! \brief Returns the local atom count array for all constraints
246 * The local atom count for a constraint, possible values 2/1/0, is needed
247 * to avoid not/double-counting contributions linked to the Lagrange
248 * multiplier, such as the virial and free-energy derivatives.
250 * \note When \p dd = nullptr, an empty reference is returned.
252 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd);
254 /* In domdec_top.c */
256 /*! \brief Print error output when interactions are missing */
257 [[noreturn]] void dd_print_missing_interactions(const gmx::MDLogger& mdlog,
258 t_commrec* cr,
259 int local_count,
260 const gmx_mtop_t* top_global,
261 const gmx_localtop_t* top_local,
262 const rvec* x,
263 const matrix box);
265 /*! \brief Generate and store the reverse topology */
266 void dd_make_reverse_top(FILE* fplog,
267 gmx_domdec_t* dd,
268 const gmx_mtop_t* mtop,
269 const gmx_vsite_t* vsite,
270 const t_inputrec* ir,
271 gmx_bool bBCheck);
273 /*! \brief Generate the local topology and virtual site data */
274 void dd_make_local_top(struct gmx_domdec_t* dd,
275 struct gmx_domdec_zones_t* zones,
276 int npbcdim,
277 matrix box,
278 rvec cellsize_min,
279 const ivec npulse,
280 t_forcerec* fr,
281 rvec* cgcm_or_x,
282 const gmx_mtop_t& top,
283 gmx_localtop_t* ltop);
285 /*! \brief Sort ltop->ilist when we are doing free energy. */
286 void dd_sort_local_top(gmx_domdec_t* dd, const t_mdatoms* mdatoms, gmx_localtop_t* ltop);
288 /*! \brief Initialize local topology
290 * \param[in] top_global Reference to global topology.
291 * \param[in,out] top Pointer to new local topology
293 void dd_init_local_top(const gmx_mtop_t& top_global, gmx_localtop_t* top);
295 /*! \brief Construct local state */
296 void dd_init_local_state(struct gmx_domdec_t* dd, const t_state* state_global, t_state* local_state);
298 /*! \brief Generate a list of links between atoms that are linked by bonded interactions
300 * Also stores whether atoms are linked in \p cginfo_mb.
302 t_blocka* makeBondedLinks(const gmx_mtop_t& mtop, gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
304 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
305 void dd_bonded_cg_distance(const gmx::MDLogger& mdlog,
306 const gmx_mtop_t* mtop,
307 const t_inputrec* ir,
308 const rvec* x,
309 const matrix box,
310 gmx_bool bBCheck,
311 real* r_2b,
312 real* r_mb);
314 #endif