2 * This file is part of the GROMACS molecular simulation package.
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.
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.
37 * \brief Implements functions for generating update groups
39 * \author Berk Hess <hess@kth.se>
40 * \ingroup module_mdlib
45 #include "updategroups.h"
49 #include <unordered_map>
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/mdlib/constr.h"
55 #include "gromacs/pbcutil/pbc.h"
56 #include "gromacs/topology/idef.h"
57 #include "gromacs/topology/ifunc.h"
58 #include "gromacs/topology/mtop_util.h"
59 #include "gromacs/topology/topology.h"
60 #include "gromacs/utility/listoflists.h"
61 #include "gromacs/utility/logger.h"
62 #include "gromacs/utility/message_string_collector.h"
67 /*! \brief Returns whether \p moltype contains flexible constraints */
68 static bool hasFlexibleConstraints(const gmx_moltype_t& moltype, gmx::ArrayRef<const t_iparams> iparams)
70 for (auto& ilist : extractILists(moltype.ilist, IF_CONSTRAINT))
72 if (ilist.functionType != F_SETTLE)
74 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
76 if (isConstraintFlexible(iparams, ilist.iatoms[i]))
87 /*! \brief Returns whether moltype has incompatible vsites.
89 * For simplicity the only compatible vsites are linear 2 or 3 atom sites
90 * that are constructed in between the 2 or 3 contructing atoms,
92 static bool hasIncompatibleVsites(const gmx_moltype_t& moltype, gmx::ArrayRef<const t_iparams> iparams)
94 bool hasIncompatibleVsites = false;
96 for (auto& ilist : extractILists(moltype.ilist, IF_VSITE))
98 if (ilist.functionType == F_VSITE2 || ilist.functionType == F_VSITE3)
100 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
102 const t_iparams& iparam = iparams[ilist.iatoms[i]];
105 if (ilist.functionType == F_VSITE2)
107 coeffMin = iparam.vsite.a;
108 coeffSum = iparam.vsite.a;
112 coeffMin = std::min(iparam.vsite.a, iparam.vsite.b);
113 coeffSum = iparam.vsite.a + iparam.vsite.b;
115 if (coeffMin < 0 || coeffSum > 1)
117 hasIncompatibleVsites = true;
124 hasIncompatibleVsites = true;
129 return hasIncompatibleVsites;
132 /*! \brief Returns a merged list with constraints of all types */
133 static InteractionList jointConstraintList(const gmx_moltype_t& moltype)
135 InteractionList ilistCombined;
136 std::vector<int>& iatoms = ilistCombined.iatoms;
138 for (auto& ilist : extractILists(moltype.ilist, IF_CONSTRAINT))
140 if (ilist.functionType == F_SETTLE)
142 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
144 iatoms.push_back(-1);
145 iatoms.push_back(ilist.iatoms[i + 1]);
146 iatoms.push_back(ilist.iatoms[i + 2]);
147 iatoms.push_back(-1);
148 iatoms.push_back(ilist.iatoms[i + 1]);
149 iatoms.push_back(ilist.iatoms[i + 3]);
150 iatoms.push_back(-1);
151 iatoms.push_back(ilist.iatoms[i + 2]);
152 iatoms.push_back(ilist.iatoms[i + 3]);
157 GMX_RELEASE_ASSERT(NRAL(ilist.functionType) == 2,
158 "Can only handle two-atom non-SETTLE constraints");
160 iatoms.insert(iatoms.end(), ilist.iatoms.begin(), ilist.iatoms.end());
164 return ilistCombined;
167 /*! \brief Struct for returning an atom range */
168 struct AtomIndexExtremes
170 int minAtom; //!< The minimum atom index
171 int maxAtom; //!< The maximum atom index
174 /*! \brief Returns the range of constructing atom for vsite with atom index \p a */
175 static AtomIndexExtremes vsiteConstructRange(int a, const gmx_moltype_t& moltype)
177 AtomIndexExtremes extremes = { -1, -1 };
179 for (auto& ilist : extractILists(moltype.ilist, IF_VSITE))
181 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
183 if (ilist.iatoms[i + 1] == a)
185 extremes.minAtom = ilist.iatoms[i + 2];
186 extremes.maxAtom = ilist.iatoms[i + 2];
187 for (size_t j = i + 3; j < i + ilistStride(ilist); j++)
189 extremes.minAtom = std::min(extremes.minAtom, ilist.iatoms[j]);
190 extremes.maxAtom = std::max(extremes.maxAtom, ilist.iatoms[j]);
197 GMX_RELEASE_ASSERT(false, "If a is a vsite, we should have found constructing atoms");
202 /*! \brief Returns the range of atoms constrained to atom \p a (including \p a itself) */
203 static AtomIndexExtremes constraintAtomRange(int a,
204 const ListOfLists<int>& at2con,
205 const InteractionList& ilistConstraints)
207 AtomIndexExtremes extremes = { a, a };
209 for (const int constraint : at2con[a])
211 for (int j = 0; j < 2; j++)
213 int atomJ = ilistConstraints.iatoms[constraint * 3 + 1 + j];
214 extremes.minAtom = std::min(extremes.minAtom, atomJ);
215 extremes.maxAtom = std::max(extremes.maxAtom, atomJ);
222 /*! \brief Returns a list that tells whether atoms in \p moltype are vsites */
223 static std::vector<bool> buildIsParticleVsite(const gmx_moltype_t& moltype)
225 std::vector<bool> isVsite(moltype.atoms.nr);
227 for (auto& ilist : extractILists(moltype.ilist, IF_VSITE))
229 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
231 int vsiteAtom = ilist.iatoms[i + 1];
232 isVsite[vsiteAtom] = true;
239 /*! \brief Returns the size of the update group starting at \p firstAtom or 0 when criteria (see updategroups.h) are not met */
240 static int detectGroup(int firstAtom,
241 const gmx_moltype_t& moltype,
242 const ListOfLists<int>& at2con,
243 const InteractionList& ilistConstraints)
245 /* We should be using moltype.atoms.atom[].ptype for checking whether
246 * a particle is a vsite. But the test code can't fill t_atoms,
247 * because it uses C pointers which get double freed.
249 std::vector<bool> isParticleVsite = buildIsParticleVsite(moltype);
251 /* A non-vsite atom without constraints is an update group by itself */
252 if (!isParticleVsite[firstAtom] && at2con[firstAtom].empty())
257 /* A (potential) update group starts at firstAtom and should consist
258 * of two or more atoms and possibly vsites. At least one atom should
259 * have constraints with all other atoms and vsites should have all
260 * constructing atoms inside the group. Here we increase lastAtom until
261 * the criteria are fulfilled or exit when criteria cannot be met.
263 int numAtomsWithConstraints = 0;
264 int maxConstraintsPerAtom = 0;
265 int lastAtom = firstAtom;
267 while (a <= lastAtom)
269 if (isParticleVsite[a])
271 AtomIndexExtremes extremes = vsiteConstructRange(a, moltype);
272 if (extremes.minAtom < firstAtom)
274 /* A constructing atom is outside the group,
275 * we can not use update groups.
279 lastAtom = std::max(lastAtom, extremes.maxAtom);
283 const int numConstraints = at2con[a].ssize();
284 if (numConstraints == 0)
286 /* We can not have unconstrained atoms in an update group */
289 /* This atom has at least one constraint.
290 * Check whether all constraints are within the group
291 * and whether we need to extend the group.
293 numAtomsWithConstraints += 1;
294 maxConstraintsPerAtom = std::max(maxConstraintsPerAtom, numConstraints);
296 AtomIndexExtremes extremes = constraintAtomRange(a, at2con, ilistConstraints);
297 if (extremes.minAtom < firstAtom)
299 /* Constraint to atom outside the "group" */
302 lastAtom = std::max(lastAtom, extremes.maxAtom);
308 /* lastAtom might be followed by a vsite that is constructed from atoms
309 * with index <= lastAtom. Handle this case.
311 if (lastAtom + 1 < moltype.atoms.nr && isParticleVsite[lastAtom + 1])
313 AtomIndexExtremes extremes = vsiteConstructRange(lastAtom + 1, moltype);
314 if (extremes.minAtom < firstAtom)
315 { // NOLINT bugprone-branch-clone
316 /* Constructing atom precedes the group */
319 else if (extremes.maxAtom <= lastAtom)
321 /* All constructing atoms are in the group, add the vsite to the group */
324 else if (extremes.minAtom <= lastAtom)
326 /* Some, but not all constructing atoms are in the group */
331 GMX_RELEASE_ASSERT(maxConstraintsPerAtom < numAtomsWithConstraints,
332 "We have checked that atoms are only constrained to atoms within the group,"
333 "so each atom should have fewer constraints than the number of atoms");
334 /* Check that at least one atom is constrained to all others */
335 if (maxConstraintsPerAtom != numAtomsWithConstraints - 1)
340 return lastAtom - firstAtom + 1;
343 /*! \brief Returns a list of update groups for \p moltype */
344 static RangePartitioning makeUpdateGroupingsPerMoleculeType(const gmx_moltype_t& moltype,
345 gmx::ArrayRef<const t_iparams> iparams)
347 RangePartitioning groups;
349 /* Update groups are not compatible with flexible constraints.
350 * Without dynamics the flexible constraints are ignored,
351 * but since performance for EM/NM is less critical, we do not
352 * use update groups to keep the code here simpler.
354 if (hasFlexibleConstraints(moltype, iparams) || hasIncompatibleVsites(moltype, iparams))
359 /* Combine all constraint ilists into a single one */
360 std::array<InteractionList, F_NRE> ilistsCombined;
361 ilistsCombined[F_CONSTR] = jointConstraintList(moltype);
362 /* We "include" flexible constraints, but none are present (checked above) */
363 const ListOfLists<int> at2con = make_at2con(
364 moltype.atoms.nr, ilistsCombined, iparams, FlexibleConstraintTreatment::Include);
366 bool satisfiesCriteria = true;
369 while (satisfiesCriteria && firstAtom < moltype.atoms.nr)
371 int numAtomsInGroup = detectGroup(firstAtom, moltype, at2con, ilistsCombined[F_CONSTR]);
373 if (numAtomsInGroup == 0)
375 satisfiesCriteria = false;
379 groups.appendBlock(numAtomsInGroup);
381 firstAtom += numAtomsInGroup;
384 if (!satisfiesCriteria)
386 /* Make groups empty, to signal not satisfying the criteria */
393 std::vector<RangePartitioning> makeUpdateGroupingsPerMoleculeType(const gmx_mtop_t& mtop)
395 std::vector<RangePartitioning> updateGroupingsPerMoleculeType;
397 bool systemSatisfiesCriteria = true;
398 for (const gmx_moltype_t& moltype : mtop.moltype)
400 updateGroupingsPerMoleculeType.push_back(
401 makeUpdateGroupingsPerMoleculeType(moltype, mtop.ffparams.iparams));
403 if (updateGroupingsPerMoleculeType.back().numBlocks() == 0)
405 systemSatisfiesCriteria = false;
409 if (!systemSatisfiesCriteria)
411 updateGroupingsPerMoleculeType.clear();
414 return updateGroupingsPerMoleculeType;
417 /*! \brief Returns a map of angles ilist.iatoms indices with the middle atom as key */
418 static std::unordered_multimap<int, int> getAngleIndices(const gmx_moltype_t& moltype)
420 const InteractionList& angles = moltype.ilist[F_ANGLES];
422 std::unordered_multimap<int, int> indices(angles.size());
424 for (int i = 0; i < angles.size(); i += 1 + NRAL(F_ANGLES))
426 indices.insert({ angles.iatoms[i + 2], i });
432 /*! \brief When possible, computes the maximum radius of constrained atom in an update group
434 * Supports groups with 2 or 3 atoms where all partner atoms are connected to
435 * each other by angle potentials. The temperature is used to compute a radius
436 * that is not exceeded with a chance of 10^-9. Note that this computation
437 * assumes there are no other strong forces working on these angular
438 * degrees of freedom.
439 * The return value is -1 when all partners are not connected to each other
440 * by one angle potential, when a potential is perturbed or when an angle
441 * could reach more than 180 degrees.
443 template<int numPartnerAtoms>
444 static real constraintGroupRadius(const gmx_moltype_t& moltype,
445 gmx::ArrayRef<const t_iparams> iparams,
446 const int centralAtom,
447 const ListOfLists<int>& at2con,
448 const std::unordered_multimap<int, int>& angleIndices,
449 const real constraintLength,
450 const real temperature)
452 const int numConstraints = at2con[centralAtom].ssize();
453 GMX_RELEASE_ASSERT(numConstraints == numPartnerAtoms,
454 "We expect as many constraints as partner atoms here");
456 std::array<int, numPartnerAtoms> partnerAtoms;
457 for (int i = 0; i < numPartnerAtoms; i++)
459 const int ind = at2con[centralAtom][i] * 3;
460 if (ind >= moltype.ilist[F_CONSTR].size())
462 /* This is a flexible constraint, we don't optimize for that */
465 const int a1 = moltype.ilist[F_CONSTR].iatoms[ind + 1];
466 const int a2 = moltype.ilist[F_CONSTR].iatoms[ind + 2];
467 partnerAtoms[i] = (a1 == centralAtom ? a2 : a1);
470 const InteractionList& angles = moltype.ilist[F_ANGLES];
471 auto range = angleIndices.equal_range(centralAtom);
473 std::array<int, numPartnerAtoms> numAngles = { 0 };
474 bool areSameType = true;
475 for (auto it = range.first; it != range.second; ++it)
477 /* Check if the outer atoms in the angle are both partner atoms */
478 int numAtomsFound = 0;
479 for (int ind = it->second + 1; ind < it->second + 4; ind += 2)
481 for (const int& partnerAtom : partnerAtoms)
483 if (angles.iatoms[ind] == partnerAtom)
489 if (numAtomsFound == 2)
491 /* Check that the angle potentials have the same type */
494 angleType = angles.iatoms[it->second];
496 else if (angles.iatoms[it->second] != angleType)
500 /* Count the number of angle interactions per atoms */
501 for (int ind = it->second + 1; ind < it->second + 4; ind += 2)
503 for (size_t i = 0; i < partnerAtoms.size(); i++)
505 if (angles.iatoms[ind] == partnerAtoms[i])
514 bool criteriaSatisfied = areSameType;
515 for (int i = 0; i < numPartnerAtoms; i++)
517 if (numAngles[i] != numPartnerAtoms - 1)
519 criteriaSatisfied = false;
523 /* We don't bother optimizing the perturbed angle case */
524 const t_iparams& angleParams = iparams[angleType];
525 if (criteriaSatisfied && angleParams.harmonic.rB == angleParams.harmonic.rA
526 && angleParams.harmonic.krB == angleParams.harmonic.krA)
528 /* Set number of stddevs such that change of exceeding < 10^-9 */
529 constexpr real c_numSigma = 6.0;
530 /* Compute the maximally stretched angle */
531 const real eqAngle = angleParams.harmonic.rA * gmx::c_deg2Rad;
532 const real fc = angleParams.harmonic.krA;
533 const real maxAngle =
534 eqAngle + c_numSigma * gmx::c_boltz * temperature / ((numPartnerAtoms - 1) * fc);
535 if (maxAngle >= M_PI)
540 if (numPartnerAtoms == 2)
542 /* With two atoms constrainted to a cental atom we have a triangle
543 * with two equal sides because the constraint type is equal.
544 * Return the distance from the COG to the farthest two corners,
545 * i.e. the partner atoms.
547 real distMidPartner = std::sin(0.5 * maxAngle) * constraintLength;
548 real distCentralMid = std::cos(0.5 * maxAngle) * constraintLength;
549 real distCogMid = distCentralMid * numPartnerAtoms / (numPartnerAtoms + 1);
550 real distCogPartner = std::sqrt(gmx::square(distMidPartner) + gmx::square(distCogMid));
552 return distCogPartner;
554 else if (numPartnerAtoms == 3)
556 /* With three atoms constrainted to a cental atom we have the
557 * largest COG-atom distance when one partner atom (the outlier)
558 * moves out by stretching its two angles by an equal amount.
559 * The math here gets a bit more involved, but it is still
560 * rather straightforward geometry.
561 * We first compute distances in the plane of the three partners.
562 * Here we have two "equilibrium" partners and one outlier.
563 * We make use of the "Mid" point between the two "Eq" partners.
564 * We project the central atom on this plane.
565 * Then we compute the distance of the central atom to the plane.
566 * The distance of the COG to the ourlier is returned.
568 real halfDistEqPartner = std::sin(0.5 * eqAngle) * constraintLength;
569 real distPartnerOutlier = 2 * std::sin(0.5 * maxAngle) * constraintLength;
570 real distMidOutlier =
571 std::sqrt(gmx::square(distPartnerOutlier) - gmx::square(halfDistEqPartner));
572 real distMidCenterInPlane =
573 0.5 * (distMidOutlier - gmx::square(halfDistEqPartner) / distMidOutlier);
574 real distCentralToPlane =
575 std::sqrt(gmx::square(constraintLength) - gmx::square(halfDistEqPartner)
576 - gmx::square(distMidCenterInPlane));
577 real distCogOutlierH = distCentralToPlane / (numPartnerAtoms + 1);
578 real distCogOutlierP =
579 distMidOutlier - (distMidOutlier + distMidCenterInPlane) / (numPartnerAtoms + 1);
580 real distCogOutlier = std::sqrt(gmx::square(distCogOutlierH) + gmx::square(distCogOutlierP));
582 return distCogOutlier;
586 GMX_RELEASE_ASSERT(false, "Only 2 or 3 constraints are supported here");
593 /*! \brief Returns the maximum update group radius for \p moltype */
594 static real computeMaxUpdateGroupRadius(const gmx_moltype_t& moltype,
595 gmx::ArrayRef<const t_iparams> iparams,
596 const RangePartitioning& updateGrouping,
599 GMX_RELEASE_ASSERT(!hasFlexibleConstraints(moltype, iparams),
600 "Flexible constraints are not supported here");
602 const InteractionList& settles = moltype.ilist[F_SETTLE];
604 const ListOfLists<int> at2con = make_at2con(moltype, iparams, FlexibleConstraintTreatment::Include);
606 const auto angleIndices = getAngleIndices(moltype);
609 for (int group = 0; group < updateGrouping.numBlocks(); group++)
611 if (updateGrouping.block(group).size() == 1)
613 /* Single atom group, radius is zero */
617 /* Find the atom maxAtom with the maximum number of constraints */
618 int maxNumConstraints = 0;
620 for (int a : updateGrouping.block(group))
622 const int numConstraints = at2con[a].ssize();
623 if (numConstraints > maxNumConstraints)
625 maxNumConstraints = numConstraints;
629 GMX_ASSERT(maxAtom >= 0 || !settles.empty(),
630 "We should have at least two atoms in the group with constraints");
636 bool allTypesAreEqual = true;
637 int constraintType = -1;
638 real maxConstraintLength = 0;
639 real sumConstraintLengths = 0;
640 bool isFirstConstraint = true;
641 for (const int constraint : at2con[maxAtom])
643 int conIndex = constraint * (1 + NRAL(F_CONSTR));
645 if (conIndex < moltype.ilist[F_CONSTR].size())
647 iparamsIndex = moltype.ilist[F_CONSTR].iatoms[conIndex];
652 moltype.ilist[F_CONSTRNC].iatoms[conIndex - moltype.ilist[F_CONSTR].size()];
654 if (isFirstConstraint)
656 constraintType = iparamsIndex;
657 isFirstConstraint = false;
659 else if (iparamsIndex != constraintType)
661 allTypesAreEqual = false;
663 /* Here we take the maximum constraint length of the A and B
664 * topology, which assumes lambda is between 0 and 1 for
667 real constraintLength =
668 std::max(iparams[iparamsIndex].constr.dA, iparams[iparamsIndex].constr.dB);
669 maxConstraintLength = std::max(maxConstraintLength, constraintLength);
670 sumConstraintLengths += constraintLength;
673 int numConstraints = at2con[maxAtom].ssize();
675 if (numConstraints == 1)
677 /* Single constraint: the radius is the distance from the midpoint */
678 radius = 0.5_real * maxConstraintLength;
684 /* With 2 constraints the maximum possible radius is the
685 * constraint length, so we can use that for the 0 K case.
687 if (numConstraints == 2 && allTypesAreEqual && temperature > 0)
689 radius = constraintGroupRadius<2>(
690 moltype, iparams, maxAtom, at2con, angleIndices, maxConstraintLength, temperature);
692 /* With 3 constraints the maximum possible radius is 1.4 times
693 * the constraint length, so it is worth computing a smaller
694 * radius to enable update groups for systems in a small box.
696 if (numConstraints == 3 && allTypesAreEqual && temperature >= 0)
698 radius = constraintGroupRadius<3>(
699 moltype, iparams, maxAtom, at2con, angleIndices, maxConstraintLength, temperature);
700 if (temperature == 0 && radius >= 0)
702 /* Add a 10% margin for deviation at 0 K */
709 /* Worst case: atom with the longest constraint on one side
710 * of the center, all others on the opposite side
712 radius = maxConstraintLength
713 + (sumConstraintLengths - 2 * maxConstraintLength) / (numConstraints + 1);
716 maxRadius = std::max(maxRadius, radius);
719 for (int i = 0; i < settles.size(); i += 1 + NRAL(F_SETTLE))
721 const real dOH = iparams[settles.iatoms[i]].settle.doh;
722 const real dHH = iparams[settles.iatoms[i]].settle.dhh;
723 /* Compute distance^2 from center of geometry to O and H */
724 const real dCO2 = (4 * dOH * dOH - dHH * dHH) / 9;
725 const real dCH2 = (dOH * dOH + 2 * dHH * dHH) / 9;
726 const real dCAny = std::sqrt(std::max(dCO2, dCH2));
727 maxRadius = std::max(maxRadius, dCAny);
733 real computeMaxUpdateGroupRadius(const gmx_mtop_t& mtop,
734 gmx::ArrayRef<const RangePartitioning> updateGroupingsPerMoleculeType,
737 if (updateGroupingsPerMoleculeType.empty())
742 GMX_RELEASE_ASSERT(updateGroupingsPerMoleculeType.size() == mtop.moltype.size(),
743 "We need one update group entry per moleculetype");
747 for (size_t moltype = 0; moltype < mtop.moltype.size(); moltype++)
749 const real radiusOfThisMoleculeType = computeMaxUpdateGroupRadius(
750 mtop.moltype[moltype], mtop.ffparams.iparams, updateGroupingsPerMoleculeType[moltype], temperature);
751 maxRadius = std::max(maxRadius, radiusOfThisMoleculeType);
757 real computeCutoffMargin(PbcType pbcType, matrix box, const real rlist)
759 return std::sqrt(max_cutoff2(pbcType, box)) - rlist;
762 UpdateGroups::UpdateGroups(std::vector<RangePartitioning>&& updateGroupingPerMoleculeType,
763 const real maxUpdateGroupRadius) :
764 useUpdateGroups_(true),
765 updateGroupingPerMoleculeType_(std::move(updateGroupingPerMoleculeType)),
766 maxUpdateGroupRadius_(maxUpdateGroupRadius)
770 bool systemHasConstraintsOrVsites(const gmx_mtop_t& mtop)
772 IListRange ilistRange(mtop);
773 return std::any_of(ilistRange.begin(), ilistRange.end(), [](const auto& ilists) {
774 return !extractILists(ilists.list(), IF_CONSTRAINT | IF_VSITE).empty();
778 UpdateGroups makeUpdateGroups(const gmx::MDLogger& mdlog,
779 std::vector<RangePartitioning>&& updateGroupingPerMoleculeType,
780 const real maxUpdateGroupRadius,
781 const bool useDomainDecomposition,
782 const bool systemHasConstraintsOrVsites,
783 const real cutoffMargin)
785 MessageStringCollector messages;
787 messages.startContext("When checking whether update groups are usable:");
789 messages.appendIf(!useDomainDecomposition,
790 "Domain decomposition is not active, so there is no need for update groups");
792 messages.appendIf(!systemHasConstraintsOrVsites,
793 "No constraints or virtual sites are in use, so it is best not to use update "
797 updateGroupingPerMoleculeType.empty(),
798 "At least one moleculetype does not conform to the requirements for using update "
802 getenv("GMX_NO_UPDATEGROUPS") != nullptr,
803 "Environment variable GMX_NO_UPDATEGROUPS prohibited the use of update groups");
805 // To use update groups, the large domain-to-domain cutoff
806 // distance should be compatible with the box size.
807 messages.appendIf(2 * maxUpdateGroupRadius >= cutoffMargin,
808 "The combination of rlist and box size prohibits the use of update groups");
810 if (!messages.isEmpty())
812 // Log why we can't use update groups
813 GMX_LOG(mdlog.info).appendText(messages.toString());
814 return UpdateGroups();
818 return UpdateGroups(std::move(updateGroupingPerMoleculeType), maxUpdateGroupRadius);