Avoid passing a temporary object to emplace_back
authorAndrey Alekseenko <al42and@gmail.com>
Fri, 2 Jul 2021 06:45:08 +0000 (06:45 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Fri, 2 Jul 2021 06:45:08 +0000 (06:45 +0000)
src/gromacs/gmxpreprocess/add_par.cpp
src/gromacs/gmxpreprocess/gen_ad.cpp
src/gromacs/gmxpreprocess/topio.cpp
src/gromacs/gmxpreprocess/toppush.cpp
src/gromacs/gmxpreprocess/vsite_parm.cpp
src/gromacs/modularsimulator/nosehooverchains.cpp
src/gromacs/pbcutil/com.cpp
src/gromacs/tables/forcetable.cpp

index a07152567859b6f7de6927f445946ba1490bc3e4..4cae7c0168c30055ad68c1e6103871b9ac97abf7 100644 (file)
@@ -4,7 +4,7 @@
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
  * Copyright (c) 2011,2014,2015,2017,2018 by the GROMACS development team.
- * Copyright (c) 2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -62,33 +62,33 @@ void add_param(InteractionsOfType* ps, int ai, int aj, gmx::ArrayRef<const real>
     std::vector<int>  atoms = { ai, aj };
     std::vector<real> forceParm(c.begin(), c.end());
 
-    ps->interactionTypes.emplace_back(InteractionOfType(atoms, forceParm, s ? s : ""));
+    ps->interactionTypes.emplace_back(atoms, forceParm, s ? s : "");
 }
 
 void add_cmap_param(InteractionsOfType* ps, int ai, int aj, int ak, int al, int am, const char* s)
 {
     std::vector<int> atoms = { ai, aj, ak, al, am };
-    ps->interactionTypes.emplace_back(InteractionOfType(atoms, {}, s ? s : ""));
+    ps->interactionTypes.emplace_back(atoms, gmx::ArrayRef<const real>{}, s ? s : "");
 }
 
 void add_vsite2_param(InteractionsOfType* ps, int ai, int aj, int ak, real c0)
 {
     std::vector<int>  atoms     = { ai, aj, ak };
     std::vector<real> forceParm = { c0 };
-    ps->interactionTypes.emplace_back(InteractionOfType(atoms, forceParm));
+    ps->interactionTypes.emplace_back(atoms, forceParm);
 }
 
 void add_vsite3_param(InteractionsOfType* ps, int ai, int aj, int ak, int al, real c0, real c1)
 {
     std::vector<int>  atoms     = { ai, aj, ak, al };
     std::vector<real> forceParm = { c0, c1 };
-    ps->interactionTypes.emplace_back(InteractionOfType(atoms, forceParm));
+    ps->interactionTypes.emplace_back(atoms, forceParm);
 }
 
 void add_vsite3_atoms(InteractionsOfType* ps, int ai, int aj, int ak, int al, bool bSwapParity)
 {
     std::vector<int> atoms = { ai, aj, ak, al };
-    ps->interactionTypes.emplace_back(InteractionOfType(atoms, {}));
+    ps->interactionTypes.emplace_back(atoms, gmx::ArrayRef<const real>{});
 
     if (bSwapParity)
     {
@@ -99,7 +99,7 @@ void add_vsite3_atoms(InteractionsOfType* ps, int ai, int aj, int ak, int al, bo
 void add_vsite4_atoms(InteractionsOfType* ps, int ai, int aj, int ak, int al, int am)
 {
     std::vector<int> atoms = { ai, aj, ak, al, am };
-    ps->interactionTypes.emplace_back(InteractionOfType(atoms, {}));
+    ps->interactionTypes.emplace_back(atoms, gmx::ArrayRef<const real>{});
 }
 
 int search_jtype(const PreprocessResidue& localPpResidue, const char* name, bool bNterm)
index a39d5c881084e5a3e2b29c89704a49c8409d8f3b..f8486762464f4aad870194e31a72c16d3a161d3b 100644 (file)
@@ -247,7 +247,7 @@ static std::vector<InteractionOfType> clean_dih(gmx::ArrayRef<const InteractionO
         int i = 0;
         for (const auto& dihedral : dih)
         {
-            newDihedrals.emplace_back(std::pair<InteractionOfType, int>(dihedral, i++));
+            newDihedrals.emplace_back(dihedral, i++);
         }
     }
     else
@@ -266,7 +266,7 @@ static std::vector<InteractionOfType> clean_dih(gmx::ArrayRef<const InteractionO
             if (was_dihedral_set_in_rtp(*dihedral) || dihedral == dih.begin()
                 || !is_dihedral_on_same_bond(*dihedral, *(dihedral - 1)))
             {
-                newDihedrals.emplace_back(std::pair<InteractionOfType, int>(*dihedral, i++));
+                newDihedrals.emplace_back(*dihedral, i++);
             }
         }
     }
@@ -372,7 +372,7 @@ static std::vector<InteractionOfType> get_impropers(t_atoms*
                 if (!bStop)
                 {
                     /* Not broken out */
-                    improper.emplace_back(InteractionOfType(ai, {}, bondeds.s));
+                    improper.emplace_back(ai, gmx::ArrayRef<const real>{}, bondeds.s);
                 }
             }
             while ((start < atoms->nr) && (atoms->atom[start].resind == i))
index 714b8680ea29446c01e8f832f1450a230c777d73..0535441ecde9fac04cf4613fbf35572844e4e8a2 100644 (file)
@@ -117,7 +117,7 @@ static void gen_pairs(const InteractionsOfType& nbs, InteractionsOfType* pairs,
         /* Copy normal and FEP parameters and multiply by fudge factor */
         gmx::ArrayRef<const real> existingParam = type.forceParam();
         GMX_RELEASE_ASSERT(2 * nrfp <= MAXFORCEPARAM,
-                           "Can't have more parameters than half of maximum p  arameter number");
+                           "Can't have more parameters than half of maximum parameter number");
         for (int j = 0; j < nrfp; j++)
         {
             /* If we are using sigma/epsilon values, only the epsilon values
@@ -137,7 +137,7 @@ static void gen_pairs(const InteractionsOfType& nbs, InteractionsOfType* pairs,
             forceParam[j]        = scaling * existingParam[j];
             forceParam[nrfp + j] = scaling * existingParam[j];
         }
-        pairs->interactionTypes.emplace_back(InteractionOfType(atomNumbers, forceParam));
+        pairs->interactionTypes.emplace_back(atomNumbers, forceParam);
         i++;
     }
 }
index 5e983a0cd34b0f7690e216e48675104bc2355c55..32ed0111a154ee48b017c68c022d2e144b0ab7a9 100644 (file)
@@ -114,7 +114,8 @@ void generate_nbparams(CombinationRule         comb,
                                 c  = std::sqrt(ci * cj);
                                 forceParam[nf] = c;
                             }
-                            interactions->interactionTypes.emplace_back(InteractionOfType({}, forceParam));
+                            interactions->interactionTypes.emplace_back(gmx::ArrayRef<const int>{},
+                                                                        forceParam);
                         }
                     }
                     break;
@@ -138,7 +139,8 @@ void generate_nbparams(CombinationRule         comb,
                                 forceParam[0] *= -1;
                             }
                             forceParam[1] = std::sqrt(ci1 * cj1);
-                            interactions->interactionTypes.emplace_back(InteractionOfType({}, forceParam));
+                            interactions->interactionTypes.emplace_back(gmx::ArrayRef<const int>{},
+                                                                        forceParam);
                         }
                     }
 
@@ -162,7 +164,8 @@ void generate_nbparams(CombinationRule         comb,
                                 forceParam[0] *= -1;
                             }
                             forceParam[1] = std::sqrt(ci1 * cj1);
-                            interactions->interactionTypes.emplace_back(InteractionOfType({}, forceParam));
+                            interactions->interactionTypes.emplace_back(gmx::ArrayRef<const int>{},
+                                                                        forceParam);
                         }
                     }
 
@@ -196,7 +199,7 @@ void generate_nbparams(CombinationRule         comb,
                         forceParam[1] = 2.0 / (1 / bi + 1 / bj);
                     }
                     forceParam[2] = std::sqrt(ci2 * cj2);
-                    interactions->interactionTypes.emplace_back(InteractionOfType({}, forceParam));
+                    interactions->interactionTypes.emplace_back(gmx::ArrayRef<const int>{}, forceParam);
                 }
             }
 
@@ -725,8 +728,7 @@ static void push_bondtype(InteractionsOfType*      bt,
     if (addBondType)
     {
         /* fill the arrays up and down */
-        bt->interactionTypes.emplace_back(
-                InteractionOfType(b.atoms(), b.forceParam(), b.interactionTypeName()));
+        bt->interactionTypes.emplace_back(b.atoms(), b.forceParam(), b.interactionTypeName());
         /* need to store force values because they might change below */
         std::vector<real> forceParam(b.forceParam().begin(), b.forceParam().end());
 
@@ -745,7 +747,7 @@ static void push_bondtype(InteractionsOfType*      bt,
         {
             atoms.emplace_back(*oldAtom);
         }
-        bt->interactionTypes.emplace_back(InteractionOfType(atoms, forceParam, b.interactionTypeName()));
+        bt->interactionTypes.emplace_back(atoms, forceParam, b.interactionTypeName());
     }
 }
 
@@ -2694,7 +2696,7 @@ static void convert_pairs_to_pairsQ(gmx::ArrayRef<InteractionsOfType> interactio
         std::vector<real> forceParam = {
             fudgeQQ, atoms->atom[param.ai()].q, atoms->atom[param.aj()].q, param.c0(), param.c1()
         };
-        paramnew.emplace_back(InteractionOfType(param.atoms(), forceParam, ""));
+        paramnew.emplace_back(param.atoms(), forceParam, "");
     }
 
     /* now assign the new data to the F_LJC14_Q structure */
index 9e6eabffd0594ee47bd5df65150a4da898df828d..655b5673a8f38e22da2578fca283749a91dbde84 100644 (file)
@@ -159,7 +159,7 @@ static int vsite_bond_nrcheck(int ftype)
 static void enter_bonded(int nratoms, std::vector<VsiteBondedInteraction>* bondeds, const InteractionOfType& type)
 {
     GMX_RELEASE_ASSERT(nratoms == type.atoms().ssize(), "Size of atom array must match");
-    bondeds->emplace_back(VsiteBondedInteraction(type.atoms(), type.c0()));
+    bondeds->emplace_back(type.atoms(), type.c0());
 }
 
 /*! \internal \brief
index 02b411553131c5b21021f7de904f9bf5d30f6e0d..72d3275e51fb81f68784cd25cbd4dd23ccc57e59 100644 (file)
@@ -175,12 +175,12 @@ NoseHooverChainsData::NoseHooverChainsData(int                  numTemperatureGr
     {
         for (auto temperatureGroup = 0; temperatureGroup < numTemperatureGroups; ++temperatureGroup)
         {
-            noseHooverGroups_.emplace_back(NoseHooverGroup(chainLength,
-                                                           referenceTemperature[temperatureGroup],
-                                                           numDegreesOfFreedom[temperatureGroup],
-                                                           couplingTime[temperatureGroup],
-                                                           couplingTimeStep,
-                                                           nhcUsage));
+            noseHooverGroups_.emplace_back(chainLength,
+                                           referenceTemperature[temperatureGroup],
+                                           numDegreesOfFreedom[temperatureGroup],
+                                           couplingTime[temperatureGroup],
+                                           couplingTimeStep,
+                                           nhcUsage);
         }
     }
     else if (nhcUsage == NhcUsage::Barostat)
@@ -189,8 +189,8 @@ NoseHooverChainsData::NoseHooverChainsData(int                  numTemperatureGr
                            "There can only be one barostat for the system");
         // Barostat has a single degree of freedom
         const int degreesOfFreedom = 1;
-        noseHooverGroups_.emplace_back(NoseHooverGroup(
-                chainLength, referenceTemperature[0], degreesOfFreedom, couplingTime[0], couplingTimeStep, nhcUsage));
+        noseHooverGroups_.emplace_back(
+                chainLength, referenceTemperature[0], degreesOfFreedom, couplingTime[0], couplingTimeStep, nhcUsage);
     }
 }
 
index 3aef25fdd569d850dac9cff8fbd0b4b1e02215f3..f1be05722543694c51a1f693c6e28392946c6fee 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -148,7 +148,7 @@ void placeCoordinatesWithCOMInBox(const PbcType&      pbcType,
             std::vector<Range<int>> atomRanges;
             if (comShiftType == COMShiftType::Molecule)
             {
-                atomRanges.emplace_back(gmx::Range<int>(0, atomsPerMolecule));
+                atomRanges.emplace_back(0, atomsPerMolecule);
             }
             else if (comShiftType == COMShiftType::Residue)
             {
index 54a08a65438067987b70e3d93d55d5f6460c12f2..d8ad10e2cc5d85854f1402b013453d46b814a211 100644 (file)
@@ -754,7 +754,7 @@ static std::vector<t_tabledata> read_tables(FILE* fp, const char* filename, int
     std::vector<t_tabledata> td;
     for (k = 0; (k < ntab); k++)
     {
-        td.emplace_back(t_tabledata(numRows, nx0, tabscale, true));
+        td.emplace_back(numRows, nx0, tabscale, true);
         for (i = 0; (i < numRows); i++)
         {
             td[k].x[i] = yy[0][i];