2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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 functionality for checking whether
38 * local topologies describe all bonded interactions.
41 * \ingroup module_domdec
44 #ifndef GMX_DOMDEC_LOCALTOPOLOGYCHECKER_H
45 #define GMX_DOMDEC_LOCALTOPOLOGYCHECKER_H
49 #include "gromacs/math/vectypes.h"
51 struct gmx_localtop_t;
66 /*! \brief Has responsibility for checking that the local topology
67 * distributed across domains describes a total number of bonded
68 * interactions that matches the system topology
70 * Because this check is not urgent, the communication that it
71 * requires is done at the next opportunity, rather than requiring
72 * extra communication. If the check fails, a fatal error stops
73 * execution. In principle, if there was a bug, GROMACS might crash in
74 * the meantime because of the wrong forces. However as a bug is
75 * unlikely, we optimize by avoiding creating extra overhead from
78 class LocalTopologyChecker
81 /*! \brief Constructor
82 * \param[in] mdlog Logger
83 * \param[in] cr Communication object
84 * \param[in] mtop Global system topology
85 * \param[in] useUpdateGroups Whether update groups are in use
87 LocalTopologyChecker(const MDLogger& mdlog, const t_commrec* cr, const gmx_mtop_t& mtop, bool useUpdateGroups);
89 ~LocalTopologyChecker();
91 LocalTopologyChecker(LocalTopologyChecker&& other) noexcept;
93 LocalTopologyChecker& operator=(LocalTopologyChecker&& other) noexcept;
95 /*! \brief Set that the local topology should be checked via
96 * observables reduction whenever that reduction is required by
98 void scheduleCheckOfLocalTopology(int numBondedInteractionsToReduce);
100 /*! \brief Return whether the total bonded interaction count across
101 * domains should be checked in observables reduction this step. */
102 bool shouldCheckNumberOfBondedInteractions() const;
104 //! Return the number of bonded interactions in this domain.
105 int numBondedInteractions() const;
107 /*! \brief Set total bonded interaction count across domains. */
108 void setNumberOfBondedInteractionsOverAllDomains(int newValue);
110 /*! \brief Check whether bonded interactions are missing from the reverse topology
111 * produced by domain decomposition.
113 * \param[in] top_local Local topology for the error message
114 * \param[in] x Position vector for the error message
115 * \param[in] box Box matrix for the error message
117 void checkNumberOfBondedInteractions(const gmx_localtop_t* top_local,
118 ArrayRef<const RVec> x,
123 std::unique_ptr<Impl> impl_;