Fix leaked memory created in make_ljpme_c6grid.
authorSebastian Kehl <sebastian.kehl@mpcdf.mpg.de>
Fri, 12 Feb 2021 19:24:47 +0000 (19:24 +0000)
committerJoe Jordan <ejjordan12@gmail.com>
Fri, 12 Feb 2021 19:24:47 +0000 (19:24 +0000)
admin/lsan-suppressions.txt
src/gromacs/gmxlib/nonbonded/nb_free_energy.cpp
src/gromacs/mdlib/forcerec.cpp
src/gromacs/mdtypes/forcerec.h

index f04a688a8f68fac9a0ca4eef3bd40b4f05458865..de386269962237d20d0478d5555cfdfb8d90b076 100644 (file)
@@ -36,7 +36,6 @@ leak:make_bondeds_zone
 leak:make_dd_indices
 leak:make_exclusions_zone
 leak:make_fep_list
-leak:make_ljpme_c6grid
 leak:make_pull_groups
 leak:make_reverse_top
 leak:make_rotation_groups
index 3c9b933758a5b4d5f77f0812f8cf01f867e69113..c377221472d1afe5b4af9b89a2e137c73934ac22 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) 2013,2014,2015,2016,2017 by the GROMACS development team.
- * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2018,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.
@@ -57,6 +57,7 @@
 #include "gromacs/mdtypes/mdatom.h"
 #include "gromacs/simd/simd.h"
 #include "gromacs/utility/fatalerror.h"
+#include "gromacs/utility/arrayref.h"
 
 
 //! Scalar (non-SIMD) data types.
@@ -233,15 +234,16 @@ static void nb_free_energy_kernel(const t_nblist* gmx_restrict nlist,
     const int* shift  = nlist->shift;
     const int* gid    = nlist->gid;
 
-    const real* shiftvec      = fr->shift_vec[0];
-    const real* chargeA       = mdatoms->chargeA;
-    const real* chargeB       = mdatoms->chargeB;
-    real*       Vc            = kernel_data->energygrp_elec;
-    const int*  typeA         = mdatoms->typeA;
-    const int*  typeB         = mdatoms->typeB;
-    const int   ntype         = fr->ntype;
-    const real* nbfp          = fr->nbfp.data();
-    const real* nbfp_grid     = fr->ljpme_c6grid;
+    const real*               shiftvec  = fr->shift_vec[0];
+    const real*               chargeA   = mdatoms->chargeA;
+    const real*               chargeB   = mdatoms->chargeB;
+    real*                     Vc        = kernel_data->energygrp_elec;
+    const int*                typeA     = mdatoms->typeA;
+    const int*                typeB     = mdatoms->typeB;
+    const int                 ntype     = fr->ntype;
+    gmx::ArrayRef<const real> nbfp      = fr->nbfp;
+    gmx::ArrayRef<const real> nbfp_grid = fr->ljpme_c6grid;
+
     real*       Vv            = kernel_data->energygrp_vdw;
     const real  lambda_coul   = kernel_data->lambda[efptCOUL];
     const real  lambda_vdw    = kernel_data->lambda[efptVDW];
index b9e8e42abf702712ced69bce0eade06d35083db0..e17946088e439e91e455611a201e705ac84a476b 100644 (file)
@@ -155,11 +155,10 @@ static std::vector<real> mk_nbfp(const gmx_ffparams_t* idef, gmx_bool bBHAM)
     return nbfp;
 }
 
-static real* make_ljpme_c6grid(const gmx_ffparams_t* idef, t_forcerec* fr)
+static std::vector<real> make_ljpme_c6grid(const gmx_ffparams_t* idef, t_forcerec* fr)
 {
-    int   i, j, k, atnr;
-    real  c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
-    real* grid;
+    int  i, j, k, atnr;
+    real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
 
     /* For LJ-PME simulations, we correct the energies with the reciprocal space
      * inside of the cut-off. To do this the non-bonded kernels needs to have
@@ -167,7 +166,7 @@ static real* make_ljpme_c6grid(const gmx_ffparams_t* idef, t_forcerec* fr)
      */
 
     atnr = idef->atnr;
-    snew(grid, 2 * atnr * atnr);
+    std::vector<real> grid(2 * atnr * atnr, 0.0);
     for (i = k = 0; (i < atnr); i++)
     {
         for (j = 0; (j < atnr); j++, k++)
index be0bf81954b1230066f0a3089fbb2cf18ae02240..8c1d1500155f0d38b0d4a33289381d8c4799c151 100644 (file)
@@ -271,7 +271,7 @@ struct t_forcerec
     int               ntype = 0; /* Number of atom types */
     gmx_bool          bBHAM = FALSE;
     std::vector<real> nbfp;
-    real*             ljpme_c6grid = nullptr; /* C6-values used on grid in LJPME */
+    std::vector<real> ljpme_c6grid; /* C6-values used on grid in LJPME */
 
     /* Energy group pair flags */
     int* egp_flags = nullptr;