Remove EmptyArrayRef
authorRoland Schulz <roland.schulz@intel.com>
Fri, 8 Feb 2019 21:45:15 +0000 (13:45 -0800)
committerDavid van der Spoel <spoel@xray.bmc.uu.se>
Mon, 11 Feb 2019 16:23:28 +0000 (17:23 +0100)
Make it more aligned with std::span.
Almost everywhere it wasn't needed anyhow.
Only exception: ternary conditional operator. But there the type
was already explict in all but one case and that one case
(src/gromacs/domdec/distribute.cpp) was confusing because of
multiple implicit conversions.

Related #2859

Change-Id: I0b61abdc2e60285a7ca17e3bdc7e53cb72c5cf6a

24 files changed:
src/gromacs/applied_forces/tests/electricfield.cpp
src/gromacs/commandline/filenm.cpp
src/gromacs/commandline/tests/pargs.cpp
src/gromacs/domdec/collect.cpp
src/gromacs/domdec/distribute.cpp
src/gromacs/domdec/domdec_constraints.cpp
src/gromacs/domdec/domdec_topology.cpp
src/gromacs/ewald/pme_only.cpp
src/gromacs/gmxana/gmx_trjcat.cpp
src/gromacs/gmxpreprocess/x2top.cpp
src/gromacs/math/arrayrefwithpadding.h
src/gromacs/math/tests/arrayrefwithpadding.cpp
src/gromacs/mdlib/constr.cpp
src/gromacs/mdlib/mdoutf.cpp
src/gromacs/mdlib/simulationsignal.cpp
src/gromacs/mdrun/minimize.cpp
src/gromacs/mdtypes/state.h
src/gromacs/selection/tests/nbsearch.cpp
src/gromacs/selection/tests/poscalc.cpp
src/gromacs/simd/simd_memory.h
src/gromacs/simd/tests/simd_memory.cpp
src/gromacs/trajectory/trajectoryframe.cpp
src/gromacs/utility/arrayref.h
src/gromacs/utility/tests/arrayref.cpp

index 82c983e0efd94529d1b67260c033f8aa13081133..ebf0f701d5b4ccab9bc5568bf6caf2b238158547 100644 (file)
@@ -120,7 +120,7 @@ class ElectricFieldTest : public ::testing::Test
             matrix                         boxDummy = { {0, 0, 0}, {0, 0, 0}, {0, 0, 0} };
             gmx_enerdata_t                 enerdDummy;
 
-            gmx::ForceProviderInput        forceProviderInput(gmx::EmptyArrayRef {}, md, 0.0, boxDummy, *cr);
+            gmx::ForceProviderInput        forceProviderInput({}, md, 0.0, boxDummy, *cr);
             gmx::ForceProviderOutput       forceProviderOutput(&forceWithVirial, &enerdDummy);
             forceProviders.calculateForces(forceProviderInput, &forceProviderOutput);
             done_commrec(cr);
index 0f0a1e2b327887f2da6281364074a3466a11c0ae..1d42e148a6edd699187bf60f7260f38d5c823292 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
@@ -102,7 +102,7 @@ opt2fns(const char *opt, int nfile, const t_filenm fnm[])
 
     GMX_RELEASE_ASSERT(false, "opt2fns should be called with a valid option");
 
-    return gmx::EmptyArrayRef();
+    return {};
 }
 
 gmx::ArrayRef<const std::string>
@@ -114,7 +114,7 @@ opt2fnsIfOptionSet(const char *opt, int nfile, const t_filenm fnm[])
     }
     else
     {
-        return gmx::EmptyArrayRef();
+        return {};
     }
 }
 
@@ -148,7 +148,7 @@ ftp2fns(int ftp, int nfile, const t_filenm fnm[])
 
     GMX_RELEASE_ASSERT(false, "ftp2fns should be called with a valid option");
 
-    return gmx::EmptyArrayRef();
+    return {};
 }
 
 gmx_bool ftp2bSet(int ftp, int nfile, const t_filenm fnm[])
index 2717afa8c589ac29da94edb3138f9fb270bcb69c..a0441754b1974a4e7c21c9cbdeaff3ba40f42572 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
@@ -152,7 +152,7 @@ TEST_F(ParseCommonArgsTest, ParsesIntegerArgs)
     const char *const cmdline[] = {
         "test", "-i1", "2", "-i2", "-3"
     };
-    parseFromArray(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    parseFromArray(cmdline, 0, {}, pa);
     EXPECT_EQ( 2, value1);
     EXPECT_EQ(-3, value2);
     EXPECT_EQ( 3, value3);
@@ -169,7 +169,7 @@ TEST_F(ParseCommonArgsTest, ParsesInt64Args)
     const char *const cmdline[] = {
         "test", "-i1", "2", "-i2", "-3"
     };
-    parseFromArray(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    parseFromArray(cmdline, 0, {}, pa);
     EXPECT_EQ( 2, value1);
     EXPECT_EQ(-3, value2);
     EXPECT_EQ( 3, value3);
@@ -186,7 +186,7 @@ TEST_F(ParseCommonArgsTest, ParsesRealArgs)
     const char *const cmdline[] = {
         "test", "-r1", "2", "-r2", "-.5"
     };
-    parseFromArray(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    parseFromArray(cmdline, 0, {}, pa);
     EXPECT_EQ( 2.0, value1);
     EXPECT_EQ(-0.5, value2);
     EXPECT_EQ( 2.5, value3);
@@ -203,7 +203,7 @@ TEST_F(ParseCommonArgsTest, ParsesStringArgs)
     const char *const cmdline[] = {
         "test", "-s1", "", "-s2", "test"
     };
-    parseFromArray(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    parseFromArray(cmdline, 0, {}, pa);
     EXPECT_STREQ("", value1);
     EXPECT_STREQ("test", value2);
     EXPECT_STREQ("default", value3);
@@ -220,7 +220,7 @@ TEST_F(ParseCommonArgsTest, ParsesBooleanArgs)
     const char *const cmdline[] = {
         "test", "-nob1", "-b2"
     };
-    parseFromArray(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    parseFromArray(cmdline, 0, {}, pa);
     EXPECT_FALSE(value1);
     EXPECT_TRUE(value2);
     EXPECT_TRUE(value3);
@@ -237,7 +237,7 @@ TEST_F(ParseCommonArgsTest, ParsesVectorArgs)
     const char *const cmdline[] = {
         "test", "-v1", "2", "1", "3", "-v2", "1"
     };
-    parseFromArray(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    parseFromArray(cmdline, 0, {}, pa);
     EXPECT_EQ(2.0, value1[XX]);
     EXPECT_EQ(1.0, value1[YY]);
     EXPECT_EQ(3.0, value1[ZZ]);
@@ -260,7 +260,7 @@ TEST_F(ParseCommonArgsTest, ParsesTimeArgs)
     const char *const cmdline[] = {
         "test", "-t1", "2", "-t2", "-.5"
     };
-    parseFromArray(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    parseFromArray(cmdline, 0, {}, pa);
     EXPECT_EQ( 2.0, value1);
     EXPECT_EQ(-0.5, value2);
     EXPECT_EQ( 2.5, value3);
@@ -277,7 +277,7 @@ TEST_F(ParseCommonArgsTest, ParsesTimeArgsWithTimeUnit)
     const char *const cmdline[] = {
         "test", "-t1", "2", "-t2", "-.5", "-tu", "ns"
     };
-    parseFromArray(cmdline, PCA_TIME_UNIT, gmx::EmptyArrayRef(), pa);
+    parseFromArray(cmdline, PCA_TIME_UNIT, {}, pa);
     EXPECT_EQ( 2000.0, value1);
     EXPECT_EQ(-500.0, value2);
     EXPECT_EQ( 2.5, value3);
@@ -296,7 +296,7 @@ TEST_F(ParseCommonArgsTest, ParsesEnumArgs)
     const char *const cmdline[] = {
         "test", "-s1", "off", "-s2", "val"
     };
-    parseFromArray(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    parseFromArray(cmdline, 0, {}, pa);
     EXPECT_STREQ("off", value1[0]);
     EXPECT_STREQ("value", value2[0]);
     EXPECT_STREQ("default", value3[0]);
@@ -320,7 +320,7 @@ TEST_F(ParseCommonArgsTest, ParsesFileArgs)
     const char *const cmdline[] = {
         "test", "-o1", "-o2", "test", "-om", "test1", "test2.xvg", "-om2"
     };
-    parseFromArray(cmdline, 0, fnm, gmx::EmptyArrayRef());
+    parseFromArray(cmdline, 0, fnm, {});
     EXPECT_STREQ("out1.xvg", opt2fn_null("-o1", nfile(), fnm));
     EXPECT_STREQ("test.xvg", opt2fn_null("-o2", nfile(), fnm));
     gmx::ArrayRef<const std::string> files = opt2fns("-om", nfile(), fnm);
@@ -345,7 +345,7 @@ TEST_F(ParseCommonArgsTest, ParsesFileArgsWithDefaults)
     const char *const cmdline[] = {
         "test"
     };
-    parseFromArray(cmdline, 0, fnm, gmx::EmptyArrayRef());
+    parseFromArray(cmdline, 0, fnm, {});
     EXPECT_STREQ("topol.tpr", ftp2fn(efTPS, nfile(), fnm));
     EXPECT_STREQ("traj.xtc", opt2fn("-f2", nfile(), fnm));
     EXPECT_EQ(nullptr, opt2fn_null("-f2", nfile(), fnm));
@@ -366,7 +366,7 @@ TEST_F(ParseCommonArgsTest, ParsesFileArgsWithDefaultFileName)
     const char *const cmdline[] = {
         "test", "-deffnm", "def", "-f2", "other", "-o"
     };
-    parseFromArray(cmdline, PCA_CAN_SET_DEFFNM, fnm, gmx::EmptyArrayRef());
+    parseFromArray(cmdline, PCA_CAN_SET_DEFFNM, fnm, {});
     EXPECT_STREQ("def.tpr", ftp2fn(efTPS, nfile(), fnm));
     EXPECT_STREQ("other.xtc", opt2fn("-f2", nfile(), fnm));
     EXPECT_STREQ("def.xtc", opt2fn("-f3", nfile(), fnm));
@@ -384,7 +384,7 @@ TEST_F(ParseCommonArgsTest, ParseFileArgsWithCustomDefaultExtension)
     const char *const cmdline[] = {
         "test", "-o2", "-o3", "test"
     };
-    parseFromArray(cmdline, PCA_CAN_SET_DEFFNM, fnm, gmx::EmptyArrayRef());
+    parseFromArray(cmdline, PCA_CAN_SET_DEFFNM, fnm, {});
     EXPECT_STREQ("conf1.gro", opt2fn("-o1", nfile(), fnm));
     EXPECT_STREQ("conf2.pdb", opt2fn("-o2", nfile(), fnm));
     EXPECT_STREQ("test.gro", opt2fn("-o3", nfile(), fnm));
@@ -408,7 +408,7 @@ TEST_F(ParseCommonArgsTest, HandlesNonExistentInputFiles)
     const char *const cmdline[] = {
         "test", "-f2", "-f3", "other", "-f4", "trj.gro", "-g2", "foo"
     };
-    parseFromArray(cmdline, PCA_DISABLE_INPUT_FILE_CHECKING, fnm, gmx::EmptyArrayRef());
+    parseFromArray(cmdline, PCA_DISABLE_INPUT_FILE_CHECKING, fnm, {});
     EXPECT_STREQ("topol.tpr", ftp2fn(efTPS, nfile(), fnm));
     EXPECT_STREQ("trj.xtc", opt2fn("-f", nfile(), fnm));
     EXPECT_STREQ("trj2.xtc", opt2fn("-f2", nfile(), fnm));
@@ -427,7 +427,7 @@ TEST_F(ParseCommonArgsTest, HandlesNonExistentOptionalInputFiles)
     const char *const cmdline[] = {
         "test"
     };
-    parseFromArray(cmdline, 0, fnm, gmx::EmptyArrayRef());
+    parseFromArray(cmdline, 0, fnm, {});
     EXPECT_STREQ("topol.tpr", ftp2fn(efTPS, nfile(), fnm));
     EXPECT_STREQ("trj.xtc", opt2fn("-f", nfile(), fnm));
 }
@@ -444,7 +444,7 @@ TEST_F(ParseCommonArgsTest, AcceptsNonExistentInputFilesIfSpecified)
     const char *const cmdline[] = {
         "test", "-c2", "-c3", "nonexistent", "-c4", "nonexistent.cpt", "-f", "nonexistent"
     };
-    parseFromArray(cmdline, 0, fnm, gmx::EmptyArrayRef());
+    parseFromArray(cmdline, 0, fnm, {});
     EXPECT_STREQ("file1.cpt", opt2fn("-c", nfile(), fnm));
     EXPECT_STREQ("file2.cpt", opt2fn("-c2", nfile(), fnm));
     EXPECT_STREQ("nonexistent.cpt", opt2fn("-c3", nfile(), fnm));
@@ -463,7 +463,7 @@ TEST_F(ParseCommonArgsTest, HandlesCompressedFiles)
     std::string       expectedG = addFileArg("-g", ".gro.Z", efFull);
     expectedF = gmx::Path::stripExtension(expectedF);
     expectedG = gmx::Path::stripExtension(expectedG);
-    parseFromArgs(0, fnm, gmx::EmptyArrayRef());
+    parseFromArgs(0, fnm, {});
     EXPECT_EQ(expectedF, opt2fn("-f", nfile(), fnm));
     EXPECT_EQ(expectedG, opt2fn("-g", nfile(), fnm));
 }
@@ -475,7 +475,7 @@ TEST_F(ParseCommonArgsTest, AcceptsUnknownTrajectoryExtension)
     };
     args_.append("test");
     std::string       expected = addFileArg("-f", ".foo", efFull);
-    parseFromArgs(0, fnm, gmx::EmptyArrayRef());
+    parseFromArgs(0, fnm, {});
     EXPECT_EQ(expected, opt2fn("-f", nfile(), fnm));
 }
 
@@ -493,7 +493,7 @@ TEST_F(ParseCommonArgsTest, CompletesExtensionFromExistingFile)
         { efTRX, "-f3", nullptr, ffREAD },
         { efTRX, "-f4", def4.c_str(), ffREAD }
     };
-    parseFromArgs(0, fnm, gmx::EmptyArrayRef());
+    parseFromArgs(0, fnm, {});
     EXPECT_EQ(expected1, opt2fn("-f1", nfile(), fnm));
     EXPECT_EQ(expected2, opt2fn("-f2", nfile(), fnm));
     EXPECT_EQ(expected3, opt2fn("-f3", nfile(), fnm));
@@ -516,7 +516,7 @@ TEST_F(ParseCommonArgsTest, CompletesExtensionFromExistingFileWithDefaultFileNam
     std::string       deffnm    = gmx::Path::stripExtension(expected3);
     args_.append("-deffnm");
     args_.append(deffnm);
-    parseFromArgs(PCA_CAN_SET_DEFFNM, fnm, gmx::EmptyArrayRef());
+    parseFromArgs(PCA_CAN_SET_DEFFNM, fnm, {});
     EXPECT_EQ(expected1, opt2fn("-f1", nfile(), fnm));
     EXPECT_EQ(expected2, opt2fn("-f2", nfile(), fnm));
     EXPECT_EQ(expected3, opt2fn("-f3", nfile(), fnm));
index 1091872a53a5bfeee588cc0e7a9f183b36b1176c..00fd425223a0feabac62d442b4bb6d2f0065185b 100644 (file)
@@ -306,17 +306,17 @@ void dd_collect_state(gmx_domdec_t *dd,
     }
     if (state_local->flags & (1 << estX))
     {
-        gmx::ArrayRef<gmx::RVec> globalXRef = state ? makeArrayRef(state->x) : gmx::EmptyArrayRef();
-        dd_collect_vec(dd, state_local, makeConstArrayRef(state_local->x), globalXRef);
+        auto globalXRef = state ? state->x : gmx::ArrayRef<gmx::RVec>();
+        dd_collect_vec(dd, state_local, state_local->x, globalXRef);
     }
     if (state_local->flags & (1 << estV))
     {
-        gmx::ArrayRef<gmx::RVec> globalVRef = state ? makeArrayRef(state->v) : gmx::EmptyArrayRef();
-        dd_collect_vec(dd, state_local, makeConstArrayRef(state_local->v), globalVRef);
+        auto globalVRef = state ? state->v : gmx::ArrayRef<gmx::RVec>();
+        dd_collect_vec(dd, state_local, state_local->v, globalVRef);
     }
     if (state_local->flags & (1 << estCGP))
     {
-        gmx::ArrayRef<gmx::RVec> globalCgpRef = state ? makeArrayRef(state->cg_p) : gmx::EmptyArrayRef();
-        dd_collect_vec(dd, state_local, makeConstArrayRef(state_local->cg_p), globalCgpRef);
+        auto globalCgpRef = state ? state->cg_p : gmx::ArrayRef<gmx::RVec>();
+        dd_collect_vec(dd, state_local, state_local->cg_p, globalCgpRef);
     }
 }
index c540117931ebed8f469e91b10b25fdbb5bdaec08..e437c5f488f9a30737c12ddfc32dda85dc8d06e6 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2018, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019, 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.
@@ -267,15 +267,15 @@ static void dd_distribute_state(gmx_domdec_t *dd,
 
     if (state_local->flags & (1 << estX))
     {
-        distributeVec(dd, DDMASTER(dd) ? makeArrayRef(state->x) : gmx::EmptyArrayRef(), state_local->x);
+        distributeVec(dd, DDMASTER(dd) ? state->x : gmx::ArrayRef<const gmx::RVec>(), state_local->x);
     }
     if (state_local->flags & (1 << estV))
     {
-        distributeVec(dd, DDMASTER(dd) ? makeArrayRef(state->v) : gmx::EmptyArrayRef(), state_local->v);
+        distributeVec(dd, DDMASTER(dd) ? state->v : gmx::ArrayRef<const gmx::RVec>(), state_local->v);
     }
     if (state_local->flags & (1 << estCGP))
     {
-        distributeVec(dd, DDMASTER(dd) ? makeArrayRef(state->cg_p) : gmx::EmptyArrayRef(), state_local->cg_p);
+        distributeVec(dd, DDMASTER(dd) ? state->cg_p : gmx::ArrayRef<const gmx::RVec>(), state_local->cg_p);
     }
 }
 
index 92ffdd85bdd0f3771e582120d89666e815935afc..aa5f7b2832b9388bc179cc158bfa568948c18501 100644 (file)
@@ -117,7 +117,7 @@ gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t *dd)
     }
     else
     {
-        return gmx::EmptyArrayRef();
+        return {};
     }
 }
 
@@ -455,7 +455,7 @@ int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
     else
     {
         // Currently unreachable
-        at2con_mt = gmx::EmptyArrayRef();
+        at2con_mt = {};
         ireq      = nullptr;
     }
 
index 23b7bce5ef8f18afd23bfb12c24ac87409bb0159..9482fc1ee5dd2bbd437c0f9a00347b568db460a8 100644 (file)
@@ -615,7 +615,7 @@ static int make_reverse_ilist(const InteractionLists &ilist,
     low_make_reverse_ilist(ilist, atoms->atom, vsitePbc,
                            count,
                            bConstr, bSettle, bBCheck,
-                           gmx::EmptyArrayRef(), gmx::EmptyArrayRef(),
+                           {}, {},
                            bLinkToAllAtoms, FALSE);
 
     ril_mt->index.push_back(0);
@@ -705,7 +705,7 @@ static gmx_reverse_top_t make_reverse_top(const gmx_mtop_t *mtop, gmx_bool bFE,
         *nint +=
             make_reverse_ilist(*mtop->intermolecular_ilist,
                                &atoms_global,
-                               gmx::EmptyArrayRef(),
+                               {},
                                rt.bConstr, rt.bSettle, rt.bBCheck, FALSE,
                                &rt.ril_intermol);
     }
@@ -2285,7 +2285,7 @@ t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
 
         make_reverse_ilist(*mtop->intermolecular_ilist,
                            &atoms,
-                           gmx::EmptyArrayRef(),
+                           {},
                            FALSE, FALSE, FALSE, TRUE, &ril_intermol);
     }
 
@@ -2313,7 +2313,7 @@ t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
          * The constraints are discarded here.
          */
         reverse_ilist_t ril;
-        make_reverse_ilist(molt.ilist, &molt.atoms, gmx::EmptyArrayRef(),
+        make_reverse_ilist(molt.ilist, &molt.atoms, {},
                            FALSE, FALSE, FALSE, TRUE, &ril);
 
         cgi_mb = &cginfo_mb[mb];
index 9fd7a138e9e359cb8c04ca2a24fe37a36634ce55..c6a7b8726a912101a7de38aa92306d8983176226 100644 (file)
@@ -617,7 +617,7 @@ int gmx_pmeonly(struct gmx_pme_t *pme,
         // from mdatoms for the other call to gmx_pme_do), so we have
         // fewer lines of code and less parameter passing.
         const int pmeFlags = GMX_PME_DO_ALL_F | (bEnerVir ? GMX_PME_CALC_ENER_VIR : 0);
-        PmeOutput output   = {gmx::EmptyArrayRef {}, 0, {{0}}, 0, {{0}}};
+        PmeOutput output   = {{}, 0, {{0}}, 0, {{0}}};
         if (useGpuForPme)
         {
             const bool boxChanged = false;
index b251dd9af954727ce8cf1dc0635fe68857f07f7e..0b11af8b86749cc046702b16c8fd0b8490b3ab19 100644 (file)
@@ -662,7 +662,7 @@ int gmx_trjcat(int argc, char *argv[])
                 else
                 {
                     trxout = trjtools_gmx_prepare_tng_writing(out_file, 'w', nullptr,
-                                                              inFilesEdited[0].c_str(), -1, nullptr, gmx::EmptyArrayRef(), nullptr);
+                                                              inFilesEdited[0].c_str(), -1, nullptr, {}, nullptr);
                 }
             }
             else
index 199f192f04236b29281667e901d6a2d191e555f1..299ac0b150e7f7ff71146490e7fc45e5622accb0 100644 (file)
@@ -540,7 +540,7 @@ int gmx_x2top(int argc, char *argv[])
     init_nnb(&nnb, atoms->nr, 4);
     gen_nnb(&nnb, plist);
     print_nnb(&nnb, "NNB");
-    gen_pad(&nnb, atoms, &rtp_header_settings, plist, excls, gmx::EmptyArrayRef(), TRUE);
+    gen_pad(&nnb, atoms, &rtp_header_settings, plist, excls, {}, TRUE);
     done_nnb(&nnb);
 
     if (!bPairs)
index c83a541da8210cdc38c2b94d93572b66d86d13cf..02269da0f1824d850afe71de2978bc7519cfa54f 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2018, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019, 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.
@@ -90,15 +90,6 @@ class ArrayRefWithPadding
          */
         ArrayRefWithPadding()
             : begin_(nullptr), end_(nullptr), paddedEnd_(nullptr) {}
-        /*! \brief
-         * Constructs an empty reference.
-         *
-         * This is provided for convenience, such that EmptyArrayRef can be
-         * used to initialize any ArrayRefWithPadding, without specifying the template
-         * type.  It is not explicit to enable that usage.
-         */
-        ArrayRefWithPadding(const EmptyArrayRef & /*unused*/)
-            : begin_(nullptr), end_(nullptr), paddedEnd_(nullptr) {}
         /*! \brief
          * Constructs a reference to a particular range.
          *
index b9b165bab004e6e138e11edb3ff26b7c6ce2f9ae..5a3907bc448fb4c90aaadb514ca630646fd17cee 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2018, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019, 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.
@@ -59,8 +59,7 @@ namespace
 
 TEST(EmptyArrayRefWithPaddingTest, IsEmpty)
 {
-    EmptyArrayRef             emptyArray = {};
-    ArrayRefWithPadding<real> empty(emptyArray);
+    ArrayRefWithPadding<real> empty;
 
     EXPECT_EQ(0U, empty.size());
     EXPECT_TRUE(empty.empty());
@@ -68,8 +67,7 @@ TEST(EmptyArrayRefWithPaddingTest, IsEmpty)
 
 TEST(EmptyConstArrayRefWithPaddingTest, IsEmpty)
 {
-    EmptyArrayRef                   emptyArray = {};
-    ArrayRefWithPadding<const real> empty(emptyArray);
+    ArrayRefWithPadding<const real> empty;
 
     EXPECT_EQ(0U, empty.size());
     EXPECT_TRUE(empty.empty());
index 785d315c99d092680e2dba575b90c468c010d755..cf66b046a7a6d6955a043bb0c4fd74329c9e4618 100644 (file)
@@ -711,7 +711,7 @@ ArrayRef<real> Constraints::rmsdData() const
     }
     else
     {
-        return EmptyArrayRef();
+        return {};
     }
 }
 
index e706c35ae9fa03505f6228266a012f88e53a1f86..225dbc1ed818be2db2d6551471d6602a38763e1b 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
@@ -265,13 +265,13 @@ void mdoutf_write_to_trajectory_files(FILE *fplog, const t_commrec *cr,
         {
             if (mdof_flags & (MDOF_X | MDOF_X_COMPRESSED))
             {
-                gmx::ArrayRef<gmx::RVec> globalXRef = MASTER(cr) ? makeArrayRef(state_global->x) : gmx::EmptyArrayRef();
-                dd_collect_vec(cr->dd, state_local, makeArrayRef(state_local->x), globalXRef);
+                auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
+                dd_collect_vec(cr->dd, state_local, state_local->x, globalXRef);
             }
             if (mdof_flags & MDOF_V)
             {
-                gmx::ArrayRef<gmx::RVec> globalVRef = MASTER(cr) ? makeArrayRef(state_global->v) : gmx::EmptyArrayRef();
-                dd_collect_vec(cr->dd, state_local, makeArrayRef(state_local->v), globalVRef);
+                auto globalVRef = MASTER(cr) ? state_global->v : gmx::ArrayRef<gmx::RVec>();
+                dd_collect_vec(cr->dd, state_local, state_local->v, globalVRef);
             }
         }
         f_global = of->f_global;
index 57cd0043fd5619a71b7bcbb8c21f8d1d33571b22..f3c3ef0688ee0cd6ec0ed4abf7bd35d6788e96c1 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
@@ -88,7 +88,7 @@ SimulationSignaller::getCommunicationBuffer()
     }
     else
     {
-        return gmx::EmptyArrayRef();
+        return {};
     }
 }
 
index 371419c54be89dd7e2f2aecb61e13b46339dc54e..88a78f329aa41037000e4e11d9c02156aa5089ec 100644 (file)
@@ -542,8 +542,8 @@ static void write_em_traj(FILE *fplog, const t_commrec *cr,
             /* If bX=true, x was collected to state_global in the call above */
             if (!bX)
             {
-                gmx::ArrayRef<gmx::RVec> globalXRef = MASTER(cr) ? makeArrayRef(state_global->x) : gmx::EmptyArrayRef();
-                dd_collect_vec(cr->dd, &state->s, makeArrayRef(state->s.x), globalXRef);
+                auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
+                dd_collect_vec(cr->dd, &state->s, state->s.x, globalXRef);
             }
         }
         else
index 52b7d11fe367fe2135d627e3f107035a4523ba2b..b666a7b83bc66ff7fcf0fe818669e75ab7bf48a3 100644 (file)
@@ -311,7 +311,7 @@ positionsFromStatePointer(const t_state *state)
     }
     else
     {
-        return gmx::EmptyArrayRef();
+        return {};
     }
 };
 
index e3afd50051d75036fd9316f6490fc8cd7636bf70..7b6d88e78a2aab54eda42db66eae4690da36dc3f 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
@@ -524,7 +524,7 @@ void NeighborhoodSearchTest::testPairSearch(
         const NeighborhoodSearchTestData &data)
 {
     testPairSearchFull(search, data, data.testPositions(), nullptr,
-                       gmx::EmptyArrayRef(), gmx::EmptyArrayRef(), false);
+                       {}, {}, false);
 }
 
 void NeighborhoodSearchTest::testPairSearchIndexed(
@@ -1049,7 +1049,7 @@ TEST_F(NeighborhoodSearchTest, SimpleSelfPairsSearch)
     ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Simple, search.mode());
 
     testPairSearchFull(&search, data, data.testPositions(), nullptr,
-                       gmx::EmptyArrayRef(), gmx::EmptyArrayRef(), true);
+                       {}, {}, true);
 }
 
 TEST_F(NeighborhoodSearchTest, GridSelfPairsSearch)
@@ -1063,7 +1063,7 @@ TEST_F(NeighborhoodSearchTest, GridSelfPairsSearch)
     ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Grid, search.mode());
 
     testPairSearchFull(&search, data, data.testPositions(), nullptr,
-                       gmx::EmptyArrayRef(), gmx::EmptyArrayRef(), true);
+                       {}, {}, true);
 }
 
 TEST_F(NeighborhoodSearchTest, HandlesConcurrentSearches)
@@ -1180,8 +1180,8 @@ TEST_F(NeighborhoodSearchTest, SimpleSearchExclusions)
 
     testPairSearchFull(&search, data,
                        data.testPositions().exclusionIds(helper.testPosIds()),
-                       helper.exclusions(), gmx::EmptyArrayRef(),
-                       gmx::EmptyArrayRef(), false);
+                       helper.exclusions(), {},
+                       {}, false);
 }
 
 TEST_F(NeighborhoodSearchTest, GridSearchExclusions)
@@ -1201,8 +1201,8 @@ TEST_F(NeighborhoodSearchTest, GridSearchExclusions)
 
     testPairSearchFull(&search, data,
                        data.testPositions().exclusionIds(helper.testPosIds()),
-                       helper.exclusions(), gmx::EmptyArrayRef(),
-                       gmx::EmptyArrayRef(), false);
+                       helper.exclusions(), {},
+                       {}, false);
 }
 
 } // namespace
index 5bb798e04e50e2eca46cdaaea086927164e15545..89348c8e5966f7f46df8708cc521e36d3989ecaa 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
@@ -87,11 +87,11 @@ class PositionCalculationTest : public ::testing::Test
 
         void testSingleStatic(e_poscalc_t type, int flags, bool bExpectTop,
                               const gmx::ArrayRef<const int> &atoms,
-                              const gmx::ArrayRef<const int> &index = gmx::EmptyArrayRef());
+                              const gmx::ArrayRef<const int> &index = {});
         void testSingleDynamic(e_poscalc_t type, int flags, bool bExpectTop,
                                const gmx::ArrayRef<const int> &initAtoms,
                                const gmx::ArrayRef<const int> &evalAtoms,
-                               const gmx::ArrayRef<const int> &index = gmx::EmptyArrayRef());
+                               const gmx::ArrayRef<const int> &index = {});
 
         gmx::test::TestReferenceData        data_;
         gmx::test::TestReferenceChecker     checker_;
index 90e67630f7809bea7f7b138ea415d20e775e9533..750c391bc35dc9f4bb18af758222a0d0b771106d 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2017,2018,2019, 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.
@@ -202,8 +202,6 @@ class SimdArrayRef
             GMX_ASSERT((reinterpret_cast<size_type>(end)/sizeof(*end))%simdWidth == 0,
                        "Size of ArrayRef needs to be divisible by type size");
         }
-        //! \copydoc ArrayRef::ArrayRef(const EmptyArrayRef&)
-        SimdArrayRef(const EmptyArrayRef & /*unused*/) : begin_(nullptr), end_(nullptr) {}
         //! \copydoc ArrayRef::ArrayRef(U)
         template<typename U,
                  typename = typename std::enable_if<
index 47529122e77a0b064853e245984e5ca919a98d0e..007cdd5ccee75fa87b1fa735bf0644470eee7091 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2015,2016,2017, by the GROMACS development team, led by
+ * Copyright (c) 2015,2016,2017,2019, 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.
@@ -78,8 +78,7 @@ namespace
 
 TEST(EmptyArrayRefTest, IsEmpty)
 {
-    EmptyArrayRef      emptyArray = {};
-    ArrayRef<SimdReal> empty(emptyArray);
+    ArrayRef<SimdReal> empty = ArrayRef<real>();
 
     EXPECT_EQ(0U, empty.size());
     EXPECT_TRUE(empty.empty());
index a22c45e257040465a45251ad685a6cce5ace8271..5cbf476c34f3ec64f3fcdfa0536656f38add1e0d 100644 (file)
@@ -169,7 +169,7 @@ ArrayRef<const RVec> TrajectoryFrame::v() const
     }
     else
     {
-        return EmptyArrayRef();
+        return {};
     }
 }
 
@@ -182,7 +182,7 @@ ArrayRef<const RVec> TrajectoryFrame::f() const
     }
     else
     {
-        return EmptyArrayRef();
+        return {};
     }
 }
 
index fdb754d06b60db74fc461c372d7e4f8b0bed8a04..ba2c041ab04ee69309c2d021511c97c95e8ba511 100644 (file)
 namespace gmx
 {
 
-/*! \brief
- * Tag type to initialize empty array references.
- *
- * This type (together with appropriate constructors in ArrayRef)
- * allows initializing any array reference to an empty value
- * without explicitly specifying its type.  This is convenient when calling
- * a function that takes an array reference, where constructing an empty
- * reference explicitly would otherwise require specifying the full array
- * reference type, including the template parameter.
- */
-struct EmptyArrayRef {};
-
 /*! \brief STL-like interface to a C array of T (or part
  * of a std container of T).
  *
@@ -142,14 +130,6 @@ class ArrayRef
          * Constructs an empty reference.
          */
         ArrayRef() : begin_(nullptr), end_(nullptr) {}
-        /*! \brief
-         * Constructs an empty reference.
-         *
-         * This is provided for convenience, such that EmptyArrayRef can be
-         * used to initialize any ArrayRef, without specifying the template
-         * type.  It is not explicit to enable that usage.
-         */
-        ArrayRef(const EmptyArrayRef & /*unused*/) : begin_(nullptr), end_(nullptr) {}
         /*! \brief
          * Constructs a reference to a container or reference
          *
index 679b4d135312263d3c195ef1c2f6683f5d074d9f..e8a8388f862b8897813af1a72ac54743aab8b4c2 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2015,2016,2017,2018,2019, 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.
@@ -57,8 +57,7 @@ namespace
 
 TEST(EmptyArrayRefTest, IsEmpty)
 {
-    EmptyArrayRef  emptyArray = {};
-    ArrayRef<real> empty(emptyArray);
+    ArrayRef<real> empty;
 
     EXPECT_EQ(0U, empty.size());
     EXPECT_TRUE(empty.empty());
@@ -66,8 +65,7 @@ TEST(EmptyArrayRefTest, IsEmpty)
 
 TEST(EmptyConstArrayRefTest, IsEmpty)
 {
-    EmptyArrayRef        emptyArray = {};
-    ArrayRef<const real> empty(emptyArray);
+    ArrayRef<const real> empty;
 
     EXPECT_EQ(0U, empty.size());
     EXPECT_TRUE(empty.empty());