5ee169a4737271a94532bc99e41145d21ed2b87d
[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, 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/block.h"
56 #include "gromacs/topology/idef.h"
57 #include "gromacs/topology/ifunc.h"
58 #include "gromacs/topology/topology.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.data(), 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, const t_blocka& at2con, const InteractionList& ilistConstraints)
200 {
201     AtomIndexExtremes extremes = { a, a };
202
203     for (int i = at2con.index[a]; i < at2con.index[a + 1]; i++)
204     {
205         for (int j = 0; j < 2; j++)
206         {
207             int atomJ        = ilistConstraints.iatoms[at2con.a[i] * 3 + 1 + j];
208             extremes.minAtom = std::min(extremes.minAtom, atomJ);
209             extremes.maxAtom = std::max(extremes.maxAtom, atomJ);
210         }
211     }
212
213     return extremes;
214 }
215
216 /*! \brief Returns a list that tells whether atoms in \p moltype are vsites */
217 static std::vector<bool> buildIsParticleVsite(const gmx_moltype_t& moltype)
218 {
219     std::vector<bool> isVsite(moltype.atoms.nr);
220
221     for (auto& ilist : extractILists(moltype.ilist, IF_VSITE))
222     {
223         for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
224         {
225             int vsiteAtom      = ilist.iatoms[i + 1];
226             isVsite[vsiteAtom] = true;
227         }
228     }
229
230     return isVsite;
231 }
232
233 /*! \brief Returns the size of the update group starting at \p firstAtom or 0 when criteria (see updategroups.h) are not met */
234 static int detectGroup(int                    firstAtom,
235                        const gmx_moltype_t&   moltype,
236                        const t_blocka&        at2con,
237                        const InteractionList& ilistConstraints)
238 {
239     /* We should be using moltype.atoms.atom[].ptype for checking whether
240      * a particle is a vsite. But the test code can't fill t_atoms,
241      * because it uses C pointers which get double freed.
242      */
243     std::vector<bool> isParticleVsite = buildIsParticleVsite(moltype);
244
245     /* A non-vsite atom without constraints is an update group by itself */
246     if (!isParticleVsite[firstAtom] && at2con.index[firstAtom + 1] - at2con.index[firstAtom] == 0)
247     {
248         return 1;
249     }
250
251     /* A (potential) update group starts at firstAtom and should consist
252      * of two or more atoms and possibly vsites. At least one atom should
253      * have constraints with all other atoms and vsites should have all
254      * constructing atoms inside the group. Here we increase lastAtom until
255      * the criteria are fulfilled or exit when criteria cannot be met.
256      */
257     int numAtomsWithConstraints = 0;
258     int maxConstraintsPerAtom   = 0;
259     int lastAtom                = firstAtom;
260     int a                       = firstAtom;
261     while (a <= lastAtom)
262     {
263         if (isParticleVsite[a])
264         {
265             AtomIndexExtremes extremes = vsiteConstructRange(a, moltype);
266             if (extremes.minAtom < firstAtom)
267             {
268                 /* A constructing atom is outside the group,
269                  * we can not use update groups.
270                  */
271                 return 0;
272             }
273             lastAtom = std::max(lastAtom, extremes.maxAtom);
274         }
275         else
276         {
277             int numConstraints = at2con.index[a + 1] - at2con.index[a];
278             if (numConstraints == 0)
279             {
280                 /* We can not have unconstrained atoms in an update group */
281                 return 0;
282             }
283             /* This atom has at least one constraint.
284              * Check whether all constraints are within the group
285              * and whether we need to extend the group.
286              */
287             numAtomsWithConstraints += 1;
288             maxConstraintsPerAtom = std::max(maxConstraintsPerAtom, numConstraints);
289
290             AtomIndexExtremes extremes = constraintAtomRange(a, at2con, ilistConstraints);
291             if (extremes.minAtom < firstAtom)
292             {
293                 /* Constraint to atom outside the "group" */
294                 return 0;
295             }
296             lastAtom = std::max(lastAtom, extremes.maxAtom);
297         }
298
299         a++;
300     }
301
302     /* lastAtom might be followed by a vsite that is constructed from atoms
303      * with index <= lastAtom. Handle this case.
304      */
305     if (lastAtom + 1 < moltype.atoms.nr && isParticleVsite[lastAtom + 1])
306     {
307         AtomIndexExtremes extremes = vsiteConstructRange(lastAtom + 1, moltype);
308         if (extremes.minAtom < firstAtom)
309         {
310             /* Constructing atom precedes the group */
311             return 0;
312         }
313         else if (extremes.maxAtom <= lastAtom)
314         {
315             /* All constructing atoms are in the group, add the vsite to the group */
316             lastAtom++;
317         }
318         else if (extremes.minAtom <= lastAtom)
319         {
320             /* Some, but not all constructing atoms are in the group */
321             return 0;
322         }
323     }
324
325     GMX_RELEASE_ASSERT(maxConstraintsPerAtom < numAtomsWithConstraints,
326                        "We have checked that atoms are only constrained to atoms within the group,"
327                        "so each atom should have fewer constraints than the number of atoms");
328     /* Check that at least one atom is constrained to all others */
329     if (maxConstraintsPerAtom != numAtomsWithConstraints - 1)
330     {
331         return 0;
332     }
333
334     return lastAtom - firstAtom + 1;
335 }
336
337 /*! \brief Returns a list of update groups for \p moltype */
338 static RangePartitioning makeUpdateGroups(const gmx_moltype_t& moltype, gmx::ArrayRef<const t_iparams> iparams)
339 {
340     RangePartitioning groups;
341
342     /* Update groups are not compatible with flexible constraints.
343      * Without dynamics the flexible constraints are ignored,
344      * but since performance for EM/NM is less critical, we do not
345      * use update groups to keep the code here simpler.
346      */
347     if (hasFlexibleConstraints(moltype, iparams) || hasIncompatibleVsites(moltype, iparams))
348     {
349         return groups;
350     }
351
352     /* Combine all constraint ilists into a single one */
353     InteractionList constraintsCombined = jointConstraintList(moltype);
354     t_ilist         ilistsCombined[F_NRE];
355     ilistsCombined[F_CONSTR].nr     = constraintsCombined.iatoms.size();
356     ilistsCombined[F_CONSTR].iatoms = constraintsCombined.iatoms.data();
357     ilistsCombined[F_CONSTRNC].nr   = 0;
358     /* We "include" flexible constraints, but none are present (checked above) */
359     t_blocka at2con = make_at2con(moltype.atoms.nr, ilistsCombined, iparams.data(),
360                                   FlexibleConstraintTreatment::Include);
361
362     bool satisfiesCriteria = true;
363
364     int firstAtom = 0;
365     while (satisfiesCriteria && firstAtom < moltype.atoms.nr)
366     {
367         int numAtomsInGroup = detectGroup(firstAtom, moltype, at2con, constraintsCombined);
368
369         if (numAtomsInGroup == 0)
370         {
371             satisfiesCriteria = false;
372         }
373         else
374         {
375             groups.appendBlock(numAtomsInGroup);
376         }
377         firstAtom += numAtomsInGroup;
378     }
379
380     if (!satisfiesCriteria)
381     {
382         /* Make groups empty, to signal not satisfying the criteria */
383         groups.clear();
384     }
385
386     done_blocka(&at2con);
387
388     return groups;
389 }
390
391 std::vector<RangePartitioning> makeUpdateGroups(const gmx_mtop_t& mtop)
392 {
393     std::vector<RangePartitioning> updateGroups;
394
395     bool systemSatisfiesCriteria = true;
396     for (const gmx_moltype_t& moltype : mtop.moltype)
397     {
398         updateGroups.push_back(makeUpdateGroups(moltype, mtop.ffparams.iparams));
399
400         if (updateGroups.back().numBlocks() == 0)
401         {
402             systemSatisfiesCriteria = false;
403         }
404     }
405
406     if (!systemSatisfiesCriteria)
407     {
408         updateGroups.clear();
409     }
410
411     return updateGroups;
412 }
413
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)
416 {
417     const InteractionList& angles = moltype.ilist[F_ANGLES];
418
419     std::unordered_multimap<int, int> indices(angles.size());
420
421     for (int i = 0; i < angles.size(); i += 1 + NRAL(F_ANGLES))
422     {
423         indices.insert({ angles.iatoms[i + 2], i });
424     }
425
426     return indices;
427 }
428
429 /*! \brief When possible, computes the maximum radius of constrained atom in an update group
430  *
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.
439  */
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 t_blocka&                          at2con,
445                                   const std::unordered_multimap<int, int>& angleIndices,
446                                   const real                               constraintLength,
447                                   const real                               temperature)
448 {
449     const int numConstraints = at2con.index[centralAtom + 1] - at2con.index[centralAtom];
450     GMX_RELEASE_ASSERT(numConstraints == numPartnerAtoms,
451                        "We expect as many constraints as partner atoms here");
452
453     std::array<int, numPartnerAtoms> partnerAtoms;
454     for (int i = 0; i < numPartnerAtoms; i++)
455     {
456         const int ind = at2con.a[at2con.index[centralAtom] + i] * 3;
457         if (ind >= moltype.ilist[F_CONSTR].size())
458         {
459             /* This is a flexible constraint, we don't optimize for that */
460             return -1;
461         }
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);
465     }
466
467     const InteractionList&           angles      = moltype.ilist[F_ANGLES];
468     auto                             range       = angleIndices.equal_range(centralAtom);
469     int                              angleType   = -1;
470     std::array<int, numPartnerAtoms> numAngles   = { 0 };
471     bool                             areSameType = true;
472     for (auto it = range.first; it != range.second; ++it)
473     {
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)
477         {
478             for (const int& partnerAtom : partnerAtoms)
479             {
480                 if (angles.iatoms[ind] == partnerAtom)
481                 {
482                     numAtomsFound++;
483                 }
484             }
485         }
486         if (numAtomsFound == 2)
487         {
488             /* Check that the angle potentials have the same type */
489             if (angleType == -1)
490             {
491                 angleType = angles.iatoms[it->second];
492             }
493             else if (angles.iatoms[it->second] != angleType)
494             {
495                 areSameType = false;
496             }
497             /* Count the number of angle interactions per atoms */
498             for (int ind = it->second + 1; ind < it->second + 4; ind += 2)
499             {
500                 for (size_t i = 0; i < partnerAtoms.size(); i++)
501                 {
502                     if (angles.iatoms[ind] == partnerAtoms[i])
503                     {
504                         numAngles[i]++;
505                     }
506                 }
507             }
508         }
509     }
510
511     bool criteriaSatisfied = areSameType;
512     for (int i = 0; i < numPartnerAtoms; i++)
513     {
514         if (numAngles[i] != numPartnerAtoms - 1)
515         {
516             criteriaSatisfied = false;
517         }
518     }
519
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)
524     {
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 * DEG2RAD;
529         const real fc      = angleParams.harmonic.krA;
530         const real maxAngle = eqAngle + c_numSigma * BOLTZ * temperature / ((numPartnerAtoms - 1) * fc);
531         if (maxAngle >= M_PI)
532         {
533             return -1;
534         }
535
536         if (numPartnerAtoms == 2)
537         {
538             /* With two atoms constrainted to a cental atom we have a triangle
539              * with two equal sides because the constraint type is equal.
540              * Return the distance from the COG to the farthest two corners,
541              * i.e. the partner atoms.
542              */
543             real distMidPartner = std::sin(0.5 * maxAngle) * constraintLength;
544             real distCentralMid = std::cos(0.5 * maxAngle) * constraintLength;
545             real distCogMid     = distCentralMid * numPartnerAtoms / (numPartnerAtoms + 1);
546             real distCogPartner = std::sqrt(gmx::square(distMidPartner) + gmx::square(distCogMid));
547
548             return distCogPartner;
549         }
550         else if (numPartnerAtoms == 3)
551         {
552             /* With three atoms constrainted to a cental atom we have the
553              * largest COG-atom distance when one partner atom (the outlier)
554              * moves out by stretching its two angles by an equal amount.
555              * The math here gets a bit more involved, but it is still
556              * rather straightforward geometry.
557              * We first compute distances in the plane of the three partners.
558              * Here we have two "equilibrium" partners and one outlier.
559              * We make use of the "Mid" point between the two "Eq" partners.
560              * We project the central atom on this plane.
561              * Then we compute the distance of the central atom to the plane.
562              * The distance of the COG to the ourlier is returned.
563              */
564             real halfDistEqPartner  = std::sin(0.5 * eqAngle) * constraintLength;
565             real distPartnerOutlier = 2 * std::sin(0.5 * maxAngle) * constraintLength;
566             real distMidOutlier =
567                     std::sqrt(gmx::square(distPartnerOutlier) - gmx::square(halfDistEqPartner));
568             real distMidCenterInPlane =
569                     0.5 * (distMidOutlier - gmx::square(halfDistEqPartner) / distMidOutlier);
570             real distCentralToPlane =
571                     std::sqrt(gmx::square(constraintLength) - gmx::square(halfDistEqPartner)
572                               - gmx::square(distMidCenterInPlane));
573             real distCogOutlierH = distCentralToPlane / (numPartnerAtoms + 1);
574             real distCogOutlierP =
575                     distMidOutlier - (distMidOutlier + distMidCenterInPlane) / (numPartnerAtoms + 1);
576             real distCogOutlier = std::sqrt(gmx::square(distCogOutlierH) + gmx::square(distCogOutlierP));
577
578             return distCogOutlier;
579         }
580         else
581         {
582             GMX_RELEASE_ASSERT(false, "Only 2 or 3 constraints are supported here");
583         }
584     }
585
586     return -1;
587 }
588
589 /*! \brief Returns the maximum update group radius for \p moltype */
590 static real computeMaxUpdateGroupRadius(const gmx_moltype_t&           moltype,
591                                         gmx::ArrayRef<const t_iparams> iparams,
592                                         const RangePartitioning&       updateGroups,
593                                         real                           temperature)
594 {
595     GMX_RELEASE_ASSERT(!hasFlexibleConstraints(moltype, iparams),
596                        "Flexible constraints are not supported here");
597
598     const InteractionList& settles = moltype.ilist[F_SETTLE];
599
600     t_blocka at2con = make_at2con(moltype, iparams, FlexibleConstraintTreatment::Include);
601
602     const auto angleIndices = getAngleIndices(moltype);
603
604     real maxRadius = 0;
605     for (int group = 0; group < updateGroups.numBlocks(); group++)
606     {
607         if (updateGroups.block(group).size() == 1)
608         {
609             /* Single atom group, radius is zero */
610             continue;
611         }
612
613         /* Find the atom maxAtom with the maximum number of constraints */
614         int maxNumConstraints = 0;
615         int maxAtom           = -1;
616         for (int a : updateGroups.block(group))
617         {
618             int numConstraints = at2con.index[a + 1] - at2con.index[a];
619             if (numConstraints > maxNumConstraints)
620             {
621                 maxNumConstraints = numConstraints;
622                 maxAtom           = a;
623             }
624         }
625         GMX_ASSERT(maxAtom >= 0 || settles.size() > 0,
626                    "We should have at least two atoms in the group with constraints");
627         if (maxAtom < 0)
628         {
629             continue;
630         }
631
632         bool allTypesAreEqual     = true;
633         int  constraintType       = -1;
634         real maxConstraintLength  = 0;
635         real sumConstraintLengths = 0;
636         for (int i = at2con.index[maxAtom]; i < at2con.index[maxAtom + 1]; i++)
637         {
638             int conIndex = at2con.a[i] * (1 + NRAL(F_CONSTR));
639             int iparamsIndex;
640             if (conIndex < moltype.ilist[F_CONSTR].size())
641             {
642                 iparamsIndex = moltype.ilist[F_CONSTR].iatoms[conIndex];
643             }
644             else
645             {
646                 iparamsIndex =
647                         moltype.ilist[F_CONSTRNC].iatoms[conIndex - moltype.ilist[F_CONSTR].size()];
648             }
649             if (i == at2con.index[maxAtom])
650             {
651                 constraintType = iparamsIndex;
652             }
653             else if (iparamsIndex != constraintType)
654             {
655                 allTypesAreEqual = false;
656             }
657             /* Here we take the maximum constraint length of the A and B
658              * topology, which assumes lambda is between 0 and 1 for
659              * free-energy runs.
660              */
661             real constraintLength =
662                     std::max(iparams[iparamsIndex].constr.dA, iparams[iparamsIndex].constr.dB);
663             maxConstraintLength = std::max(maxConstraintLength, constraintLength);
664             sumConstraintLengths += constraintLength;
665         }
666
667         int  numConstraints = at2con.index[maxAtom + 1] - at2con.index[maxAtom];
668         real radius;
669         if (numConstraints == 1)
670         {
671             /* Single constraint: the radius is the distance from the midpoint */
672             radius = 0.5_real * maxConstraintLength;
673         }
674         else
675         {
676             radius = -1;
677
678             /* With 2 constraints the maximum possible radius is the
679              * constraint length, so we can use that for the 0 K case.
680              */
681             if (numConstraints == 2 && allTypesAreEqual && temperature > 0)
682             {
683                 radius = constraintGroupRadius<2>(moltype, iparams, maxAtom, at2con, angleIndices,
684                                                   maxConstraintLength, temperature);
685             }
686             /* With 3 constraints the maximum possible radius is 1.4 times
687              * the constraint length, so it is worth computing a smaller
688              * radius to enable update groups for systems in a small box.
689              */
690             if (numConstraints == 3 && allTypesAreEqual && temperature >= 0)
691             {
692                 radius = constraintGroupRadius<3>(moltype, iparams, maxAtom, at2con, angleIndices,
693                                                   maxConstraintLength, temperature);
694                 if (temperature == 0 && radius >= 0)
695                 {
696                     /* Add a 10% margin for deviation at 0 K */
697                     radius *= 1.1_real;
698                 }
699             }
700
701             if (radius < 0)
702             {
703                 /* Worst case: atom with the longest constraint on one side
704                  * of the center, all others on the opposite side
705                  */
706                 radius = maxConstraintLength
707                          + (sumConstraintLengths - 2 * maxConstraintLength) / (numConstraints + 1);
708             }
709         }
710         maxRadius = std::max(maxRadius, radius);
711     }
712
713     for (int i = 0; i < settles.size(); i += 1 + NRAL(F_SETTLE))
714     {
715         const real dOH = iparams[settles.iatoms[i]].settle.doh;
716         const real dHH = iparams[settles.iatoms[i]].settle.dhh;
717         /* Compute distance^2 from center of geometry to O and H */
718         const real dCO2  = (4 * dOH * dOH - dHH * dHH) / 9;
719         const real dCH2  = (dOH * dOH + 2 * dHH * dHH) / 9;
720         const real dCAny = std::sqrt(std::max(dCO2, dCH2));
721         maxRadius        = std::max(maxRadius, dCAny);
722     }
723
724     done_blocka(&at2con);
725
726     return maxRadius;
727 }
728
729 real computeMaxUpdateGroupRadius(const gmx_mtop_t&                      mtop,
730                                  gmx::ArrayRef<const RangePartitioning> updateGroups,
731                                  real                                   temperature)
732 {
733     if (updateGroups.empty())
734     {
735         return 0;
736     }
737
738     GMX_RELEASE_ASSERT(updateGroups.size() == mtop.moltype.size(),
739                        "We need one update group entry per moleculetype");
740
741     real maxRadius = 0;
742
743     for (size_t moltype = 0; moltype < mtop.moltype.size(); moltype++)
744     {
745         maxRadius = std::max(
746                 maxRadius, computeMaxUpdateGroupRadius(mtop.moltype[moltype], mtop.ffparams.iparams,
747                                                        updateGroups[moltype], temperature));
748     }
749
750     return maxRadius;
751 }
752
753 } // namespace gmx