Use unique_ptr for tables
authorJoe Jordan <ejjordan12@gmail.com>
Thu, 11 Feb 2021 18:21:13 +0000 (18:21 +0000)
committerArtem Zhmurov <zhmurov@gmail.com>
Thu, 11 Feb 2021 18:21:13 +0000 (18:21 +0000)
Instead of returning a raw pointer, make_tables now returns a unique_ptr.

src/gromacs/mdlib/forcerec.cpp
src/gromacs/mdlib/wall.cpp
src/gromacs/mdtypes/forcerec.h
src/gromacs/mdtypes/nblist.h
src/gromacs/tables/forcetable.cpp
src/gromacs/tables/forcetable.h

index a0af9e6e7f2be0e60736dfcfbc4464295649403d..b9e8e42abf702712ced69bce0eade06d35083db0 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-2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2013-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.
@@ -83,6 +83,7 @@
 #include "gromacs/mdtypes/interaction_const.h"
 #include "gromacs/mdtypes/md_enums.h"
 #include "gromacs/mdtypes/multipletimestepping.h"
+#include "gromacs/mdtypes/nblist.h"
 #include "gromacs/nbnxm/nbnxm.h"
 #include "gromacs/pbcutil/ishift.h"
 #include "gromacs/pbcutil/pbc.h"
index e4d49cc83df15c90c9066e2ab03b6c4483b7f189..ecc1f44ba6b0bfa1d5691f69103b446e04ee9761 100644 (file)
@@ -4,7 +4,7 @@
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2008, 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.
@@ -78,10 +78,10 @@ void make_wall_tables(FILE*                   fplog,
         fprintf(fplog, "Reading user tables for %d energy groups with %d walls\n", negp_pp, ir->nwall);
     }
 
-    snew(fr->wall_tab, ir->nwall);
+    fr->wall_tab.resize(ir->nwall);
     for (int w = 0; w < ir->nwall; w++)
     {
-        snew(fr->wall_tab[w], negp_pp);
+        fr->wall_tab[w].resize(negp_pp);
         for (int egp = 0; egp < negp_pp; egp++)
         {
             /* If the energy group pair is excluded, we don't need a table */
index 98ad9ff1aba668fb5ef72c08fe0a20500244e7cc..be0bf81954b1230066f0a3089fbb2cf18ae02240 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.
@@ -231,7 +231,7 @@ struct t_forcerec
     gmx_bool bcoultab = FALSE;
     gmx_bool bvdwtab  = FALSE;
 
-    t_forcetable* pairsTable = nullptr; /* for 1-4 interactions, [pairs] and [pairs_nb] */
+    std::unique_ptr<t_forcetable> pairsTable; /* for 1-4 interactions, [pairs] and [pairs_nb] */
 
     /* Free energy */
     int efep = 0;
@@ -249,8 +249,8 @@ struct t_forcerec
     std::unique_ptr<nonbonded_verlet_t> nbv;
 
     /* The wall tables (if used) */
-    int             nwall    = 0;
-    t_forcetable*** wall_tab = nullptr;
+    int                                                     nwall = 0;
+    std::vector<std::vector<std::unique_ptr<t_forcetable>>> wall_tab;
 
     /* The number of atoms participating in do_force_lowlevel */
     int natoms_force = 0;
index 9ac98c5aaab459d7110ff2ace87c5fde8a4870be..45e9ce11edab0a5a0462af4eaab4a247dd62c3ab 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) 2012,2014,2015,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2012,2014,2015,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.
@@ -173,6 +173,8 @@ struct t_forcetable
 {
     t_forcetable(enum gmx_table_interaction interaction, enum gmx_table_format format);
 
+    ~t_forcetable();
+
     enum gmx_table_interaction interaction; /* Types of interactions stored in this table */
     enum gmx_table_format      format;      /* Interpolation type and data format */
 
index b8036525557d95a9debcd142c9591e25ceb82c9f..ebd24fe6bb378179217615a50f44d9890842f752 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.
@@ -1259,15 +1259,16 @@ static void set_table_type(int tabsel[], const interaction_const_t* ic, gmx_bool
     }
 }
 
-t_forcetable* make_tables(FILE* out, const interaction_const_t* ic, const char* fn, real rtab, int flags)
+std::unique_ptr<t_forcetable>
+make_tables(FILE* fp, const interaction_const_t* ic, const char* fn, real rtab, int flags)
 {
     t_tabledata* td;
     gmx_bool     b14only, useUserTable;
     int          nx0, tabsel[etiNR];
     real         scalefactor;
 
-    t_forcetable* table = new t_forcetable(GMX_TABLE_INTERACTION_ELEC_VDWREP_VDWDISP,
-                                           GMX_TABLE_FORMAT_CUBICSPLINE_YFGH);
+    auto table = std::make_unique<t_forcetable>(GMX_TABLE_INTERACTION_ELEC_VDWREP_VDWDISP,
+                                                GMX_TABLE_FORMAT_CUBICSPLINE_YFGH);
 
     b14only = ((flags & GMX_MAKETABLES_14ONLY) != 0);
 
@@ -1301,7 +1302,7 @@ t_forcetable* make_tables(FILE* out, const interaction_const_t* ic, const char*
     }
     if (useUserTable)
     {
-        read_tables(out, fn, etiNR, 0, td);
+        read_tables(fp, fn, etiNR, 0, td);
         if (rtab == 0 || (flags & GMX_MAKETABLES_14ONLY))
         {
             table->n = td[0].nx;
@@ -1354,9 +1355,9 @@ t_forcetable* make_tables(FILE* out, const interaction_const_t* ic, const char*
             init_table(table->n, nx0, scale, &(td[k]), !useUserTable);
 
             fill_table(&(td[k]), tabsel[k], ic, b14only);
-            if (out)
+            if (fp)
             {
-                fprintf(out,
+                fprintf(fp,
                         "Generated table with %d data points for %s%s.\n"
                         "Tabscale = %g points/nm\n",
                         td[k].nx,
@@ -1434,7 +1435,7 @@ makeDispersionCorrectionTable(FILE* fp, const interaction_const_t* ic, real rtab
     GMX_RELEASE_ASSERT(ic->vdwtype != evdwUSER || tabfn,
                        "With VdW user tables we need a table file name");
 
-    t_forcetable* fullTable = make_tables(fp, ic, tabfn, rtab, 0);
+    std::unique_ptr<t_forcetable> fullTable = make_tables(fp, ic, tabfn, rtab, 0);
     /* Copy the contents of the table to one that has just dispersion
      * and repulsion, to improve cache performance. We want the table
      * data to be aligned to 32-byte boundaries.
@@ -1458,7 +1459,6 @@ makeDispersionCorrectionTable(FILE* fp, const interaction_const_t* ic, real rtab
             dispersionCorrectionTable->data[8 * i + j] = fullTable->data[12 * i + 4 + j];
         }
     }
-    delete fullTable;
 
     return dispersionCorrectionTable;
 }
@@ -1474,3 +1474,5 @@ t_forcetable::t_forcetable(enum gmx_table_interaction interaction, enum gmx_tabl
     stride(0)
 {
 }
+
+t_forcetable::~t_forcetable() = default;
index a889cca9536bee87466e8e559fd0e4e0b690c434..14cbfd738412f1217fb92c65ca7ec4a9ae548c08 100644 (file)
@@ -2,7 +2,7 @@
  * This file is part of the GROMACS molecular simulation package.
  *
  * Copyright (c) 2012,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.
@@ -133,7 +133,8 @@ double v_lj_ewald_lr(double beta, double r);
  *
  * \return Pointer to inner loop table structure
  */
-t_forcetable* make_tables(FILE* fp, const interaction_const_t* ic, const char* fn, real rtab, int flags);
+std::unique_ptr<t_forcetable>
+make_tables(FILE* fp, const interaction_const_t* ic, const char* fn, real rtab, int flags);
 
 /*! \brief Return a table for bonded interactions,
  *