2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019,2020, 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.
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.
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.
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.
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.
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.
35 /*! \libinternal \file
37 * \brief This file declares functions for mdrun to call to make a new
38 * domain decomposition, and check it.
40 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_domdec
45 #ifndef GMX_DOMDEC_PARTITION_H
46 #define GMX_DOMDEC_PARTITION_H
50 #include "gromacs/math/vectypes.h"
51 #include "gromacs/utility/basedefinitions.h"
52 #include "gromacs/utility/real.h"
56 struct gmx_localtop_t
;
75 class VirtualSitesHandler
;
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
);
81 /*! \brief Print statistics for domain decomposition communication */
82 void print_dd_statistics(const t_commrec
* cr
, const t_inputrec
* ir
, FILE* fplog
);
84 /*! \brief Partition the system over the nodes.
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.
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] ir 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
113 void dd_partition_system(FILE* fplog
,
114 const gmx::MDLogger
& mdlog
,
117 gmx_bool bMasterState
,
119 t_state
* state_global
,
120 const gmx_mtop_t
& top_global
,
121 const t_inputrec
* ir
,
122 gmx::ImdSession
* imdSession
,
124 t_state
* state_local
,
125 gmx::ForceBuffers
* f
,
126 gmx::MDAtoms
* mdAtoms
,
127 gmx_localtop_t
* top_local
,
129 gmx::VirtualSitesHandler
* vsite
,
130 gmx::Constraints
* constr
,
132 gmx_wallcycle
* wcycle
,
135 /*! \brief Check whether bonded interactions are missing, if appropriate
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.
146 void checkNumberOfBondedInteractions(const gmx::MDLogger
& mdlog
,
148 int totalNumberOfBondedInteractions
,
149 const gmx_mtop_t
* top_global
,
150 const gmx_localtop_t
* top_local
,
151 gmx::ArrayRef
<const gmx::RVec
> x
,
153 bool* shouldCheckNumberOfBondedInteractions
);