Remove gmx custom fixed int (e.g. gmx_int64_t) types
[alexxy/gromacs.git] / src / gromacs / domdec / domdec.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2005,2006,2007,2008,2009,2010,2012,2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \defgroup module_domdec Spatial domain decomposition (for parallelization over MPI)
36  * \ingroup group_mdrun
37  *
38  * \brief Manages the decomposition of the simulation volume over MPI
39  * ranks to try to distribute work evenly with minimal communication
40  * overheads.
41  *
42  * \todo Get domdec stuff out of mdtypes/commrec.h
43  *
44  * \author Berk Hess <hess@kth.se>
45  *
46  */
47
48 /*! \libinternal \file
49  *
50  * \brief This file declares functions for mdrun to call to manage the
51  * details of its domain decomposition.
52  *
53  * \author Berk Hess <hess@kth.se>
54  * \inlibraryapi
55  * \ingroup module_domdec
56  */
57
58 #ifndef GMX_DOMDEC_DOMDEC_H
59 #define GMX_DOMDEC_DOMDEC_H
60
61 #include <vector>
62
63 #include "gromacs/math/paddedvector.h"
64 #include "gromacs/math/vectypes.h"
65 #include "gromacs/utility/arrayref.h"
66 #include "gromacs/utility/basedefinitions.h"
67 #include "gromacs/utility/real.h"
68
69 struct cginfo_mb_t;
70 struct gmx_domdec_t;
71 struct gmx_ddbox_t;
72 struct gmx_domdec_zones_t;
73 struct gmx_enfrot;
74 struct gmx_localtop_t;
75 struct gmx_mtop_t;
76 struct gmx_vsite_t;
77 struct MdrunOptions;
78 struct t_block;
79 struct t_blocka;
80 struct t_commrec;
81 struct t_forcerec;
82 struct t_inputrec;
83 struct t_mdatoms;
84 struct t_nrnb;
85 struct gmx_wallcycle;
86 class t_state;
87
88 namespace gmx
89 {
90 class Constraints;
91 class MDAtoms;
92 class LocalAtomSetManager;
93 } // namespace
94
95 /*! \brief Returns the global topology atom number belonging to local atom index i.
96  *
97  * This function is intended for writing ASCII output
98  * and returns atom numbers starting at 1.
99  * When dd=NULL returns i+1.
100  */
101 int ddglatnr(const gmx_domdec_t *dd, int i);
102
103 /*! \brief Return a block struct for the charge groups of the whole system */
104 t_block *dd_charge_groups_global(struct gmx_domdec_t *dd);
105
106 /*! \brief Store the global cg indices of the home cgs in state,
107  *
108  * This means it can be reset, even after a new DD partitioning.
109  */
110 void dd_store_state(struct gmx_domdec_t *dd, t_state *state);
111
112 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
113 struct gmx_domdec_zones_t *domdec_zones(struct gmx_domdec_t *dd);
114
115 /*! \brief Sets the j-charge-group range for i-charge-group \p icg */
116 void dd_get_ns_ranges(const gmx_domdec_t *dd, int icg,
117                       int *jcg0, int *jcg1, ivec shift0, ivec shift1);
118
119 /*! \brief Returns the number of home atoms */
120 int dd_numHomeAtoms(const gmx_domdec_t &dd);
121
122 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
123 int dd_natoms_mdatoms(const gmx_domdec_t *dd);
124
125 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
126 int dd_natoms_vsite(const gmx_domdec_t *dd);
127
128 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
129 void dd_get_constraint_range(const gmx_domdec_t *dd,
130                              int *at_start, int *at_end);
131
132 /*! \libinternal \brief Struct for passing around the number of PME domains */
133 struct NumPmeDomains
134 {
135     int x; //!< The number of PME domains along dimension x
136     int y; //!< The number of PME domains along dimension y
137 };
138
139 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
140 NumPmeDomains getNumPmeDomains(const gmx_domdec_t *dd);
141
142 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
143 std::vector<int> get_pme_ddranks(const t_commrec *cr, int pmenodeid);
144
145 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
146 int dd_pme_maxshift_x(const gmx_domdec_t *dd);
147
148 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
149 int dd_pme_maxshift_y(const gmx_domdec_t *dd);
150
151 /*! \brief The options for the domain decomposition MPI task ordering. */
152 enum class DdRankOrder
153 {
154     select,     //!< First value (needed to cope with command-line parsing)
155     interleave, //!< Interleave the PP and PME ranks
156     pp_pme,     //!< First all PP ranks, all PME rank at the end
157     cartesian,  //!< Use Cartesian communicators for PP, PME and PP-PME
158     nr          //!< The number of options
159 };
160
161 /*! \brief The options for the dynamic load balancing. */
162 enum class DlbOption
163 {
164     select,           //!< First value (needed to cope with command-line parsing)
165     turnOnWhenUseful, //!< Turn on DLB when we think it would improve performance
166     no,               //!< Never turn on DLB
167     yes,              //!< Turn on DLB from the start and keep it on
168     nr                //!< The number of options
169 };
170
171 /*! \libinternal \brief Structure containing all (command line) options for the domain decomposition */
172 struct DomdecOptions
173 {
174     /*! \brief Constructor */
175     DomdecOptions();
176
177     //! If true, check that all bonded interactions have been assigned to exactly one domain/rank.
178     gmx_bool          checkBondedInteractions;
179     //! If true, don't communicate all atoms between the non-bonded cut-off and the larger bonded cut-off, but only those that have non-local bonded interactions. This significantly reduces the communication volume.
180     gmx_bool          useBondedCommunication;
181     //! The domain decomposition grid cell count, 0 means let domdec choose based on the number of ranks.
182     ivec              numCells;
183     //! The number of separate PME ranks requested, -1 = auto.
184     int               numPmeRanks;
185     //! Ordering of the PP and PME ranks, values from enum above.
186     DdRankOrder       rankOrder;
187     //! The minimum communication range, used for extended the communication range for bonded interactions (nm).
188     real              minimumCommunicationRange;
189     //! Communication range for atom involved in constraints (P-LINCS) (nm).
190     real              constraintCommunicationRange;
191     //! Dynamic load balancing option, values from enum above.
192     DlbOption         dlbOption;
193     /*! \brief Fraction in (0,1) by whose reciprocal the initial
194      * DD cell size will be increased in order to provide a margin
195      * in which dynamic load balancing can act, while preserving
196      * the minimum cell size. */
197     real              dlbScaling;
198     //! String containing a vector of the relative sizes in the x direction of the corresponding DD cells.
199     const char       *cellSizeX;
200     //! String containing a vector of the relative sizes in the y direction of the corresponding DD cells.
201     const char       *cellSizeY;
202     //! String containing a vector of the relative sizes in the z direction of the corresponding DD cells.
203     const char       *cellSizeZ;
204 };
205
206 /*! \brief Initialized the domain decomposition, chooses the DD grid and PME ranks, return the DD struct */
207 gmx_domdec_t *
208 init_domain_decomposition(FILE                           *fplog,
209                           t_commrec                      *cr,
210                           const DomdecOptions            &options,
211                           const MdrunOptions             &mdrunOptions,
212                           const gmx_mtop_t               *mtop,
213                           const t_inputrec               *ir,
214                           const matrix                    box,
215                           gmx::ArrayRef<const gmx::RVec>  xGlobal,
216                           gmx::LocalAtomSetManager       *atomSets);
217
218 /*! \brief Initialize data structures for bonded interactions */
219 void dd_init_bondeds(FILE              *fplog,
220                      gmx_domdec_t      *dd,
221                      const gmx_mtop_t  *mtop,
222                      const gmx_vsite_t *vsite,
223                      const t_inputrec  *ir,
224                      gmx_bool           bBCheck,
225                      cginfo_mb_t       *cginfo_mb);
226
227 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
228 gmx_bool dd_bonded_molpbc(const gmx_domdec_t *dd, int ePBC);
229
230 /*! \brief Change the DD non-bonded communication cut-off.
231  *
232  * This could fail when trying to increase the cut-off,
233  * then FALSE will be returned and the cut-off is not modified.
234  */
235 gmx_bool change_dd_cutoff(struct t_commrec *cr,
236                           t_state *state, const t_inputrec *ir,
237                           real cutoff_req );
238
239 /*! \brief Limit DLB to preserve the option of returning to the current cut-off.
240  *
241  * Domain boundary changes due to the DD dynamic load balancing can limit
242  * the cut-off distance that can be set in change_dd_cutoff. This function
243  * sets/changes the DLB limit such that using the passed (pair-list) cut-off
244  * should still be possible after subsequently setting a shorter cut-off
245  * with change_dd_cutoff.
246  */
247 void set_dd_dlb_max_cutoff(struct t_commrec *cr, real cutoff);
248
249 /*! \brief Return if we are currently using dynamic load balancing */
250 gmx_bool dd_dlb_is_on(const struct gmx_domdec_t *dd);
251
252 /*! \brief Return if the DLB lock is set */
253 gmx_bool dd_dlb_is_locked(const struct gmx_domdec_t *dd);
254
255 /*! \brief Set a lock such that with DLB=auto DLB cannot get turned on */
256 void dd_dlb_lock(struct gmx_domdec_t *dd);
257
258 /*! \brief Clear a lock such that with DLB=auto DLB may get turned on later */
259 void dd_dlb_unlock(struct gmx_domdec_t *dd);
260
261 /*! \brief Set up communication for averaging GPU wait times over domains
262  *
263  * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
264  * are meaningless, as it depends on the order in which tasks on the same
265  * GPU finish. Therefore there wait times need to be averaged over the ranks
266  * sharing the same GPU. This function sets up the communication for that.
267  */
268 void dd_setup_dlb_resource_sharing(t_commrec           *cr,
269                                    int                  gpu_id);
270
271 /*! \brief Cycle counter indices used internally in the domain decomposition */
272 enum {
273     ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr
274 };
275
276 /*! \brief Add the wallcycle count to the DD counter */
277 void dd_cycles_add(const gmx_domdec_t *dd, float cycles, int ddCycl);
278
279 /*! \brief Start the force flop count */
280 void dd_force_flop_start(struct gmx_domdec_t *dd, t_nrnb *nrnb);
281
282 /*! \brief Stop the force flop count */
283 void dd_force_flop_stop(struct gmx_domdec_t *dd, t_nrnb *nrnb);
284
285 /*! \brief Return the PME/PP force load ratio, or -1 if nothing was measured.
286  *
287  * Should only be called on the DD master node.
288  */
289 float dd_pme_f_ratio(struct gmx_domdec_t *dd);
290
291 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
292 void dd_move_x(struct gmx_domdec_t      *dd,
293                matrix                    box,
294                gmx::ArrayRef<gmx::RVec>  x,
295                gmx_wallcycle            *wcycle);
296
297 /*! \brief Sum the forces over the neighboring cells.
298  *
299  * When fshift!=NULL the shift forces are updated to obtain
300  * the correct virial from the single sum including f.
301  */
302 void dd_move_f(struct gmx_domdec_t      *dd,
303                gmx::ArrayRef<gmx::RVec>  f,
304                rvec                     *fshift,
305                gmx_wallcycle            *wcycle);
306
307 /*! \brief Communicate a real for each atom to the neighboring cells. */
308 void dd_atom_spread_real(struct gmx_domdec_t *dd, real v[]);
309
310 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
311 void dd_atom_sum_real(struct gmx_domdec_t *dd, real v[]);
312
313 /*! \brief Partition the system over the nodes.
314  *
315  * step is only used for printing error messages.
316  * If bMasterState==TRUE then state_global from the master node is used,
317  * else state_local is redistributed between the nodes.
318  * When f!=NULL, *f will be reallocated to the size of state_local.
319  */
320 void dd_partition_system(FILE                *fplog,
321                          int64_t              step,
322                          const t_commrec     *cr,
323                          gmx_bool             bMasterState,
324                          int                  nstglobalcomm,
325                          t_state             *state_global,
326                          const gmx_mtop_t    *top_global,
327                          const t_inputrec    *ir,
328                          gmx_enfrot          *enforcedRotation,
329                          t_state             *state_local,
330                          PaddedRVecVector    *f,
331                          gmx::MDAtoms        *mdatoms,
332                          gmx_localtop_t      *top_local,
333                          t_forcerec          *fr,
334                          gmx_vsite_t         *vsite,
335                          gmx::Constraints    *constr,
336                          t_nrnb              *nrnb,
337                          gmx_wallcycle       *wcycle,
338                          gmx_bool             bVerbose);
339
340 /*! \brief Reset all the statistics and counters for total run counting */
341 void reset_dd_statistics_counters(struct gmx_domdec_t *dd);
342
343 /*! \brief Print statistics for domain decomposition communication */
344 void print_dd_statistics(const t_commrec *cr, const t_inputrec *ir, FILE *fplog);
345
346 /*! \brief Check whether bonded interactions are missing, if appropriate
347  *
348  * \param[in]    fplog                                  Log file pointer
349  * \param[in]    cr                                     Communication object
350  * \param[in]    totalNumberOfBondedInteractions        Result of the global reduction over the number of bonds treated in each domain
351  * \param[in]    top_global                             Global topology for the error message
352  * \param[in]    top_local                              Local topology for the error message
353  * \param[in]    state                                  Global state for the error message
354  * \param[in,out] shouldCheckNumberOfBondedInteractions Whether we should do the check. Always set to false.
355  */
356 void checkNumberOfBondedInteractions(FILE                 *fplog,
357                                      t_commrec            *cr,
358                                      int                   totalNumberOfBondedInteractions,
359                                      const gmx_mtop_t     *top_global,
360                                      const gmx_localtop_t *top_local,
361                                      const t_state        *state,
362                                      bool                 *shouldCheckNumberOfBondedInteractions);
363
364 /* In domdec_con.c */
365
366 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
367 void dd_move_f_vsites(struct gmx_domdec_t *dd, rvec *f, rvec *fshift);
368
369 /*! \brief Clears the forces for virtual sites */
370 void dd_clear_f_vsites(struct gmx_domdec_t *dd, rvec *f);
371
372 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
373 void dd_move_x_constraints(struct gmx_domdec_t *dd, matrix box,
374                            rvec *x0, rvec *x1, gmx_bool bX1IsCoord);
375
376 /*! \brief Communicates the coordinates involved in virtual sites */
377 void dd_move_x_vsites(struct gmx_domdec_t *dd, const matrix box, rvec *x);
378
379 /*! \brief Returns the local atom count array for all constraints
380  *
381  * The local atom count for a constraint, possible values 2/1/0, is needed
382  * to avoid not/double-counting contributions linked to the Lagrange
383  * multiplier, such as the virial and free-energy derivatives.
384  *
385  * \note When \p dd = nullptr, an empty reference is returned.
386  */
387 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t *dd);
388
389 /* In domdec_top.c */
390
391 /*! \brief Print error output when interactions are missing */
392 [[ noreturn ]]
393 void dd_print_missing_interactions(FILE *fplog, struct t_commrec *cr,
394                                    int local_count,
395                                    const gmx_mtop_t *top_global,
396                                    const gmx_localtop_t *top_local,
397                                    const t_state *state_local);
398
399 /*! \brief Generate and store the reverse topology */
400 void dd_make_reverse_top(FILE *fplog,
401                          gmx_domdec_t *dd, const gmx_mtop_t *mtop,
402                          const gmx_vsite_t *vsite,
403                          const t_inputrec *ir, gmx_bool bBCheck);
404
405 /*! \brief Store the local charge group index in \p lcgs */
406 void dd_make_local_cgs(struct gmx_domdec_t *dd, t_block *lcgs);
407
408 /*! \brief Generate the local topology and virtual site data */
409 void dd_make_local_top(struct gmx_domdec_t *dd, struct gmx_domdec_zones_t *zones,
410                        int npbcdim, matrix box,
411                        rvec cellsize_min, const ivec npulse,
412                        t_forcerec *fr,
413                        rvec *cgcm_or_x,
414                        gmx_vsite_t *vsite,
415                        const gmx_mtop_t *top, gmx_localtop_t *ltop);
416
417 /*! \brief Sort ltop->ilist when we are doing free energy. */
418 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
419                        gmx_localtop_t *ltop);
420
421 /*! \brief Construct local topology */
422 gmx_localtop_t *dd_init_local_top(const gmx_mtop_t *top_global);
423
424 /*! \brief Construct local state */
425 void dd_init_local_state(struct gmx_domdec_t *dd,
426                          t_state *state_global, t_state *local_state);
427
428 /*! \brief Generate a list of links between charge groups that are linked by bonded interactions */
429 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
430                                   cginfo_mb_t *cginfo_mb);
431
432 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
433 void dd_bonded_cg_distance(FILE *fplog, const gmx_mtop_t *mtop,
434                            const t_inputrec *ir,
435                            const rvec *x, const matrix box,
436                            gmx_bool bBCheck,
437                            real *r_2b, real *r_mb);
438
439 /*! \brief Dump a pdb file with the current DD home + communicated atoms.
440  *
441  * When natoms=-1, dump all known atoms.
442  */
443 void write_dd_pdb(const char *fn, int64_t step, const char *title,
444                   const gmx_mtop_t *mtop,
445                   const t_commrec *cr,
446                   int natoms, const rvec x[], const matrix box);
447
448
449 /* In domdec_setup.c */
450
451 /*! \brief Returns the volume fraction of the system that is communicated */
452 real comm_box_frac(const ivec dd_nc, real cutoff, const gmx_ddbox_t *ddbox);
453
454 /*! \brief Determines the optimal DD cell setup dd->nc and possibly npmenodes
455  * for the system.
456  *
457  * On the master node returns the actual cellsize limit used.
458  */
459 real dd_choose_grid(FILE *fplog,
460                     t_commrec *cr, gmx_domdec_t *dd,
461                     const t_inputrec *ir,
462                     const gmx_mtop_t *mtop,
463                     const matrix box, const gmx_ddbox_t *ddbox,
464                     int nPmeRanks,
465                     gmx_bool bDynLoadBal, real dlb_scale,
466                     real cellsize_limit, real cutoff_dd,
467                     gmx_bool bInterCGBondeds);
468
469
470 /* In domdec_box.c */
471
472 /*! \brief Set the box and PBC data in \p ddbox */
473 void set_ddbox(gmx_domdec_t *dd, bool masterRankHasTheSystemState,
474                const t_inputrec *ir, const matrix box,
475                bool calculateUnboundedSize,
476                gmx::ArrayRef<const gmx::RVec> x,
477                gmx_ddbox_t *ddbox);
478
479 /*! \brief Set the box and PBC data in \p ddbox */
480 void set_ddbox_cr(const t_commrec *cr, const ivec *dd_nc,
481                   const t_inputrec *ir, const matrix box,
482                   gmx::ArrayRef<const gmx::RVec> x,
483                   gmx_ddbox_t *ddbox);
484
485 #endif