Convert t_forcetable to C++
authorBerk Hess <hess@kth.se>
Wed, 17 Jun 2020 21:28:59 +0000 (23:28 +0200)
committerBerk Hess <hess@kth.se>
Thu, 18 Jun 2020 07:52:40 +0000 (09:52 +0200)
Also forward declared t_nblist in nbnxm.h to reduce dependencies.

src/gromacs/listed_forces/pairs.cpp
src/gromacs/mdlib/dispersioncorrection.cpp
src/gromacs/mdlib/wall.cpp
src/gromacs/mdtypes/nblist.h
src/gromacs/nbnxm/pairlist.cpp
src/gromacs/nbnxm/pairlist.h
src/gromacs/nbnxm/pairlistset.cpp
src/gromacs/nbnxm/pairsearch.cpp
src/gromacs/tables/forcetable.cpp

index 1f758598b63bbae9a8ae7247f5a1ee910f858405..7945d689f0a3ce6272874f850ddde84011e1bae8 100644 (file)
@@ -522,14 +522,14 @@ static real do_pairs_general(int                 ftype,
 
             fscal = free_energy_evaluate_single(
                     r2, fr->sc_r_power, fr->sc_alphacoul, fr->sc_alphavdw, fr->pairsTable->scale,
-                    fr->pairsTable->data, fr->pairsTable->stride, qq, c6, c12, qqB, c6B, c12B, LFC,
-                    LFV, DLF, lfac_coul, lfac_vdw, dlfac_coul, dlfac_vdw, fr->sc_sigma6_def,
+                    fr->pairsTable->data.data(), fr->pairsTable->stride, qq, c6, c12, qqB, c6B, c12B,
+                    LFC, LFV, DLF, lfac_coul, lfac_vdw, dlfac_coul, dlfac_vdw, fr->sc_sigma6_def,
                     fr->sc_sigma6_min, sigma2_def, sigma2_min, &velec, &vvdw, dvdl);
         }
         else
         {
             /* Evaluate tabulated interaction without free energy */
-            fscal = evaluate_single(r2, fr->pairsTable->scale, fr->pairsTable->data,
+            fscal = evaluate_single(r2, fr->pairsTable->scale, fr->pairsTable->data.data(),
                                     fr->pairsTable->stride, qq, c6, c12, &velec, &vvdw);
         }
 
index d53a7e6b383638666780781bee3e9ad6cc601911..b4871765ced15e39c3909915758cb70f15cafe28 100644 (file)
@@ -407,7 +407,7 @@ void DispersionCorrection::setInteractionParameters(InteractionParams*         i
                    "Dispersion-correction code needs a table with both repulsion and dispersion "
                    "terms");
         const real  scale  = iParams->dispersionCorrectionTable_->scale;
-        const real* vdwtab = iParams->dispersionCorrectionTable_->data;
+        const real* vdwtab = iParams->dispersionCorrectionTable_->data.data();
 
         /* Round the cut-offs to exact table values for precision */
         int ri0 = static_cast<int>(std::floor(ic.rvdw_switch * scale));
index 2944b5ec8678e921b13fa8508191f21c9ce58c1b..0b838990c54cf1e8457ceebacbaa32ca33296d74 100644 (file)
@@ -117,7 +117,7 @@ void make_wall_tables(FILE*                   fplog,
 static void tableForce(real r, const t_forcetable& tab, real Cd, real Cr, real* V, real* F)
 {
     const real  tabscale = tab.scale;
-    const real* VFtab    = tab.data;
+    const real* VFtab    = tab.data.data();
 
     real rt = r * tabscale;
     int  n0 = static_cast<int>(rt);
index 1f542d7974da109af82c96648cde6cc7b3864a95..9ac98c5aaab459d7110ff2ace87c5fde8a4870be 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, by the GROMACS development team, led by
+ * Copyright (c) 2012,2014,2015,2019,2020, 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.
@@ -37,6 +37,9 @@
 #ifndef GMX_MDTYPES_NBLIST_H
 #define GMX_MDTYPES_NBLIST_H
 
+#include <vector>
+
+#include "gromacs/utility/alignedallocator.h"
 #include "gromacs/utility/basedefinitions.h"
 #include "gromacs/utility/real.h"
 
@@ -170,15 +173,13 @@ 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 */
 
-    real  r;     /* range of the table */
-    int   n;     /* n+1 is the number of table points */
-    real  scale; /* distance (nm) between two table points */
-    real* data;  /* the actual table data */
+    real r;     /* range of the table */
+    int  n;     /* n+1 is the number of table points */
+    real scale; /* distance (nm) between two table points */
+    std::vector<real, gmx::AlignedAllocator<real>> data; /* the actual table data */
 
     /* Some information about the table layout. This can also be derived from the interpolation
      * type and the table interactions, but it is convenient to have here for sanity checks, and it
@@ -190,19 +191,6 @@ struct t_forcetable
     int stride; /* Distance to next table point (number of fp variables per table point in total) */
 };
 
-struct t_nblists
-{
-    struct t_forcetable* table_elec;
-    struct t_forcetable* table_vdw;
-    struct t_forcetable* table_elec_vdw;
-
-    /* The actual neighbor lists, short and long range, see enum above
-     * for definition of neighborlist indices.
-     */
-    struct t_nblist nlist_sr[eNL_NR];
-    struct t_nblist nlist_lr[eNL_NR];
-};
-
 struct gmx_ns_t
 {
     gmx_bool       bCGlist;
index 4365b02a0842075a86ef511d63aac1226f22c7d8..253a4cd97ae770259d47d1b839b2afaf83c0bab2 100644 (file)
@@ -54,6 +54,7 @@
 #include "gromacs/mdlib/gmx_omp_nthreads.h"
 #include "gromacs/mdtypes/group.h"
 #include "gromacs/mdtypes/md_enums.h"
+#include "gromacs/mdtypes/nblist.h"
 #include "gromacs/nbnxm/atomdata.h"
 #include "gromacs/nbnxm/gpu_data_mgmt.h"
 #include "gromacs/pbcutil/ishift.h"
index 923a40082da885e445f8ee39f81840ac256ed8af..5bdfba33a83a672f98b6ef8dff6db958ab1b5683 100644 (file)
@@ -42,7 +42,6 @@
 #include "gromacs/gpu_utils/hostallocator.h"
 #include "gromacs/math/vectypes.h"
 #include "gromacs/mdtypes/locality.h"
-#include "gromacs/mdtypes/nblist.h"
 #include "gromacs/utility/basedefinitions.h"
 #include "gromacs/utility/defaultinitializationallocator.h"
 #include "gromacs/utility/enumerationhelpers.h"
@@ -52,6 +51,7 @@
 
 struct NbnxnPairlistCpuWork;
 struct NbnxnPairlistGpuWork;
+struct t_nblist;
 
 
 //! Convenience type for vector with aligned memory
index 2ebaa0fcc4a9a9d5924898f7be9afec9525bd4e0..82c344b191f6aa73340d896bba5d4857f4fb9aa4 100644 (file)
@@ -46,6 +46,8 @@
 
 #include "pairlistset.h"
 
+#include "gromacs/mdtypes/nblist.h"
+
 #include "pairlistwork.h"
 
 PairlistSet::~PairlistSet() = default;
index d8ffbcb5aeb041274941fc4c5f8b0e509a95ca07..be981232e2182344757d52db878a3f0da8920e0b 100644 (file)
@@ -45,6 +45,7 @@
 
 #include "pairsearch.h"
 
+#include "gromacs/mdtypes/nblist.h"
 #include "gromacs/utility/smalloc.h"
 
 #include "pairlist.h"
index 5b75deb2e488057c3db33d6c66da0c798d9a5d2d..af3382855a7c63be889aa72358b4c18956266244 100644 (file)
@@ -1319,9 +1319,9 @@ t_forcetable* make_tables(FILE* out, const interaction_const_t* ic, const char*
 
     /* Each table type (e.g. coul,lj6,lj12) requires four
      * numbers per table->n+1 data points. For performance reasons we want
-     * the table data to be aligned to a 32-byte boundary.
+     * the table data to be aligned to (at least) a 32-byte boundary.
      */
-    snew_aligned(table->data, table->stride * (table->n + 1) * sizeof(real), 32);
+    table->data.resize(table->stride * (table->n + 1) * sizeof(real));
 
     for (int k = 0; (k < etiNR); k++)
     {
@@ -1367,7 +1367,7 @@ t_forcetable* make_tables(FILE* out, const interaction_const_t* ic, const char*
         }
 
         copy2table(table->n, k * table->formatsize, table->stride, td[k].x, td[k].v, td[k].f,
-                   scalefactor, table->data);
+                   scalefactor, table->data.data());
 
         done_tabledata(&(td[k]));
     }
@@ -1428,8 +1428,8 @@ makeDispersionCorrectionTable(FILE* fp, const interaction_const_t* ic, real rtab
     dispersionCorrectionTable->ninteractions = 2;
     dispersionCorrectionTable->stride =
             dispersionCorrectionTable->formatsize * dispersionCorrectionTable->ninteractions;
-    snew_aligned(dispersionCorrectionTable->data,
-                 dispersionCorrectionTable->stride * (dispersionCorrectionTable->n + 1), 32);
+    dispersionCorrectionTable->data.resize(dispersionCorrectionTable->stride
+                                           * (dispersionCorrectionTable->n + 1));
 
     for (int i = 0; i <= fullTable->n; i++)
     {
@@ -1449,14 +1449,8 @@ t_forcetable::t_forcetable(enum gmx_table_interaction interaction, enum gmx_tabl
     r(0),
     n(0),
     scale(0),
-    data(nullptr),
     formatsize(0),
     ninteractions(0),
     stride(0)
 {
 }
-
-t_forcetable::~t_forcetable()
-{
-    sfree_aligned(data);
-}