Change t_extmass pointers to std::vector
authorKevin Boyd <kevin.boyd@uconn.edu>
Tue, 18 Dec 2018 03:18:16 +0000 (22:18 -0500)
committerDavid van der Spoel <spoel@xray.bmc.uu.se>
Fri, 28 Dec 2018 10:20:02 +0000 (11:20 +0100)
Using std::vector in t_extmass fixes a leak in do_md.

Used arrayRefs to view arrays in NHC_trotter

Refs #2693

Change-Id: I6bd84952fb7a470879ddd665183680c865e22fc2

src/gromacs/mdlib/coupling.cpp
src/gromacs/mdtypes/state.h

index fbfd884ad7562e04cd5cbef8cbb5137ead505565..a4b59020eb5f57420b2ecac3ad1d3eba20813041 100644 (file)
@@ -113,7 +113,6 @@ static void NHC_trotter(const t_grpopts *opts, int nvar, const gmx_ekindata_t *e
     double        dt;
     t_grp_tcstat *tcstat;
     double       *ivxi, *ixi;
-    double       *iQinv;
     double       *GQ;
     gmx_bool      bBarostat;
     int           mstepsi, mstepsj;
@@ -139,16 +138,17 @@ static void NHC_trotter(const t_grpopts *opts, int nvar, const gmx_ekindata_t *e
 
         ivxi = &vxi[i*nh];
         ixi  = &xi[i*nh];
+        gmx::ArrayRef<const double> iQinv;
         if (bBarostat)
         {
-            iQinv = &(MassQ->QPinv[i*nh]);
+            iQinv = gmx::arrayRefFromArray(&MassQ->QPinv[i*nh], nh);
             nd    = 1.0; /* THIS WILL CHANGE IF NOT ISOTROPIC */
             reft  = std::max<real>(0, opts->ref_t[0]);
             Ekin  = gmx::square(*veta)/MassQ->Winv;
         }
         else
         {
-            iQinv  = &(MassQ->Qinv[i*nh]);
+            iQinv  = gmx::arrayRefFromArray(&MassQ->Qinv[i*nh], nh);
             tcstat = &ekind->tcstat[i];
             nd     = opts->nrdf[i];
             reft   = std::max<real>(0, opts->ref_t[i]);
@@ -988,7 +988,7 @@ extern void init_npt_masses(const t_inputrec *ir, t_state *state, t_extmass *Mas
     {
         if (bInit)
         {
-            snew(MassQ->Qinv, ngtc);
+            MassQ->Qinv.resize(ngtc);
         }
         for (i = 0; (i < ngtc); i++)
         {
@@ -1034,7 +1034,7 @@ extern void init_npt_masses(const t_inputrec *ir, t_state *state, t_extmass *Mas
         /* Allocate space for thermostat variables */
         if (bInit)
         {
-            snew(MassQ->Qinv, ngtc*nh);
+            MassQ->Qinv.resize(ngtc * nh);
         }
 
         /* now, set temperature variables */
@@ -1217,7 +1217,7 @@ std::array < std::vector < int>, ettTSEQMAX> init_npt_vars(const t_inputrec *ir,
             bmass = DIM*DIM; /* recommended mass parameters for isotropic barostat */
     }
 
-    snew(MassQ->QPinv, nnhpres*opts->nhchainlength);
+    MassQ->QPinv.resize(nnhpres*opts->nhchainlength);
 
     /* barostat temperature */
     if ((ir->tau_p > 0) && (opts->ref_t[0] > 0))
index 1516daac57f390b6a9eb3ce6d77f82def066e41f..a0885cf996e29a0d25177f5251b14438756cd924 100644 (file)
@@ -233,13 +233,13 @@ class t_state
  * TODO: Move the next two structs out of state.h.
  */
 
-typedef struct t_extmass
+struct t_extmass
 {
-    double *Qinv;  /* inverse mass of thermostat -- computed from inputs, but a good place to store */
-    double *QPinv; /* inverse mass of thermostat for barostat -- computed from inputs, but a good place to store */
-    double  Winv;  /* Pressure mass inverse -- computed, not input, but a good place to store. Need to make a matrix later */
-    tensor  Winvm; /* inverse pressure mass tensor, computed       */
-} t_extmass;
+    std::vector<double> Qinv;  /* inverse mass of thermostat -- computed from inputs, but a good place to store */
+    std::vector<double> QPinv; /* inverse mass of thermostat for barostat -- computed from inputs, but a good place to store */
+    double              Winv;  /* Pressure mass inverse -- computed, not input, but a good place to store. Need to make a matrix later */
+    tensor              Winvm; /* inverse pressure mass tensor, computed       */
+};
 
 
 typedef struct