2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019,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 Declares the functions for generating update groups
39 * \author Berk Hess <hess@kth.se>
40 * \ingroup module_mdlib
43 #ifndef GMX_MDLIB_UPDATEGROUPS
44 #define GMX_MDLIB_UPDATEGROUPS
48 #include "gromacs/math/vectypes.h"
49 #include "gromacs/utility/arrayref.h"
50 #include "gromacs/utility/real.h"
52 enum class PbcType : int;
58 class RangePartitioning;
60 /*! \brief Returns a vector with update groups for each moleculetype in \p mtop
61 * or an empty vector when the criteria (see below) are not satisfied.
63 * An empty vector is returned when at least one moleculetype does not obey
64 * the restrictions of update groups, e.g. more than two constraints in a row.
66 * Currently valid update groups are:
67 * - a single atom which is not a virtual site and does not have constraints;
68 * - or a group of atoms where all virtual sites are constructed from atoms
69 * within the group and at least one non-vsite atom is constrained to
70 * all other non-vsite atoms.
72 * To have update groups, all virtual sites should be linear 2 or 3 atom
73 * constructions with coefficients >= 0 and sum of coefficients <= 1.
75 * This vector is generally consumed in constructing an UpdateGroups object.
77 * \param[in] mtop The system topology
79 std::vector<RangePartitioning> makeUpdateGroupingsPerMoleculeType(const gmx_mtop_t& mtop);
81 /*! \brief Returns the maximum update group radius
83 * \note When \p updateGroups is empty, 0 is returned.
85 * \param[in] mtop The system topology
86 * \param[in] updateGroupingPerMoleculeType List of update group, size should match the
87 * number of moltypes in \p mtop or be 0
88 * \param[in] temperature The maximum reference temperature, pass -1
89 * when unknown or not applicable
91 real computeMaxUpdateGroupRadius(const gmx_mtop_t& mtop,
92 ArrayRef<const RangePartitioning> updateGroupingPerMoleculeType,
95 /*! \brief Return the margin required for successful domain decomposition
97 * \param[in] pbcType The PBC type in use
98 * \param[in] box The box in use
99 * \param[in] rlist The list size in use
101 real computeCutoffMargin(PbcType pbcType, matrix box, real rlist);
103 /*! \brief Returns whether mtop contains any constraints and/or vsites
105 * When we have constraints and/or vsites, it is beneficial to use
106 * update groups (when possible) to allow independent update of
108 bool systemHasConstraintsOrVsites(const gmx_mtop_t& mtop);
111 * \brief Owns the update grouping and related data */
115 //! Default constructor
116 UpdateGroups() = default;
117 //! Constructor when update groups are active
118 UpdateGroups(std::vector<RangePartitioning>&& updateGroupingPerMoleculeType, real maxUpdateGroupRadius);
120 bool useUpdateGroups() const { return useUpdateGroups_; }
121 real maxUpdateGroupRadius() const { return maxUpdateGroupRadius_; }
122 ArrayRef<const RangePartitioning> updateGroupingPerMoleculeType() const
124 return updateGroupingPerMoleculeType_;
128 //! Whether update groups are in use
129 const bool useUpdateGroups_ = false;
130 //! The update groupings within each respective molecule type, empty when not in use
131 const std::vector<RangePartitioning> updateGroupingPerMoleculeType_ = {};
132 //! The maximum radius of any update group, 0 when not in use
133 const real maxUpdateGroupRadius_ = 0.0_real;
136 /*! \brief Builder for update groups.
138 * Checks the conditions for using update groups, and logs a message
139 * if they cannot be used, along with the reason why not.
141 * If PP domain decomposition is not in use, there is no reason to use
144 * All molecule types in the system topology must be conform to the
145 * requirements, such that makeUpdateGroupingsPerMoleculeType()
146 * returns a non-empty vector.
148 * When we have constraints and/or vsites, it is beneficial to use
149 * update groups (when possible) to allow independent update of
150 * groups. But if there are no constraints or vsites, then there is no
151 * need to use update groups at all.
153 * To use update groups, the large domain-to-domain cutoff distance
154 * should be compatible with the box size.
156 UpdateGroups makeUpdateGroups(const gmx::MDLogger& mdlog,
157 std::vector<RangePartitioning>&& updateGroupingPerMoleculeType,
158 real maxUpdateGroupRadius,
159 bool useDomainDecomposition,
160 bool systemHasConstraintsOrVsites,
165 #endif // GMX_MDLIB_UPDATEGROUPS