cf4856825c432c8599f5847f3997e62cb00fae5e
[alexxy/gromacs.git] / src / gromacs / mdlib / updategroups.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019,2020, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /* \internal \file
36  *
37  * \brief Implements functions for generating update groups
38  *
39  * \author Berk Hess <hess@kth.se>
40  * \ingroup module_mdlib
41  */
42
43 #include "gmxpre.h"
44
45 #include "updategroups.h"
46
47 #include <cmath>
48
49 #include <unordered_map>
50
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/idef.h"
56 #include "gromacs/topology/ifunc.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/utility/listoflists.h"
59
60 namespace gmx
61 {
62
63 /*! \brief Returns whether \p moltype contains flexible constraints */
64 static bool hasFlexibleConstraints(const gmx_moltype_t& moltype, gmx::ArrayRef<const t_iparams> iparams)
65 {
66     for (auto& ilist : extractILists(moltype.ilist, IF_CONSTRAINT))
67     {
68         if (ilist.functionType != F_SETTLE)
69         {
70             for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
71             {
72                 if (isConstraintFlexible(iparams, ilist.iatoms[i]))
73                 {
74                     return true;
75                 }
76             }
77         }
78     }
79
80     return false;
81 }
82
83 /*! \brief Returns whether moltype has incompatible vsites.
84  *
85  * For simplicity the only compatible vsites are linear 2 or 3 atom sites
86  * that are constructed in between the 2 or 3 contructing atoms,
87  */
88 static bool hasIncompatibleVsites(const gmx_moltype_t& moltype, gmx::ArrayRef<const t_iparams> iparams)
89 {
90     bool hasIncompatibleVsites = false;
91
92     for (auto& ilist : extractILists(moltype.ilist, IF_VSITE))
93     {
94         if (ilist.functionType == F_VSITE2 || ilist.functionType == F_VSITE3)
95         {
96             for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
97             {
98                 const t_iparams& iparam = iparams[ilist.iatoms[i]];
99                 real             coeffMin;
100                 real             coeffSum;
101                 if (ilist.functionType == F_VSITE2)
102                 {
103                     coeffMin = iparam.vsite.a;
104                     coeffSum = iparam.vsite.a;
105                 }
106                 else
107                 {
108                     coeffMin = std::min(iparam.vsite.a, iparam.vsite.b);
109                     coeffSum = iparam.vsite.a + iparam.vsite.b;
110                 }
111                 if (coeffMin < 0 || coeffSum > 1)
112                 {
113                     hasIncompatibleVsites = true;
114                     break;
115                 }
116             }
117         }
118         else
119         {
120             hasIncompatibleVsites = true;
121             break;
122         }
123     }
124
125     return hasIncompatibleVsites;
126 }
127
128 /*! \brief Returns a merged list with constraints of all types */
129 static InteractionList jointConstraintList(const gmx_moltype_t& moltype)
130 {
131     InteractionList   ilistCombined;
132     std::vector<int>& iatoms = ilistCombined.iatoms;
133
134     for (auto& ilist : extractILists(moltype.ilist, IF_CONSTRAINT))
135     {
136         if (ilist.functionType == F_SETTLE)
137         {
138             for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
139             {
140                 iatoms.push_back(-1);
141                 iatoms.push_back(ilist.iatoms[i + 1]);
142                 iatoms.push_back(ilist.iatoms[i + 2]);
143                 iatoms.push_back(-1);
144                 iatoms.push_back(ilist.iatoms[i + 1]);
145                 iatoms.push_back(ilist.iatoms[i + 3]);
146                 iatoms.push_back(-1);
147                 iatoms.push_back(ilist.iatoms[i + 2]);
148                 iatoms.push_back(ilist.iatoms[i + 3]);
149             }
150         }
151         else
152         {
153             GMX_RELEASE_ASSERT(NRAL(ilist.functionType) == 2,
154                                "Can only handle two-atom non-SETTLE constraints");
155
156             iatoms.insert(iatoms.end(), ilist.iatoms.begin(), ilist.iatoms.end());
157         }
158     }
159
160     return ilistCombined;
161 }
162
163 /*! \brief Struct for returning an atom range */
164 struct AtomIndexExtremes
165 {
166     int minAtom; //!< The minimum atom index
167     int maxAtom; //!< The maximum atom index
168 };
169
170 /*! \brief Returns the range of constructing atom for vsite with atom index \p a */
171 static AtomIndexExtremes vsiteConstructRange(int a, const gmx_moltype_t& moltype)
172 {
173     AtomIndexExtremes extremes = { -1, -1 };
174
175     for (auto& ilist : extractILists(moltype.ilist, IF_VSITE))
176     {
177         for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
178         {
179             if (ilist.iatoms[i + 1] == a)
180             {
181                 extremes.minAtom = ilist.iatoms[i + 2];
182                 extremes.maxAtom = ilist.iatoms[i + 2];
183                 for (size_t j = i + 3; j < i + ilistStride(ilist); j++)
184                 {
185                     extremes.minAtom = std::min(extremes.minAtom, ilist.iatoms[j]);
186                     extremes.maxAtom = std::max(extremes.maxAtom, ilist.iatoms[j]);
187                 }
188                 return extremes;
189             }
190         }
191     }
192
193     GMX_RELEASE_ASSERT(false, "If a is a vsite, we should have found constructing atoms");
194
195     return extremes;
196 }
197
198 /*! \brief Returns the range of atoms constrained to atom \p a (including \p a itself) */
199 static AtomIndexExtremes constraintAtomRange(int                     a,
200                                              const ListOfLists<int>& at2con,
201                                              const InteractionList&  ilistConstraints)
202 {
203     AtomIndexExtremes extremes = { a, a };
204
205     for (const int constraint : at2con[a])
206     {
207         for (int j = 0; j < 2; j++)
208         {
209             int atomJ        = ilistConstraints.iatoms[constraint * 3 + 1 + j];
210             extremes.minAtom = std::min(extremes.minAtom, atomJ);
211             extremes.maxAtom = std::max(extremes.maxAtom, atomJ);
212         }
213     }
214
215     return extremes;
216 }
217
218 /*! \brief Returns a list that tells whether atoms in \p moltype are vsites */
219 static std::vector<bool> buildIsParticleVsite(const gmx_moltype_t& moltype)
220 {
221     std::vector<bool> isVsite(moltype.atoms.nr);
222
223     for (auto& ilist : extractILists(moltype.ilist, IF_VSITE))
224     {
225         for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
226         {
227             int vsiteAtom      = ilist.iatoms[i + 1];
228             isVsite[vsiteAtom] = true;
229         }
230     }
231
232     return isVsite;
233 }
234
235 /*! \brief Returns the size of the update group starting at \p firstAtom or 0 when criteria (see updategroups.h) are not met */
236 static int detectGroup(int                     firstAtom,
237                        const gmx_moltype_t&    moltype,
238                        const ListOfLists<int>& at2con,
239                        const InteractionList&  ilistConstraints)
240 {
241     /* We should be using moltype.atoms.atom[].ptype for checking whether
242      * a particle is a vsite. But the test code can't fill t_atoms,
243      * because it uses C pointers which get double freed.
244      */
245     std::vector<bool> isParticleVsite = buildIsParticleVsite(moltype);
246
247     /* A non-vsite atom without constraints is an update group by itself */
248     if (!isParticleVsite[firstAtom] && at2con[firstAtom].empty())
249     {
250         return 1;
251     }
252
253     /* A (potential) update group starts at firstAtom and should consist
254      * of two or more atoms and possibly vsites. At least one atom should
255      * have constraints with all other atoms and vsites should have all
256      * constructing atoms inside the group. Here we increase lastAtom until
257      * the criteria are fulfilled or exit when criteria cannot be met.
258      */
259     int numAtomsWithConstraints = 0;
260     int maxConstraintsPerAtom   = 0;
261     int lastAtom                = firstAtom;
262     int a                       = firstAtom;
263     while (a <= lastAtom)
264     {
265         if (isParticleVsite[a])
266         {
267             AtomIndexExtremes extremes = vsiteConstructRange(a, moltype);
268             if (extremes.minAtom < firstAtom)
269             {
270                 /* A constructing atom is outside the group,
271                  * we can not use update groups.
272                  */
273                 return 0;
274             }
275             lastAtom = std::max(lastAtom, extremes.maxAtom);
276         }
277         else
278         {
279             const int numConstraints = at2con[a].ssize();
280             if (numConstraints == 0)
281             {
282                 /* We can not have unconstrained atoms in an update group */
283                 return 0;
284             }
285             /* This atom has at least one constraint.
286              * Check whether all constraints are within the group
287              * and whether we need to extend the group.
288              */
289             numAtomsWithConstraints += 1;
290             maxConstraintsPerAtom = std::max(maxConstraintsPerAtom, numConstraints);
291
292             AtomIndexExtremes extremes = constraintAtomRange(a, at2con, ilistConstraints);
293             if (extremes.minAtom < firstAtom)
294             {
295                 /* Constraint to atom outside the "group" */
296                 return 0;
297             }
298             lastAtom = std::max(lastAtom, extremes.maxAtom);
299         }
300
301         a++;
302     }
303
304     /* lastAtom might be followed by a vsite that is constructed from atoms
305      * with index <= lastAtom. Handle this case.
306      */
307     if (lastAtom + 1 < moltype.atoms.nr && isParticleVsite[lastAtom + 1])
308     {
309         AtomIndexExtremes extremes = vsiteConstructRange(lastAtom + 1, moltype);
310         if (extremes.minAtom < firstAtom)
311         { // NOLINT bugprone-branch-clone
312             /* Constructing atom precedes the group */
313             return 0;
314         }
315         else if (extremes.maxAtom <= lastAtom)
316         {
317             /* All constructing atoms are in the group, add the vsite to the group */
318             lastAtom++;
319         }
320         else if (extremes.minAtom <= lastAtom)
321         {
322             /* Some, but not all constructing atoms are in the group */
323             return 0;
324         }
325     }
326
327     GMX_RELEASE_ASSERT(maxConstraintsPerAtom < numAtomsWithConstraints,
328                        "We have checked that atoms are only constrained to atoms within the group,"
329                        "so each atom should have fewer constraints than the number of atoms");
330     /* Check that at least one atom is constrained to all others */
331     if (maxConstraintsPerAtom != numAtomsWithConstraints - 1)
332     {
333         return 0;
334     }
335
336     return lastAtom - firstAtom + 1;
337 }
338
339 /*! \brief Returns a list of update groups for \p moltype */
340 static RangePartitioning makeUpdateGroups(const gmx_moltype_t& moltype, gmx::ArrayRef<const t_iparams> iparams)
341 {
342     RangePartitioning groups;
343
344     /* Update groups are not compatible with flexible constraints.
345      * Without dynamics the flexible constraints are ignored,
346      * but since performance for EM/NM is less critical, we do not
347      * use update groups to keep the code here simpler.
348      */
349     if (hasFlexibleConstraints(moltype, iparams) || hasIncompatibleVsites(moltype, iparams))
350     {
351         return groups;
352     }
353
354     /* Combine all constraint ilists into a single one */
355     std::array<InteractionList, F_NRE> ilistsCombined;
356     ilistsCombined[F_CONSTR] = jointConstraintList(moltype);
357     /* We "include" flexible constraints, but none are present (checked above) */
358     const ListOfLists<int> at2con = make_at2con(moltype.atoms.nr, ilistsCombined, iparams,
359                                                 FlexibleConstraintTreatment::Include);
360
361     bool satisfiesCriteria = true;
362
363     int firstAtom = 0;
364     while (satisfiesCriteria && firstAtom < moltype.atoms.nr)
365     {
366         int numAtomsInGroup = detectGroup(firstAtom, moltype, at2con, ilistsCombined[F_CONSTR]);
367
368         if (numAtomsInGroup == 0)
369         {
370             satisfiesCriteria = false;
371         }
372         else
373         {
374             groups.appendBlock(numAtomsInGroup);
375         }
376         firstAtom += numAtomsInGroup;
377     }
378
379     if (!satisfiesCriteria)
380     {
381         /* Make groups empty, to signal not satisfying the criteria */
382         groups.clear();
383     }
384
385     return groups;
386 }
387
388 std::vector<RangePartitioning> makeUpdateGroups(const gmx_mtop_t& mtop)
389 {
390     std::vector<RangePartitioning> updateGroups;
391
392     bool systemSatisfiesCriteria = true;
393     for (const gmx_moltype_t& moltype : mtop.moltype)
394     {
395         updateGroups.push_back(makeUpdateGroups(moltype, mtop.ffparams.iparams));
396
397         if (updateGroups.back().numBlocks() == 0)
398         {
399             systemSatisfiesCriteria = false;
400         }
401     }
402
403     if (!systemSatisfiesCriteria)
404     {
405         updateGroups.clear();
406     }
407
408     return updateGroups;
409 }
410
411 /*! \brief Returns a map of angles ilist.iatoms indices with the middle atom as key */
412 static std::unordered_multimap<int, int> getAngleIndices(const gmx_moltype_t& moltype)
413 {
414     const InteractionList& angles = moltype.ilist[F_ANGLES];
415
416     std::unordered_multimap<int, int> indices(angles.size());
417
418     for (int i = 0; i < angles.size(); i += 1 + NRAL(F_ANGLES))
419     {
420         indices.insert({ angles.iatoms[i + 2], i });
421     }
422
423     return indices;
424 }
425
426 /*! \brief When possible, computes the maximum radius of constrained atom in an update group
427  *
428  * Supports groups with 2 or 3 atoms where all partner atoms are connected to
429  * each other by angle potentials. The temperature is used to compute a radius
430  * that is not exceeded with a chance of 10^-9. Note that this computation
431  * assumes there are no other strong forces working on these angular
432  * degrees of freedom.
433  * The return value is -1 when all partners are not connected to each other
434  * by one angle potential, when a potential is perturbed or when an angle
435  * could reach more than 180 degrees.
436  */
437 template<int numPartnerAtoms>
438 static real constraintGroupRadius(const gmx_moltype_t&                     moltype,
439                                   gmx::ArrayRef<const t_iparams>           iparams,
440                                   const int                                centralAtom,
441                                   const ListOfLists<int>&                  at2con,
442                                   const std::unordered_multimap<int, int>& angleIndices,
443                                   const real                               constraintLength,
444                                   const real                               temperature)
445 {
446     const int numConstraints = at2con[centralAtom].ssize();
447     GMX_RELEASE_ASSERT(numConstraints == numPartnerAtoms,
448                        "We expect as many constraints as partner atoms here");
449
450     std::array<int, numPartnerAtoms> partnerAtoms;
451     for (int i = 0; i < numPartnerAtoms; i++)
452     {
453         const int ind = at2con[centralAtom][i] * 3;
454         if (ind >= moltype.ilist[F_CONSTR].size())
455         {
456             /* This is a flexible constraint, we don't optimize for that */
457             return -1;
458         }
459         const int a1    = moltype.ilist[F_CONSTR].iatoms[ind + 1];
460         const int a2    = moltype.ilist[F_CONSTR].iatoms[ind + 2];
461         partnerAtoms[i] = (a1 == centralAtom ? a2 : a1);
462     }
463
464     const InteractionList&           angles      = moltype.ilist[F_ANGLES];
465     auto                             range       = angleIndices.equal_range(centralAtom);
466     int                              angleType   = -1;
467     std::array<int, numPartnerAtoms> numAngles   = { 0 };
468     bool                             areSameType = true;
469     for (auto it = range.first; it != range.second; ++it)
470     {
471         /* Check if the outer atoms in the angle are both partner atoms */
472         int numAtomsFound = 0;
473         for (int ind = it->second + 1; ind < it->second + 4; ind += 2)
474         {
475             for (const int& partnerAtom : partnerAtoms)
476             {
477                 if (angles.iatoms[ind] == partnerAtom)
478                 {
479                     numAtomsFound++;
480                 }
481             }
482         }
483         if (numAtomsFound == 2)
484         {
485             /* Check that the angle potentials have the same type */
486             if (angleType == -1)
487             {
488                 angleType = angles.iatoms[it->second];
489             }
490             else if (angles.iatoms[it->second] != angleType)
491             {
492                 areSameType = false;
493             }
494             /* Count the number of angle interactions per atoms */
495             for (int ind = it->second + 1; ind < it->second + 4; ind += 2)
496             {
497                 for (size_t i = 0; i < partnerAtoms.size(); i++)
498                 {
499                     if (angles.iatoms[ind] == partnerAtoms[i])
500                     {
501                         numAngles[i]++;
502                     }
503                 }
504             }
505         }
506     }
507
508     bool criteriaSatisfied = areSameType;
509     for (int i = 0; i < numPartnerAtoms; i++)
510     {
511         if (numAngles[i] != numPartnerAtoms - 1)
512         {
513             criteriaSatisfied = false;
514         }
515     }
516
517     /* We don't bother optimizing the perturbed angle case */
518     const t_iparams& angleParams = iparams[angleType];
519     if (criteriaSatisfied && angleParams.harmonic.rB == angleParams.harmonic.rA
520         && angleParams.harmonic.krB == angleParams.harmonic.krA)
521     {
522         /* Set number of stddevs such that change of exceeding < 10^-9 */
523         constexpr real c_numSigma = 6.0;
524         /* Compute the maximally stretched angle */
525         const real eqAngle = angleParams.harmonic.rA * DEG2RAD;
526         const real fc      = angleParams.harmonic.krA;
527         const real maxAngle = eqAngle + c_numSigma * BOLTZ * temperature / ((numPartnerAtoms - 1) * fc);
528         if (maxAngle >= M_PI)
529         {
530             return -1;
531         }
532
533         if (numPartnerAtoms == 2)
534         {
535             /* With two atoms constrainted to a cental atom we have a triangle
536              * with two equal sides because the constraint type is equal.
537              * Return the distance from the COG to the farthest two corners,
538              * i.e. the partner atoms.
539              */
540             real distMidPartner = std::sin(0.5 * maxAngle) * constraintLength;
541             real distCentralMid = std::cos(0.5 * maxAngle) * constraintLength;
542             real distCogMid     = distCentralMid * numPartnerAtoms / (numPartnerAtoms + 1);
543             real distCogPartner = std::sqrt(gmx::square(distMidPartner) + gmx::square(distCogMid));
544
545             return distCogPartner;
546         }
547         else if (numPartnerAtoms == 3)
548         {
549             /* With three atoms constrainted to a cental atom we have the
550              * largest COG-atom distance when one partner atom (the outlier)
551              * moves out by stretching its two angles by an equal amount.
552              * The math here gets a bit more involved, but it is still
553              * rather straightforward geometry.
554              * We first compute distances in the plane of the three partners.
555              * Here we have two "equilibrium" partners and one outlier.
556              * We make use of the "Mid" point between the two "Eq" partners.
557              * We project the central atom on this plane.
558              * Then we compute the distance of the central atom to the plane.
559              * The distance of the COG to the ourlier is returned.
560              */
561             real halfDistEqPartner  = std::sin(0.5 * eqAngle) * constraintLength;
562             real distPartnerOutlier = 2 * std::sin(0.5 * maxAngle) * constraintLength;
563             real distMidOutlier =
564                     std::sqrt(gmx::square(distPartnerOutlier) - gmx::square(halfDistEqPartner));
565             real distMidCenterInPlane =
566                     0.5 * (distMidOutlier - gmx::square(halfDistEqPartner) / distMidOutlier);
567             real distCentralToPlane =
568                     std::sqrt(gmx::square(constraintLength) - gmx::square(halfDistEqPartner)
569                               - gmx::square(distMidCenterInPlane));
570             real distCogOutlierH = distCentralToPlane / (numPartnerAtoms + 1);
571             real distCogOutlierP =
572                     distMidOutlier - (distMidOutlier + distMidCenterInPlane) / (numPartnerAtoms + 1);
573             real distCogOutlier = std::sqrt(gmx::square(distCogOutlierH) + gmx::square(distCogOutlierP));
574
575             return distCogOutlier;
576         }
577         else
578         {
579             GMX_RELEASE_ASSERT(false, "Only 2 or 3 constraints are supported here");
580         }
581     }
582
583     return -1;
584 }
585
586 /*! \brief Returns the maximum update group radius for \p moltype */
587 static real computeMaxUpdateGroupRadius(const gmx_moltype_t&           moltype,
588                                         gmx::ArrayRef<const t_iparams> iparams,
589                                         const RangePartitioning&       updateGroups,
590                                         real                           temperature)
591 {
592     GMX_RELEASE_ASSERT(!hasFlexibleConstraints(moltype, iparams),
593                        "Flexible constraints are not supported here");
594
595     const InteractionList& settles = moltype.ilist[F_SETTLE];
596
597     const ListOfLists<int> at2con = make_at2con(moltype, iparams, FlexibleConstraintTreatment::Include);
598
599     const auto angleIndices = getAngleIndices(moltype);
600
601     real maxRadius = 0;
602     for (int group = 0; group < updateGroups.numBlocks(); group++)
603     {
604         if (updateGroups.block(group).size() == 1)
605         {
606             /* Single atom group, radius is zero */
607             continue;
608         }
609
610         /* Find the atom maxAtom with the maximum number of constraints */
611         int maxNumConstraints = 0;
612         int maxAtom           = -1;
613         for (int a : updateGroups.block(group))
614         {
615             const int numConstraints = at2con[a].ssize();
616             if (numConstraints > maxNumConstraints)
617             {
618                 maxNumConstraints = numConstraints;
619                 maxAtom           = a;
620             }
621         }
622         GMX_ASSERT(maxAtom >= 0 || !settles.empty(),
623                    "We should have at least two atoms in the group with constraints");
624         if (maxAtom < 0)
625         {
626             continue;
627         }
628
629         bool allTypesAreEqual     = true;
630         int  constraintType       = -1;
631         real maxConstraintLength  = 0;
632         real sumConstraintLengths = 0;
633         bool isFirstConstraint    = true;
634         for (const int constraint : at2con[maxAtom])
635         {
636             int conIndex = constraint * (1 + NRAL(F_CONSTR));
637             int iparamsIndex;
638             if (conIndex < moltype.ilist[F_CONSTR].size())
639             {
640                 iparamsIndex = moltype.ilist[F_CONSTR].iatoms[conIndex];
641             }
642             else
643             {
644                 iparamsIndex =
645                         moltype.ilist[F_CONSTRNC].iatoms[conIndex - moltype.ilist[F_CONSTR].size()];
646             }
647             if (isFirstConstraint)
648             {
649                 constraintType    = iparamsIndex;
650                 isFirstConstraint = false;
651             }
652             else if (iparamsIndex != constraintType)
653             {
654                 allTypesAreEqual = false;
655             }
656             /* Here we take the maximum constraint length of the A and B
657              * topology, which assumes lambda is between 0 and 1 for
658              * free-energy runs.
659              */
660             real constraintLength =
661                     std::max(iparams[iparamsIndex].constr.dA, iparams[iparamsIndex].constr.dB);
662             maxConstraintLength = std::max(maxConstraintLength, constraintLength);
663             sumConstraintLengths += constraintLength;
664         }
665
666         int  numConstraints = at2con[maxAtom].ssize();
667         real radius;
668         if (numConstraints == 1)
669         {
670             /* Single constraint: the radius is the distance from the midpoint */
671             radius = 0.5_real * maxConstraintLength;
672         }
673         else
674         {
675             radius = -1;
676
677             /* With 2 constraints the maximum possible radius is the
678              * constraint length, so we can use that for the 0 K case.
679              */
680             if (numConstraints == 2 && allTypesAreEqual && temperature > 0)
681             {
682                 radius = constraintGroupRadius<2>(moltype, iparams, maxAtom, at2con, angleIndices,
683                                                   maxConstraintLength, temperature);
684             }
685             /* With 3 constraints the maximum possible radius is 1.4 times
686              * the constraint length, so it is worth computing a smaller
687              * radius to enable update groups for systems in a small box.
688              */
689             if (numConstraints == 3 && allTypesAreEqual && temperature >= 0)
690             {
691                 radius = constraintGroupRadius<3>(moltype, iparams, maxAtom, at2con, angleIndices,
692                                                   maxConstraintLength, temperature);
693                 if (temperature == 0 && radius >= 0)
694                 {
695                     /* Add a 10% margin for deviation at 0 K */
696                     radius *= 1.1_real;
697                 }
698             }
699
700             if (radius < 0)
701             {
702                 /* Worst case: atom with the longest constraint on one side
703                  * of the center, all others on the opposite side
704                  */
705                 radius = maxConstraintLength
706                          + (sumConstraintLengths - 2 * maxConstraintLength) / (numConstraints + 1);
707             }
708         }
709         maxRadius = std::max(maxRadius, radius);
710     }
711
712     for (int i = 0; i < settles.size(); i += 1 + NRAL(F_SETTLE))
713     {
714         const real dOH = iparams[settles.iatoms[i]].settle.doh;
715         const real dHH = iparams[settles.iatoms[i]].settle.dhh;
716         /* Compute distance^2 from center of geometry to O and H */
717         const real dCO2  = (4 * dOH * dOH - dHH * dHH) / 9;
718         const real dCH2  = (dOH * dOH + 2 * dHH * dHH) / 9;
719         const real dCAny = std::sqrt(std::max(dCO2, dCH2));
720         maxRadius        = std::max(maxRadius, dCAny);
721     }
722
723     return maxRadius;
724 }
725
726 real computeMaxUpdateGroupRadius(const gmx_mtop_t&                      mtop,
727                                  gmx::ArrayRef<const RangePartitioning> updateGroups,
728                                  real                                   temperature)
729 {
730     if (updateGroups.empty())
731     {
732         return 0;
733     }
734
735     GMX_RELEASE_ASSERT(updateGroups.size() == mtop.moltype.size(),
736                        "We need one update group entry per moleculetype");
737
738     real maxRadius = 0;
739
740     for (size_t moltype = 0; moltype < mtop.moltype.size(); moltype++)
741     {
742         maxRadius = std::max(
743                 maxRadius, computeMaxUpdateGroupRadius(mtop.moltype[moltype], mtop.ffparams.iparams,
744                                                        updateGroups[moltype], temperature));
745     }
746
747     return maxRadius;
748 }
749
750 } // namespace gmx