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/mdtypes/inputrec.h"
56 #include "gromacs/topology/idef.h"
57 #include "gromacs/topology/ifunc.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/listoflists.h"
64 /*! \brief Returns whether \p moltype contains flexible constraints */
65 static bool hasFlexibleConstraints(const gmx_moltype_t& moltype, gmx::ArrayRef<const t_iparams> iparams)
67 for (auto& ilist : extractILists(moltype.ilist, IF_CONSTRAINT))
69 if (ilist.functionType != F_SETTLE)
71 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
73 if (isConstraintFlexible(iparams, ilist.iatoms[i]))
84 /*! \brief Returns whether moltype has incompatible vsites.
86 * For simplicity the only compatible vsites are linear 2 or 3 atom sites
87 * that are constructed in between the 2 or 3 contructing atoms,
89 static bool hasIncompatibleVsites(const gmx_moltype_t& moltype, gmx::ArrayRef<const t_iparams> iparams)
91 bool hasIncompatibleVsites = false;
93 for (auto& ilist : extractILists(moltype.ilist, IF_VSITE))
95 if (ilist.functionType == F_VSITE2 || ilist.functionType == F_VSITE3)
97 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
99 const t_iparams& iparam = iparams[ilist.iatoms[i]];
102 if (ilist.functionType == F_VSITE2)
104 coeffMin = iparam.vsite.a;
105 coeffSum = iparam.vsite.a;
109 coeffMin = std::min(iparam.vsite.a, iparam.vsite.b);
110 coeffSum = iparam.vsite.a + iparam.vsite.b;
112 if (coeffMin < 0 || coeffSum > 1)
114 hasIncompatibleVsites = true;
121 hasIncompatibleVsites = true;
126 return hasIncompatibleVsites;
129 /*! \brief Returns a merged list with constraints of all types */
130 static InteractionList jointConstraintList(const gmx_moltype_t& moltype)
132 InteractionList ilistCombined;
133 std::vector<int>& iatoms = ilistCombined.iatoms;
135 for (auto& ilist : extractILists(moltype.ilist, IF_CONSTRAINT))
137 if (ilist.functionType == F_SETTLE)
139 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
141 iatoms.push_back(-1);
142 iatoms.push_back(ilist.iatoms[i + 1]);
143 iatoms.push_back(ilist.iatoms[i + 2]);
144 iatoms.push_back(-1);
145 iatoms.push_back(ilist.iatoms[i + 1]);
146 iatoms.push_back(ilist.iatoms[i + 3]);
147 iatoms.push_back(-1);
148 iatoms.push_back(ilist.iatoms[i + 2]);
149 iatoms.push_back(ilist.iatoms[i + 3]);
154 GMX_RELEASE_ASSERT(NRAL(ilist.functionType) == 2,
155 "Can only handle two-atom non-SETTLE constraints");
157 iatoms.insert(iatoms.end(), ilist.iatoms.begin(), ilist.iatoms.end());
161 return ilistCombined;
164 /*! \brief Struct for returning an atom range */
165 struct AtomIndexExtremes
167 int minAtom; //!< The minimum atom index
168 int maxAtom; //!< The maximum atom index
171 /*! \brief Returns the range of constructing atom for vsite with atom index \p a */
172 static AtomIndexExtremes vsiteConstructRange(int a, const gmx_moltype_t& moltype)
174 AtomIndexExtremes extremes = { -1, -1 };
176 for (auto& ilist : extractILists(moltype.ilist, IF_VSITE))
178 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
180 if (ilist.iatoms[i + 1] == a)
182 extremes.minAtom = ilist.iatoms[i + 2];
183 extremes.maxAtom = ilist.iatoms[i + 2];
184 for (size_t j = i + 3; j < i + ilistStride(ilist); j++)
186 extremes.minAtom = std::min(extremes.minAtom, ilist.iatoms[j]);
187 extremes.maxAtom = std::max(extremes.maxAtom, ilist.iatoms[j]);
194 GMX_RELEASE_ASSERT(false, "If a is a vsite, we should have found constructing atoms");
199 /*! \brief Returns the range of atoms constrained to atom \p a (including \p a itself) */
200 static AtomIndexExtremes constraintAtomRange(int a,
201 const ListOfLists<int>& at2con,
202 const InteractionList& ilistConstraints)
204 AtomIndexExtremes extremes = { a, a };
206 for (const int constraint : at2con[a])
208 for (int j = 0; j < 2; j++)
210 int atomJ = ilistConstraints.iatoms[constraint * 3 + 1 + j];
211 extremes.minAtom = std::min(extremes.minAtom, atomJ);
212 extremes.maxAtom = std::max(extremes.maxAtom, atomJ);
219 /*! \brief Returns a list that tells whether atoms in \p moltype are vsites */
220 static std::vector<bool> buildIsParticleVsite(const gmx_moltype_t& moltype)
222 std::vector<bool> isVsite(moltype.atoms.nr);
224 for (auto& ilist : extractILists(moltype.ilist, IF_VSITE))
226 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
228 int vsiteAtom = ilist.iatoms[i + 1];
229 isVsite[vsiteAtom] = true;
236 /*! \brief Returns the size of the update group starting at \p firstAtom or 0 when criteria (see updategroups.h) are not met */
237 static int detectGroup(int firstAtom,
238 const gmx_moltype_t& moltype,
239 const ListOfLists<int>& at2con,
240 const InteractionList& ilistConstraints)
242 /* We should be using moltype.atoms.atom[].ptype for checking whether
243 * a particle is a vsite. But the test code can't fill t_atoms,
244 * because it uses C pointers which get double freed.
246 std::vector<bool> isParticleVsite = buildIsParticleVsite(moltype);
248 /* A non-vsite atom without constraints is an update group by itself */
249 if (!isParticleVsite[firstAtom] && at2con[firstAtom].empty())
254 /* A (potential) update group starts at firstAtom and should consist
255 * of two or more atoms and possibly vsites. At least one atom should
256 * have constraints with all other atoms and vsites should have all
257 * constructing atoms inside the group. Here we increase lastAtom until
258 * the criteria are fulfilled or exit when criteria cannot be met.
260 int numAtomsWithConstraints = 0;
261 int maxConstraintsPerAtom = 0;
262 int lastAtom = firstAtom;
264 while (a <= lastAtom)
266 if (isParticleVsite[a])
268 AtomIndexExtremes extremes = vsiteConstructRange(a, moltype);
269 if (extremes.minAtom < firstAtom)
271 /* A constructing atom is outside the group,
272 * we can not use update groups.
276 lastAtom = std::max(lastAtom, extremes.maxAtom);
280 const int numConstraints = at2con[a].ssize();
281 if (numConstraints == 0)
283 /* We can not have unconstrained atoms in an update group */
286 /* This atom has at least one constraint.
287 * Check whether all constraints are within the group
288 * and whether we need to extend the group.
290 numAtomsWithConstraints += 1;
291 maxConstraintsPerAtom = std::max(maxConstraintsPerAtom, numConstraints);
293 AtomIndexExtremes extremes = constraintAtomRange(a, at2con, ilistConstraints);
294 if (extremes.minAtom < firstAtom)
296 /* Constraint to atom outside the "group" */
299 lastAtom = std::max(lastAtom, extremes.maxAtom);
305 /* lastAtom might be followed by a vsite that is constructed from atoms
306 * with index <= lastAtom. Handle this case.
308 if (lastAtom + 1 < moltype.atoms.nr && isParticleVsite[lastAtom + 1])
310 AtomIndexExtremes extremes = vsiteConstructRange(lastAtom + 1, moltype);
311 if (extremes.minAtom < firstAtom)
312 { // NOLINT bugprone-branch-clone
313 /* Constructing atom precedes the group */
316 else if (extremes.maxAtom <= lastAtom)
318 /* All constructing atoms are in the group, add the vsite to the group */
321 else if (extremes.minAtom <= lastAtom)
323 /* Some, but not all constructing atoms are in the group */
328 GMX_RELEASE_ASSERT(maxConstraintsPerAtom < numAtomsWithConstraints,
329 "We have checked that atoms are only constrained to atoms within the group,"
330 "so each atom should have fewer constraints than the number of atoms");
331 /* Check that at least one atom is constrained to all others */
332 if (maxConstraintsPerAtom != numAtomsWithConstraints - 1)
337 return lastAtom - firstAtom + 1;
340 /*! \brief Returns a list of update groups for \p moltype */
341 static RangePartitioning makeUpdateGroupingsPerMoleculeType(const gmx_moltype_t& moltype,
342 gmx::ArrayRef<const t_iparams> iparams)
344 RangePartitioning groups;
346 /* Update groups are not compatible with flexible constraints.
347 * Without dynamics the flexible constraints are ignored,
348 * but since performance for EM/NM is less critical, we do not
349 * use update groups to keep the code here simpler.
351 if (hasFlexibleConstraints(moltype, iparams) || hasIncompatibleVsites(moltype, iparams))
356 /* Combine all constraint ilists into a single one */
357 std::array<InteractionList, F_NRE> ilistsCombined;
358 ilistsCombined[F_CONSTR] = jointConstraintList(moltype);
359 /* We "include" flexible constraints, but none are present (checked above) */
360 const ListOfLists<int> at2con = make_at2con(
361 moltype.atoms.nr, ilistsCombined, iparams, FlexibleConstraintTreatment::Include);
363 bool satisfiesCriteria = true;
366 while (satisfiesCriteria && firstAtom < moltype.atoms.nr)
368 int numAtomsInGroup = detectGroup(firstAtom, moltype, at2con, ilistsCombined[F_CONSTR]);
370 if (numAtomsInGroup == 0)
372 satisfiesCriteria = false;
376 groups.appendBlock(numAtomsInGroup);
378 firstAtom += numAtomsInGroup;
381 if (!satisfiesCriteria)
383 /* Make groups empty, to signal not satisfying the criteria */
390 std::vector<RangePartitioning> makeUpdateGroupingsPerMoleculeType(const gmx_mtop_t& mtop)
392 std::vector<RangePartitioning> updateGroupingsPerMoleculeType;
394 bool systemSatisfiesCriteria = true;
395 for (const gmx_moltype_t& moltype : mtop.moltype)
397 updateGroupingsPerMoleculeType.push_back(
398 makeUpdateGroupingsPerMoleculeType(moltype, mtop.ffparams.iparams));
400 if (updateGroupingsPerMoleculeType.back().numBlocks() == 0)
402 systemSatisfiesCriteria = false;
406 if (!systemSatisfiesCriteria)
408 updateGroupingsPerMoleculeType.clear();
411 return updateGroupingsPerMoleculeType;
414 /*! \brief Returns a map of angles ilist.iatoms indices with the middle atom as key */
415 static std::unordered_multimap<int, int> getAngleIndices(const gmx_moltype_t& moltype)
417 const InteractionList& angles = moltype.ilist[F_ANGLES];
419 std::unordered_multimap<int, int> indices(angles.size());
421 for (int i = 0; i < angles.size(); i += 1 + NRAL(F_ANGLES))
423 indices.insert({ angles.iatoms[i + 2], i });
429 /*! \brief When possible, computes the maximum radius of constrained atom in an update group
431 * Supports groups with 2 or 3 atoms where all partner atoms are connected to
432 * each other by angle potentials. The temperature is used to compute a radius
433 * that is not exceeded with a chance of 10^-9. Note that this computation
434 * assumes there are no other strong forces working on these angular
435 * degrees of freedom.
436 * The return value is -1 when all partners are not connected to each other
437 * by one angle potential, when a potential is perturbed or when an angle
438 * could reach more than 180 degrees.
440 template<int numPartnerAtoms>
441 static real constraintGroupRadius(const gmx_moltype_t& moltype,
442 gmx::ArrayRef<const t_iparams> iparams,
443 const int centralAtom,
444 const ListOfLists<int>& at2con,
445 const std::unordered_multimap<int, int>& angleIndices,
446 const real constraintLength,
447 const real temperature)
449 const int numConstraints = at2con[centralAtom].ssize();
450 GMX_RELEASE_ASSERT(numConstraints == numPartnerAtoms,
451 "We expect as many constraints as partner atoms here");
453 std::array<int, numPartnerAtoms> partnerAtoms;
454 for (int i = 0; i < numPartnerAtoms; i++)
456 const int ind = at2con[centralAtom][i] * 3;
457 if (ind >= moltype.ilist[F_CONSTR].size())
459 /* This is a flexible constraint, we don't optimize for that */
462 const int a1 = moltype.ilist[F_CONSTR].iatoms[ind + 1];
463 const int a2 = moltype.ilist[F_CONSTR].iatoms[ind + 2];
464 partnerAtoms[i] = (a1 == centralAtom ? a2 : a1);
467 const InteractionList& angles = moltype.ilist[F_ANGLES];
468 auto range = angleIndices.equal_range(centralAtom);
470 std::array<int, numPartnerAtoms> numAngles = { 0 };
471 bool areSameType = true;
472 for (auto it = range.first; it != range.second; ++it)
474 /* Check if the outer atoms in the angle are both partner atoms */
475 int numAtomsFound = 0;
476 for (int ind = it->second + 1; ind < it->second + 4; ind += 2)
478 for (const int& partnerAtom : partnerAtoms)
480 if (angles.iatoms[ind] == partnerAtom)
486 if (numAtomsFound == 2)
488 /* Check that the angle potentials have the same type */
491 angleType = angles.iatoms[it->second];
493 else if (angles.iatoms[it->second] != angleType)
497 /* Count the number of angle interactions per atoms */
498 for (int ind = it->second + 1; ind < it->second + 4; ind += 2)
500 for (size_t i = 0; i < partnerAtoms.size(); i++)
502 if (angles.iatoms[ind] == partnerAtoms[i])
511 bool criteriaSatisfied = areSameType;
512 for (int i = 0; i < numPartnerAtoms; i++)
514 if (numAngles[i] != numPartnerAtoms - 1)
516 criteriaSatisfied = false;
520 /* We don't bother optimizing the perturbed angle case */
521 const t_iparams& angleParams = iparams[angleType];
522 if (criteriaSatisfied && angleParams.harmonic.rB == angleParams.harmonic.rA
523 && angleParams.harmonic.krB == angleParams.harmonic.krA)
525 /* Set number of stddevs such that change of exceeding < 10^-9 */
526 constexpr real c_numSigma = 6.0;
527 /* Compute the maximally stretched angle */
528 const real eqAngle = angleParams.harmonic.rA * gmx::c_deg2Rad;
529 const real fc = angleParams.harmonic.krA;
530 const real maxAngle =
531 eqAngle + c_numSigma * gmx::c_boltz * temperature / ((numPartnerAtoms - 1) * fc);
532 if (maxAngle >= M_PI)
537 if (numPartnerAtoms == 2)
539 /* With two atoms constrainted to a cental atom we have a triangle
540 * with two equal sides because the constraint type is equal.
541 * Return the distance from the COG to the farthest two corners,
542 * i.e. the partner atoms.
544 real distMidPartner = std::sin(0.5 * maxAngle) * constraintLength;
545 real distCentralMid = std::cos(0.5 * maxAngle) * constraintLength;
546 real distCogMid = distCentralMid * numPartnerAtoms / (numPartnerAtoms + 1);
547 real distCogPartner = std::sqrt(gmx::square(distMidPartner) + gmx::square(distCogMid));
549 return distCogPartner;
551 else if (numPartnerAtoms == 3)
553 /* With three atoms constrainted to a cental atom we have the
554 * largest COG-atom distance when one partner atom (the outlier)
555 * moves out by stretching its two angles by an equal amount.
556 * The math here gets a bit more involved, but it is still
557 * rather straightforward geometry.
558 * We first compute distances in the plane of the three partners.
559 * Here we have two "equilibrium" partners and one outlier.
560 * We make use of the "Mid" point between the two "Eq" partners.
561 * We project the central atom on this plane.
562 * Then we compute the distance of the central atom to the plane.
563 * The distance of the COG to the ourlier is returned.
565 real halfDistEqPartner = std::sin(0.5 * eqAngle) * constraintLength;
566 real distPartnerOutlier = 2 * std::sin(0.5 * maxAngle) * constraintLength;
567 real distMidOutlier =
568 std::sqrt(gmx::square(distPartnerOutlier) - gmx::square(halfDistEqPartner));
569 real distMidCenterInPlane =
570 0.5 * (distMidOutlier - gmx::square(halfDistEqPartner) / distMidOutlier);
571 real distCentralToPlane =
572 std::sqrt(gmx::square(constraintLength) - gmx::square(halfDistEqPartner)
573 - gmx::square(distMidCenterInPlane));
574 real distCogOutlierH = distCentralToPlane / (numPartnerAtoms + 1);
575 real distCogOutlierP =
576 distMidOutlier - (distMidOutlier + distMidCenterInPlane) / (numPartnerAtoms + 1);
577 real distCogOutlier = std::sqrt(gmx::square(distCogOutlierH) + gmx::square(distCogOutlierP));
579 return distCogOutlier;
583 GMX_RELEASE_ASSERT(false, "Only 2 or 3 constraints are supported here");
590 /*! \brief Returns the maximum update group radius for \p moltype */
591 static real computeMaxUpdateGroupRadius(const gmx_moltype_t& moltype,
592 gmx::ArrayRef<const t_iparams> iparams,
593 const RangePartitioning& updateGrouping,
596 GMX_RELEASE_ASSERT(!hasFlexibleConstraints(moltype, iparams),
597 "Flexible constraints are not supported here");
599 const InteractionList& settles = moltype.ilist[F_SETTLE];
601 const ListOfLists<int> at2con = make_at2con(moltype, iparams, FlexibleConstraintTreatment::Include);
603 const auto angleIndices = getAngleIndices(moltype);
606 for (int group = 0; group < updateGrouping.numBlocks(); group++)
608 if (updateGrouping.block(group).size() == 1)
610 /* Single atom group, radius is zero */
614 /* Find the atom maxAtom with the maximum number of constraints */
615 int maxNumConstraints = 0;
617 for (int a : updateGrouping.block(group))
619 const int numConstraints = at2con[a].ssize();
620 if (numConstraints > maxNumConstraints)
622 maxNumConstraints = numConstraints;
626 GMX_ASSERT(maxAtom >= 0 || !settles.empty(),
627 "We should have at least two atoms in the group with constraints");
633 bool allTypesAreEqual = true;
634 int constraintType = -1;
635 real maxConstraintLength = 0;
636 real sumConstraintLengths = 0;
637 bool isFirstConstraint = true;
638 for (const int constraint : at2con[maxAtom])
640 int conIndex = constraint * (1 + NRAL(F_CONSTR));
642 if (conIndex < moltype.ilist[F_CONSTR].size())
644 iparamsIndex = moltype.ilist[F_CONSTR].iatoms[conIndex];
649 moltype.ilist[F_CONSTRNC].iatoms[conIndex - moltype.ilist[F_CONSTR].size()];
651 if (isFirstConstraint)
653 constraintType = iparamsIndex;
654 isFirstConstraint = false;
656 else if (iparamsIndex != constraintType)
658 allTypesAreEqual = false;
660 /* Here we take the maximum constraint length of the A and B
661 * topology, which assumes lambda is between 0 and 1 for
664 real constraintLength =
665 std::max(iparams[iparamsIndex].constr.dA, iparams[iparamsIndex].constr.dB);
666 maxConstraintLength = std::max(maxConstraintLength, constraintLength);
667 sumConstraintLengths += constraintLength;
670 int numConstraints = at2con[maxAtom].ssize();
672 if (numConstraints == 1)
674 /* Single constraint: the radius is the distance from the midpoint */
675 radius = 0.5_real * maxConstraintLength;
681 /* With 2 constraints the maximum possible radius is the
682 * constraint length, so we can use that for the 0 K case.
684 if (numConstraints == 2 && allTypesAreEqual && temperature > 0)
686 radius = constraintGroupRadius<2>(
687 moltype, iparams, maxAtom, at2con, angleIndices, maxConstraintLength, temperature);
689 /* With 3 constraints the maximum possible radius is 1.4 times
690 * the constraint length, so it is worth computing a smaller
691 * radius to enable update groups for systems in a small box.
693 if (numConstraints == 3 && allTypesAreEqual && temperature >= 0)
695 radius = constraintGroupRadius<3>(
696 moltype, iparams, maxAtom, at2con, angleIndices, maxConstraintLength, temperature);
697 if (temperature == 0 && radius >= 0)
699 /* Add a 10% margin for deviation at 0 K */
706 /* Worst case: atom with the longest constraint on one side
707 * of the center, all others on the opposite side
709 radius = maxConstraintLength
710 + (sumConstraintLengths - 2 * maxConstraintLength) / (numConstraints + 1);
713 maxRadius = std::max(maxRadius, radius);
716 for (int i = 0; i < settles.size(); i += 1 + NRAL(F_SETTLE))
718 const real dOH = iparams[settles.iatoms[i]].settle.doh;
719 const real dHH = iparams[settles.iatoms[i]].settle.dhh;
720 /* Compute distance^2 from center of geometry to O and H */
721 const real dCO2 = (4 * dOH * dOH - dHH * dHH) / 9;
722 const real dCH2 = (dOH * dOH + 2 * dHH * dHH) / 9;
723 const real dCAny = std::sqrt(std::max(dCO2, dCH2));
724 maxRadius = std::max(maxRadius, dCAny);
730 real computeMaxUpdateGroupRadius(const gmx_mtop_t& mtop,
731 gmx::ArrayRef<const RangePartitioning> updateGroupingsPerMoleculeType,
734 if (updateGroupingsPerMoleculeType.empty())
739 GMX_RELEASE_ASSERT(updateGroupingsPerMoleculeType.size() == mtop.moltype.size(),
740 "We need one update group entry per moleculetype");
744 for (size_t moltype = 0; moltype < mtop.moltype.size(); moltype++)
746 const real radiusOfThisMoleculeType = computeMaxUpdateGroupRadius(
747 mtop.moltype[moltype], mtop.ffparams.iparams, updateGroupingsPerMoleculeType[moltype], temperature);
748 maxRadius = std::max(maxRadius, radiusOfThisMoleculeType);