Change nbnxm FEP pointer to unique_ptr
authorBerk Hess <hess@kth.se>
Tue, 30 Apr 2019 07:07:24 +0000 (09:07 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Thu, 2 May 2019 14:23:08 +0000 (16:23 +0200)
This fixes an incorrect delete issue in free energy runs.

Change-Id: Ib247c6eaa98807fe3abe2fcefd09d8509d6b8598

src/gromacs/nbnxm/kerneldispatch.cpp
src/gromacs/nbnxm/pairlist.cpp
src/gromacs/nbnxm/pairlistset.h

index d9b15c8412acdf6d4d80a397d122ef22388f56ab..23be585cbc9cd41c39c241cbfe4430c1c4d4481a 100644 (file)
@@ -516,7 +516,7 @@ nonbonded_verlet_t::dispatchFreeEnergyKernel(Nbnxm::InteractionLocality  iLocali
                                              const int                   forceFlags,
                                              t_nrnb                     *nrnb)
 {
-    const gmx::ArrayRef<t_nblist const * const > nbl_fep = pairlistSets().pairlistSet(iLocality).fepLists();
+    const auto nbl_fep = pairlistSets().pairlistSet(iLocality).fepLists();
 
     /* When the first list is empty, all are empty and there is nothing to do */
     if (!pairlistSets().params().haveFep || nbl_fep[0]->nrj == 0)
@@ -558,7 +558,7 @@ nonbonded_verlet_t::dispatchFreeEnergyKernel(Nbnxm::InteractionLocality  iLocali
     {
         try
         {
-            gmx_nb_free_energy_kernel(nbl_fep[th],
+            gmx_nb_free_energy_kernel(nbl_fep[th].get(),
                                       x, f, fr, &mdatoms, &kernel_data, nrnb);
         }
         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
@@ -599,7 +599,7 @@ nonbonded_verlet_t::dispatchFreeEnergyKernel(Nbnxm::InteractionLocality  iLocali
             {
                 try
                 {
-                    gmx_nb_free_energy_kernel(nbl_fep[th],
+                    gmx_nb_free_energy_kernel(nbl_fep[th].get(),
                                               x, f, fr, &mdatoms, &kernel_data, nrnb);
                 }
                 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
index 6e0fde9908bccc65139923a742f0a4b8e5da2fec..ee890ba7a4367913cf86180efbde8d697f35d90c 100644 (file)
@@ -738,8 +738,8 @@ PairlistSet::PairlistSet(const Nbnxm::InteractionLocality  locality,
                  * master thread (but all contained list memory thread local)
                  * impacts performance.
                  */
-                snew(fepLists_[i], 1);
-                nbnxn_init_pairlist_fep(fepLists_[i]);
+                fepLists_[i] = std::make_unique<t_nblist>();
+                nbnxn_init_pairlist_fep(fepLists_[i].get());
             }
             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
         }
@@ -2740,8 +2740,8 @@ static void combine_nblists(gmx::ArrayRef<const NbnxnPairlistGpu>  nbls,
     }
 }
 
-static void balance_fep_lists(gmx::ArrayRef<t_nblist *>      fepLists,
-                              gmx::ArrayRef<PairsearchWork>  work)
+static void balance_fep_lists(gmx::ArrayRef < std::unique_ptr < t_nblist>> fepLists,
+                              gmx::ArrayRef<PairsearchWork>                work)
 {
     const int numLists = fepLists.ssize();
 
@@ -2754,7 +2754,7 @@ static void balance_fep_lists(gmx::ArrayRef<t_nblist *>      fepLists,
     /* Count the total i-lists and pairs */
     int nri_tot = 0;
     int nrj_tot = 0;
-    for (auto list : fepLists)
+    for (const auto &list : fepLists)
     {
         nri_tot += list->nri;
         nrj_tot += list->nrj;
@@ -2797,7 +2797,7 @@ static void balance_fep_lists(gmx::ArrayRef<t_nblist *>      fepLists,
     t_nblist *nbld    = work[th_dest].nbl_fep.get();
     for (int th = 0; th < numLists; th++)
     {
-        t_nblist *nbls = fepLists[th];
+        const t_nblist *nbls = fepLists[th].get();
 
         for (int i = 0; i < nbls->nri; i++)
         {
@@ -2834,9 +2834,7 @@ static void balance_fep_lists(gmx::ArrayRef<t_nblist *>      fepLists,
     /* Swap the list pointers */
     for (int th = 0; th < numLists; th++)
     {
-        t_nblist *nbl_tmp = work[th].nbl_fep.release();
-        work[th].nbl_fep.reset(fepLists[th]);
-        fepLists[th]      = nbl_tmp;
+        fepLists[th].swap(work[th].nbl_fep);
 
         if (debug)
         {
@@ -4012,7 +4010,7 @@ PairlistSet::constructPairlists(const Nbnxm::GridSet          &gridSet,
 
         if (params_.haveFep)
         {
-            clear_pairlist_fep(fepLists_[th]);
+            clear_pairlist_fep(fepLists_[th].get());
         }
     }
 
@@ -4080,7 +4078,7 @@ PairlistSet::constructPairlists(const Nbnxm::GridSet          &gridSet,
 
                     work.cycleCounter.start();
 
-                    t_nblist *fepListPtr = (fepLists_.empty() ? nullptr : fepLists_[th]);
+                    t_nblist *fepListPtr = (fepLists_.empty() ? nullptr : fepLists_[th].get());
 
                     /* Divide the i cells equally over the pairlists */
                     if (isCpuType_)
index 758d10f19660e2d505584c5249d71af465b169a9..6fbe2b17928c4ef0efe674a1723f8d55d4672715 100644 (file)
@@ -49,6 +49,8 @@
 #ifndef GMX_NBNXM_PAIRLISTSET_H
 #define GMX_NBNXM_PAIRLISTSET_H
 
+#include <memory>
+
 #include "gromacs/math/vectypes.h"
 #include "gromacs/nbnxm/pairlist.h"
 #include "gromacs/utility/basedefinitions.h"
@@ -119,7 +121,7 @@ class PairlistSet
         }
 
         //! Returns the lists of free-energy pairlists, empty when nonbonded interactions are not perturbed
-        gmx::ArrayRef<t_nblist const * const> fepLists() const
+        gmx::ArrayRef < const std::unique_ptr < t_nblist>> fepLists() const
         {
             return fepLists_;
         }
@@ -140,7 +142,7 @@ class PairlistSet
         //! Tells whether the lists is of CPU type, otherwise GPU type
         gmx_bool                       isCpuType_;
         //! Lists for perturbed interactions in simple atom-atom layout
-        std::vector<t_nblist *>        fepLists_;
+        std::vector < std::unique_ptr < t_nblist>> fepLists_;
 
     public:
         /* Pair counts for flop counting */