da3c066737a86ade3a8ce1c8866d5a4abd9c9891
[alexxy/gromacs.git] / src / gromacs / domdec / partition.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019,2020,2021, 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 /*! \libinternal \file
36  *
37  * \brief This file declares functions for mdrun to call to make a new
38  * domain decomposition, and check it.
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \inlibraryapi
42  * \ingroup module_domdec
43  */
44
45 #ifndef GMX_DOMDEC_PARTITION_H
46 #define GMX_DOMDEC_PARTITION_H
47
48 #include <cstdio>
49
50 #include "gromacs/math/vectypes.h"
51 #include "gromacs/utility/basedefinitions.h"
52 #include "gromacs/utility/real.h"
53
54 struct gmx_ddbox_t;
55 struct gmx_domdec_t;
56 struct gmx_localtop_t;
57 struct gmx_mtop_t;
58 struct gmx_wallcycle;
59 struct pull_t;
60 struct t_commrec;
61 struct t_forcerec;
62 struct t_inputrec;
63 struct t_nrnb;
64 class t_state;
65
66 namespace gmx
67 {
68 template<typename>
69 class ArrayRef;
70 class Constraints;
71 class ForceBuffers;
72 class ImdSession;
73 class MDAtoms;
74 class MDLogger;
75 class VirtualSitesHandler;
76 } // namespace gmx
77
78 //! Check whether the DD grid has moved too far for correctness.
79 bool check_grid_jump(int64_t step, const gmx_domdec_t* dd, real cutoff, const gmx_ddbox_t* ddbox, gmx_bool bFatal);
80
81 /*! \brief Print statistics for domain decomposition communication */
82 void print_dd_statistics(const t_commrec* cr, const t_inputrec& inputrec, FILE* fplog);
83
84 /*! \brief Partition the system over the nodes.
85  *
86  * step is only used for printing error messages.
87  * If bMasterState==TRUE then state_global from the master node is used,
88  * else state_local is redistributed between the nodes.
89  * When f!=NULL, *f will be reallocated to the size of state_local.
90  *
91  * \param[in] fplog         Pointer to the log file
92  * \param[in] mdlog         MD file logger
93  * \param[in] step          Current step
94  * \param[in] cr            Communication record
95  * \param[in] bMasterState  Is it a master state
96  * \param[in] nstglobalcomm Will globals be computed on this step
97  * \param[in] state_global  Global state
98  * \param[in] top_global    Global topology
99  * \param[in] inputrec      Input record
100  * \param[in] imdSession    IMD handle
101  * \param[in] pull_work     Pulling data
102  * \param[in] state_local   Local state
103  * \param[in] f             Force buffer
104  * \param[in] mdAtoms       MD atoms
105  * \param[in] top_local     Local topology
106  * \param[in] fr            Force record
107  * \param[in] vsite         Virtual sites handler
108  * \param[in] constr        Constraints
109  * \param[in] nrnb          Cycle counters
110  * \param[in] wcycle        Timers
111  * \param[in] bVerbose      Be verbose
112  */
113 void dd_partition_system(FILE*                     fplog,
114                          const gmx::MDLogger&      mdlog,
115                          int64_t                   step,
116                          const t_commrec*          cr,
117                          gmx_bool                  bMasterState,
118                          int                       nstglobalcomm,
119                          t_state*                  state_global,
120                          const gmx_mtop_t&         top_global,
121                          const t_inputrec&         inputrec,
122                          gmx::ImdSession*          imdSession,
123                          pull_t*                   pull_work,
124                          t_state*                  state_local,
125                          gmx::ForceBuffers*        f,
126                          gmx::MDAtoms*             mdAtoms,
127                          gmx_localtop_t*           top_local,
128                          t_forcerec*               fr,
129                          gmx::VirtualSitesHandler* vsite,
130                          gmx::Constraints*         constr,
131                          t_nrnb*                   nrnb,
132                          gmx_wallcycle*            wcycle,
133                          gmx_bool                  bVerbose);
134
135 /*! \brief Check whether bonded interactions are missing, if appropriate
136  *
137  * \param[in]    mdlog                                  Logger
138  * \param[in]    cr                                     Communication object
139  * \param[in]    totalNumberOfBondedInteractions        Result of the global reduction over the number of bonds treated in each domain
140  * \param[in]    top_global                             Global topology for the error message
141  * \param[in]    top_local                              Local topology for the error message
142  * \param[in]    x                                      Position vector for the error message
143  * \param[in]    box                                    Box matrix for the error message
144  * \param[in,out] shouldCheckNumberOfBondedInteractions Whether we should do the check. Always set to false.
145  */
146 void checkNumberOfBondedInteractions(const gmx::MDLogger&           mdlog,
147                                      t_commrec*                     cr,
148                                      int                            totalNumberOfBondedInteractions,
149                                      const gmx_mtop_t&              top_global,
150                                      const gmx_localtop_t*          top_local,
151                                      gmx::ArrayRef<const gmx::RVec> x,
152                                      const matrix                   box,
153                                      bool* shouldCheckNumberOfBondedInteractions);
154
155 #endif