Avoid expensive atom-type lookups in grompp
authorMark Abraham <mark.j.abraham@gmail.com>
Sat, 7 Aug 2021 14:38:20 +0000 (14:38 +0000)
committerAndrey Alekseenko <al42and@gmail.com>
Sat, 7 Aug 2021 14:38:20 +0000 (14:38 +0000)
For each dihedral interaction for which parameters needed to be looked
up, many non-inline function calls were made to see if dihedral types
matched. Instead they are now compared as two ranges with std::equal.

src/gromacs/gmxpreprocess/toppush.cpp

index c73e821844e9b36a2a72a97ad6f73ae9c355b9b7..087695001b58e036ee755d07fde2825b2f7e4116 100644 (file)
@@ -1727,18 +1727,31 @@ static int findNumberOfDihedralAtomMatches(const InteractionOfType&       bondTy
                                            const gmx::ArrayRef<const int> atomTypes)
 {
     GMX_RELEASE_ASSERT(atomTypes.size() == 4, "Dihedrals have 4 atom types");
-    if ((bondType.ai() == -1 || atomTypes[0] == bondType.ai())
-        && (bondType.aj() == -1 || atomTypes[1] == bondType.aj())
-        && (bondType.ak() == -1 || atomTypes[2] == bondType.ak())
-        && (bondType.al() == -1 || atomTypes[3] == bondType.al()))
-    {
-        return (bondType.ai() == -1 ? 0 : 1) + (bondType.aj() == -1 ? 0 : 1)
-               + (bondType.ak() == -1 ? 0 : 1) + (bondType.al() == -1 ? 0 : 1);
-    }
-    else
-    {
-        return -1;
-    }
+    const gmx::ArrayRef<const int> bondTypeAtomTypes = bondType.atoms();
+    GMX_RELEASE_ASSERT(bondTypeAtomTypes.size() == 4, "Dihedral types have 4 atom types");
+    int numExactMatches = 0;
+    if (std::equal(bondTypeAtomTypes.begin(),
+                   bondTypeAtomTypes.end(),
+                   atomTypes.begin(),
+                   [&numExactMatches](int bondTypeAtomType, int atomType) {
+                       if (bondTypeAtomType == atomType)
+                       {
+                           // Found an exact atom type match
+                           ++numExactMatches;
+                           return true;
+                       }
+                       else if (bondTypeAtomType == -1)
+                       {
+                           // Found a wildcard atom type match
+                           return true;
+                       }
+                       // Atom types do not match
+                       return false;
+                   }))
+    {
+        return numExactMatches;
+    }
+    return -1;
 }
 
 static std::vector<InteractionOfType>::iterator
@@ -1754,7 +1767,7 @@ defaultInteractionsOfType(int                               ftype,
         int nmatch_max = -1;
 
         /* For dihedrals we allow wildcards. We choose the first type
-         * that has the most real matches, i.e. non-wildcard matches.
+         * that has the most exact matches, i.e. non-wildcard matches.
          */
         auto prevPos = bondType[ftype].interactionTypes.end();
         auto pos     = bondType[ftype].interactionTypes.begin();