2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018, 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/mdlib/constr.h"
54 #include "gromacs/mdtypes/inputrec.h"
55 #include "gromacs/topology/block.h"
56 #include "gromacs/topology/idef.h"
57 #include "gromacs/topology/ifunc.h"
58 #include "gromacs/topology/topology.h"
63 /*! \brief Returns whether \p moltype contains flexible constraints */
64 static bool hasFlexibleConstraints(const gmx_moltype_t &moltype,
65 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.data(), 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,
90 gmx::ArrayRef<const t_iparams> iparams)
92 bool hasIncompatibleVsites = false;
94 for (auto &ilist : extractILists(moltype.ilist, IF_VSITE))
96 if (ilist.functionType == F_VSITE2 || ilist.functionType == F_VSITE3)
98 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
100 const t_iparams &iparam = iparams[ilist.iatoms[i]];
103 if (ilist.functionType == F_VSITE2)
105 coeffMin = iparam.vsite.a;
106 coeffSum = iparam.vsite.a;
110 coeffMin = std::min(iparam.vsite.a, iparam.vsite.b);
111 coeffSum = iparam.vsite.a + iparam.vsite.b;
113 if (coeffMin < 0 || coeffSum > 1)
115 hasIncompatibleVsites = true;
122 hasIncompatibleVsites = true;
127 return hasIncompatibleVsites;
130 /*! \brief Returns a merged list with constraints of all types */
131 static InteractionList jointConstraintList(const gmx_moltype_t &moltype)
133 InteractionList ilistCombined;
134 std::vector<int> &iatoms = ilistCombined.iatoms;
136 for (auto &ilist : extractILists(moltype.ilist, IF_CONSTRAINT))
138 if (ilist.functionType == F_SETTLE)
140 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
142 iatoms.push_back(-1);
143 iatoms.push_back(ilist.iatoms[i + 1]);
144 iatoms.push_back(ilist.iatoms[i + 2]);
145 iatoms.push_back(-1);
146 iatoms.push_back(ilist.iatoms[i + 1]);
147 iatoms.push_back(ilist.iatoms[i + 3]);
148 iatoms.push_back(-1);
149 iatoms.push_back(ilist.iatoms[i + 2]);
150 iatoms.push_back(ilist.iatoms[i + 3]);
155 GMX_RELEASE_ASSERT(NRAL(ilist.functionType) == 2, "Can only handle two-atom non-SETTLE constraints");
157 iatoms.insert(iatoms.end(),
158 ilist.iatoms.begin(), ilist.iatoms.end());
162 return ilistCombined;
165 /*! \brief Struct for returning an atom range */
166 struct AtomIndexExtremes
168 int minAtom; //!< The minimum atom index
169 int maxAtom; //!< The maximum atom index
172 /*! \brief Returns the range of constructing atom for vsite with atom index \p a */
173 static AtomIndexExtremes
174 vsiteConstructRange(int a,
175 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
204 constraintAtomRange(int a,
205 const t_blocka &at2con,
206 const InteractionList &ilistConstraints)
208 AtomIndexExtremes extremes = { a, a };
210 for (int i = at2con.index[a]; i < at2con.index[a + 1]; i++)
212 for (int j = 0; j < 2; j++)
214 int atomJ = ilistConstraints.iatoms[at2con.a[i]*3 + 1 + j];
215 extremes.minAtom = std::min(extremes.minAtom, atomJ);
216 extremes.maxAtom = std::max(extremes.maxAtom, atomJ);
223 /*! \brief Returns a list that tells whether atoms in \p moltype are vsites */
224 static std::vector<bool> buildIsParticleVsite(const gmx_moltype_t &moltype)
226 std::vector<bool> isVsite(moltype.atoms.nr);
228 for (auto &ilist : extractILists(moltype.ilist, IF_VSITE))
230 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
232 int vsiteAtom = ilist.iatoms[i + 1];
233 isVsite[vsiteAtom] = true;
240 /*! \brief Returns the size of the update group starting at \p firstAtom or 0 when criteria (see updategroups.h) are not met */
241 static int detectGroup(int firstAtom,
242 const gmx_moltype_t &moltype,
243 const t_blocka &at2con,
244 const InteractionList &ilistConstraints)
246 /* We should be using moltype.atoms.atom[].ptype for checking whether
247 * a particle is a vsite. But the test code can't fill t_atoms,
248 * because it uses C pointers which get double freed.
250 std::vector<bool> isParticleVsite = buildIsParticleVsite(moltype);
252 /* A non-vsite atom without constraints is an update group by itself */
253 if (!isParticleVsite[firstAtom] &&
254 at2con.index[firstAtom + 1] - at2con.index[firstAtom] == 0)
259 /* A (potential) update group starts at firstAtom and should consist
260 * of two or more atoms and possibly vsites. At least one atom should
261 * have constraints with all other atoms and vsites should have all
262 * constructing atoms inside the group. Here we increase lastAtom until
263 * the criteria are fulfilled or exit when criteria cannot be met.
265 int numAtomsWithConstraints = 0;
266 int maxConstraintsPerAtom = 0;
267 int lastAtom = firstAtom;
269 while (a <= lastAtom)
271 if (isParticleVsite[a])
273 AtomIndexExtremes extremes = vsiteConstructRange(a, moltype);
274 if (extremes.minAtom < firstAtom)
276 /* A constructing atom is outside the group,
277 * we can not use update groups.
281 lastAtom = std::max(lastAtom, extremes.maxAtom);
285 int numConstraints = at2con.index[a + 1] - at2con.index[a];
286 if (numConstraints == 0)
288 /* We can not have unconstrained atoms in an update group */
291 /* This atom has at least one constraint.
292 * Check whether all constraints are within the group
293 * and whether we need to extend the group.
295 numAtomsWithConstraints += 1;
296 maxConstraintsPerAtom = std::max(maxConstraintsPerAtom, numConstraints);
298 AtomIndexExtremes extremes =
299 constraintAtomRange(a, at2con, ilistConstraints);
300 if (extremes.minAtom < firstAtom)
302 /* Constraint to atom outside the "group" */
305 lastAtom = std::max(lastAtom, extremes.maxAtom);
311 /* lastAtom might be followed by a vsite that is constructed from atoms
312 * with index <= lastAtom. Handle this case.
314 if (lastAtom + 1 < moltype.atoms.nr &&
315 isParticleVsite[lastAtom + 1])
317 AtomIndexExtremes extremes = vsiteConstructRange(lastAtom + 1, moltype);
318 if (extremes.minAtom < firstAtom)
320 /* Constructing atom precedes the group */
323 else if (extremes.maxAtom <= lastAtom)
325 /* All constructing atoms are in the group, add the vsite to the group */
328 else if (extremes.minAtom <= lastAtom)
330 /* Some, but not all constructing atoms are in the group */
335 GMX_RELEASE_ASSERT(maxConstraintsPerAtom < numAtomsWithConstraints,
336 "We have checked that atoms are only constrained to atoms within the group,"
337 "so each atom should have fewer constraints than the number of atoms");
338 /* Check that at least one atom is constrained to all others */
339 if (maxConstraintsPerAtom != numAtomsWithConstraints - 1)
344 return lastAtom - firstAtom + 1;
347 /*! \brief Returns a list of update groups for \p moltype */
348 static RangePartitioning
349 makeUpdateGroups(const gmx_moltype_t &moltype,
350 gmx::ArrayRef<const t_iparams> iparams)
352 RangePartitioning groups;
354 /* Update groups are not compatible with flexible constraints.
355 * Without dynamics the flexible constraints are ignored,
356 * but since performance for EM/NM is less critical, we do not
357 * use update groups to keep the code here simpler.
359 if (hasFlexibleConstraints(moltype, iparams) ||
360 hasIncompatibleVsites(moltype, iparams))
365 /* Combine all constraint ilists into a single one */
366 InteractionList constraintsCombined = jointConstraintList(moltype);
367 t_ilist ilistsCombined[F_NRE];
368 ilistsCombined[F_CONSTR].nr = constraintsCombined.iatoms.size();
369 ilistsCombined[F_CONSTR].iatoms = constraintsCombined.iatoms.data();
370 ilistsCombined[F_CONSTRNC].nr = 0;
371 /* We "include" flexible constraints, but none are present (checked above) */
372 t_blocka at2con = make_at2con(moltype.atoms.nr,
373 ilistsCombined, iparams.data(),
374 FlexibleConstraintTreatment::Include);
376 bool satisfiesCriteria = true;
379 while (satisfiesCriteria && firstAtom < moltype.atoms.nr)
381 int numAtomsInGroup =
382 detectGroup(firstAtom, moltype, at2con, constraintsCombined);
384 if (numAtomsInGroup == 0)
386 satisfiesCriteria = false;
390 groups.appendBlock(numAtomsInGroup);
392 firstAtom += numAtomsInGroup;
395 if (!satisfiesCriteria)
397 /* Make groups empty, to signal not satisfying the criteria */
401 done_blocka(&at2con);
406 std::vector<RangePartitioning> makeUpdateGroups(const gmx_mtop_t &mtop)
408 std::vector<RangePartitioning> updateGroups;
410 bool systemSatisfiesCriteria = true;
411 for (const gmx_moltype_t &moltype : mtop.moltype)
413 updateGroups.push_back(makeUpdateGroups(moltype, mtop.ffparams.iparams));
415 if (updateGroups.back().numBlocks() == 0)
417 systemSatisfiesCriteria = false;
421 if (!systemSatisfiesCriteria)
423 updateGroups.clear();
429 /*! \brief Returns a map of angles ilist.iatoms indices with the middle atom as key */
430 static std::unordered_multimap<int, int>
431 getAngleIndices(const gmx_moltype_t &moltype)
433 const InteractionList &angles = moltype.ilist[F_ANGLES];
435 std::unordered_multimap<int, int> indices(angles.size());
437 for (int i = 0; i < angles.size(); i += 1 + NRAL(F_ANGLES))
439 indices.insert({ angles.iatoms[i + 2], i });
445 /*! \brief When possible, computes the maximum radius of constrained atom in an update group
447 * Supports groups with 2 or 3 atoms where all partner atoms are connected to
448 * each other by angle potentials. The temperature is used to compute a radius
449 * that is not exceeded with a chance of 10^-9. Note that this computation
450 * assumes there are no other strong forces working on these angular
451 * degrees of freedom.
452 * The return value is -1 when all partners are not connected to each other
453 * by one angle potential, when a potential is perturbed or when an angle
454 * could reach more than 180 degrees.
456 template <int numPartnerAtoms> static real
457 constraintGroupRadius(const gmx_moltype_t &moltype,
458 gmx::ArrayRef<const t_iparams> iparams,
459 const int centralAtom,
460 const t_blocka &at2con,
461 const std::unordered_multimap<int, int> &angleIndices,
462 const real constraintLength,
463 const real temperature)
465 const int numConstraints = at2con.index[centralAtom + 1] - at2con.index[centralAtom];
466 GMX_RELEASE_ASSERT(numConstraints == numPartnerAtoms, "We expect as many constraints as partner atoms here");
468 std::array<int, numPartnerAtoms> partnerAtoms;
469 for (int i = 0; i < numPartnerAtoms; i++)
471 const int ind = at2con.a[at2con.index[centralAtom] + i]*3;
472 if (ind >= moltype.ilist[F_CONSTR].size())
474 /* This is a flexible constraint, we don't optimize for that */
477 const int a1 = moltype.ilist[F_CONSTR].iatoms[ind + 1];
478 const int a2 = moltype.ilist[F_CONSTR].iatoms[ind + 2];
479 partnerAtoms[i] = (a1 == centralAtom ? a2 : a1);
482 const InteractionList &angles = moltype.ilist[F_ANGLES];
483 auto range = angleIndices.equal_range(centralAtom);
485 std::array<int, numPartnerAtoms> numAngles = { 0 };
486 bool areSameType = true;
487 for (auto it = range.first; it != range.second; ++it)
489 /* Check if the outer atoms in the angle are both partner atoms */
490 int numAtomsFound = 0;
491 for (int ind = it->second + 1; ind < it->second + 4; ind += 2)
493 for (const int &partnerAtom : partnerAtoms)
495 if (angles.iatoms[ind] == partnerAtom)
501 if (numAtomsFound == 2)
503 /* Check that the angle potentials have the same type */
506 angleType = angles.iatoms[it->second];
508 else if (angles.iatoms[it->second] != angleType)
512 /* Count the number of angle interactions per atoms */
513 for (int ind = it->second + 1; ind < it->second + 4; ind += 2)
515 for (size_t i = 0; i < partnerAtoms.size(); i++)
517 if (angles.iatoms[ind] == partnerAtoms[i])
526 bool criteriaSatisfied = areSameType;
527 for (int i = 0; i < numPartnerAtoms; i++)
529 if (numAngles[i] != numPartnerAtoms - 1)
531 criteriaSatisfied = false;
535 /* We don't bother optimizing the perturbed angle case */
536 const t_iparams &angleParams = iparams[angleType];
537 if (criteriaSatisfied &&
538 angleParams.harmonic.rB == angleParams.harmonic.rA &&
539 angleParams.harmonic.krB == angleParams.harmonic.krA)
541 /* Set number of stddevs such that change of exceeding < 10^-9 */
542 constexpr real c_numSigma = 6.0;
543 /* Compute the maximally stretched angle */
544 const real eqAngle = angleParams.harmonic.rA*DEG2RAD;
545 const real fc = angleParams.harmonic.krA;
546 const real maxAngle = eqAngle + c_numSigma*BOLTZ*temperature/((numPartnerAtoms - 1)*fc);
547 if (maxAngle >= M_PI)
552 if (numPartnerAtoms == 2)
554 /* With two atoms constrainted to a cental atom we have a triangle
555 * with two equal sides because the constraint type is equal.
556 * Return the distance from the COG to the farthest two corners,
557 * i.e. the partner atoms.
559 real distMidPartner = std::sin(0.5*maxAngle)*constraintLength;
560 real distCentralMid = std::cos(0.5*maxAngle)*constraintLength;
561 real distCogMid = distCentralMid*numPartnerAtoms/(numPartnerAtoms + 1);
562 real distCogPartner = std::sqrt(gmx::square(distMidPartner) + gmx::square(distCogMid));
564 return distCogPartner;
566 else if (numPartnerAtoms == 3)
568 /* With three atoms constrainted to a cental atom we have the
569 * largest COG-atom distance when one partner atom (the outlier)
570 * moves out by stretching its two angles by an equal amount.
571 * The math here gets a bit more involved, but it is still
572 * rather straightforward geometry.
573 * We first compute distances in the plane of the three partners.
574 * Here we have two "equilibrium" partners and one outlier.
575 * We make use of the "Mid" point between the two "Eq" partners.
576 * We project the central atom on this plane.
577 * Then we compute the distance of the central atom to the plane.
578 * The distance of the COG to the ourlier is returned.
580 real halfDistEqPartner = std::sin(0.5*eqAngle)*constraintLength;
581 real distPartnerOutlier = 2*std::sin(0.5*maxAngle)*constraintLength;
582 real distMidOutlier = std::sqrt(gmx::square(distPartnerOutlier) - gmx::square(halfDistEqPartner));
583 real distMidCenterInPlane = 0.5*(distMidOutlier - gmx::square(halfDistEqPartner)/distMidOutlier);
584 real distCentralToPlane = std::sqrt(gmx::square(constraintLength) - gmx::square(halfDistEqPartner) - gmx::square(distMidCenterInPlane));
585 real distCogOutlierH = distCentralToPlane/(numPartnerAtoms + 1);
586 real distCogOutlierP = distMidOutlier - (distMidOutlier + distMidCenterInPlane)/(numPartnerAtoms + 1);
587 real distCogOutlier = std::sqrt(gmx::square(distCogOutlierH) + gmx::square(distCogOutlierP));
589 return distCogOutlier;
593 GMX_RELEASE_ASSERT(false, "Only 2 or 3 constraints are supported here");
600 /*! \brief Returns the maximum update group radius for \p moltype */
602 computeMaxUpdateGroupRadius(const gmx_moltype_t &moltype,
603 gmx::ArrayRef<const t_iparams> iparams,
604 const RangePartitioning &updateGroups,
607 GMX_RELEASE_ASSERT(!hasFlexibleConstraints(moltype, iparams),
608 "Flexible constraints are not supported here");
610 const InteractionList &settles = moltype.ilist[F_SETTLE];
613 make_at2con(moltype, iparams, FlexibleConstraintTreatment::Include);
615 const auto angleIndices = getAngleIndices(moltype);
618 for (int group = 0; group < updateGroups.numBlocks(); group++)
620 if (updateGroups.block(group).size() == 1)
622 /* Single atom group, radius is zero */
626 /* Find the atom maxAtom with the maximum number of constraints */
627 int maxNumConstraints = 0;
629 for (int a : updateGroups.block(group))
631 int numConstraints = at2con.index[a + 1] - at2con.index[a];
632 if (numConstraints > maxNumConstraints)
634 maxNumConstraints = numConstraints;
638 GMX_ASSERT(maxAtom >= 0 || settles.size() > 0,
639 "We should have at least two atoms in the group with constraints");
645 bool allTypesAreEqual = true;
646 int constraintType = -1;
647 real maxConstraintLength = 0;
648 real sumConstraintLengths = 0;
649 for (int i = at2con.index[maxAtom]; i < at2con.index[maxAtom + 1]; i++)
651 int conIndex = at2con.a[i]*(1 + NRAL(F_CONSTR));
653 if (conIndex < moltype.ilist[F_CONSTR].size())
655 iparamsIndex = moltype.ilist[F_CONSTR].iatoms[conIndex];
659 iparamsIndex = moltype.ilist[F_CONSTRNC].iatoms[conIndex - moltype.ilist[F_CONSTR].size()];
661 if (i == at2con.index[maxAtom])
663 constraintType = iparamsIndex;
665 else if (iparamsIndex != constraintType)
667 allTypesAreEqual = false;
669 /* Here we take the maximum constraint length of the A and B
670 * topology, which assumes lambda is between 0 and 1 for
673 real constraintLength = std::max(iparams[iparamsIndex].constr.dA,
674 iparams[iparamsIndex].constr.dB);
675 maxConstraintLength = std::max(maxConstraintLength,
677 sumConstraintLengths += constraintLength;
680 int numConstraints = at2con.index[maxAtom + 1] - at2con.index[maxAtom];
682 if (numConstraints == 1)
684 /* Single constraint: the radius is the distance from the midpoint */
685 radius = 0.5_real*maxConstraintLength;
691 /* With 2 constraints the maximum possible radius is the
692 * constraint length, so we can use that for the 0 K case.
694 if (numConstraints == 2 && allTypesAreEqual && temperature > 0)
696 radius = constraintGroupRadius<2>(moltype, iparams,
702 /* With 3 constraints the maximum possible radius is 1.4 times
703 * the constraint length, so it is worth computing a smaller
704 * radius to enable update groups for systems in a small box.
706 if (numConstraints == 3 && allTypesAreEqual && temperature >= 0)
708 radius = constraintGroupRadius<3>(moltype, iparams,
713 if (temperature == 0 && radius >= 0)
715 /* Add a 10% margin for deviation at 0 K */
722 /* Worst case: atom with the longest constraint on one side
723 * of the center, all others on the opposite side
725 radius = maxConstraintLength + (sumConstraintLengths - 2*maxConstraintLength)/(numConstraints + 1);
728 maxRadius = std::max(maxRadius, radius);
731 for (int i = 0; i < settles.size(); i += 1 + NRAL(F_SETTLE))
733 const real dOH = iparams[settles.iatoms[i]].settle.doh;
734 const real dHH = iparams[settles.iatoms[i]].settle.dhh;
735 /* Compute distance^2 from center of geometry to O and H */
736 const real dCO2 = (4*dOH*dOH - dHH*dHH)/9;
737 const real dCH2 = (dOH*dOH + 2*dHH*dHH)/9;
738 const real dCAny = std::sqrt(std::max(dCO2, dCH2));
739 maxRadius = std::max(maxRadius, dCAny);
742 done_blocka(&at2con);
748 computeMaxUpdateGroupRadius(const gmx_mtop_t &mtop,
749 gmx::ArrayRef<const RangePartitioning> updateGroups,
752 if (updateGroups.empty())
757 GMX_RELEASE_ASSERT(static_cast<size_t>(updateGroups.size()) == mtop.moltype.size(),
758 "We need one update group entry per moleculetype");
762 for (size_t moltype = 0; moltype < mtop.moltype.size(); moltype++)
765 = std::max(maxRadius,
766 computeMaxUpdateGroupRadius(mtop.moltype[moltype],
767 mtop.ffparams.iparams,
768 updateGroups[moltype],