Use index instead of pointer in t_mde_delta_h_coll
authorejjordan <ejjordan@kth.se>
Tue, 20 Apr 2021 13:24:18 +0000 (15:24 +0200)
committerJoe Jordan <ejjordan12@gmail.com>
Wed, 21 Apr 2021 08:57:41 +0000 (08:57 +0000)
Several members of t_mde_delta_h_coll were just pointers into a
master array of t_mde_delta_h that are held by t_mde_delta_h_coll.
This change makes it so only an int instead of a pointer is stored.
After this it should be much easier to make t_mde_delta_h_coll
default constructible, which is planned for a followup.

Also removed 2 unused members.

A renaming of the members can also be left for a followup, once
the whole t_mde_delta_h_coll struct is modernized.

src/gromacs/mdlib/mdebin_bar.cpp
src/gromacs/mdlib/mdebin_bar.h

index a1242a039e1853ff085a474d9b3d03de8dbe1cf6..144526e60d8d09b03c4b4862df0570897615ce08 100644 (file)
@@ -443,11 +443,11 @@ void mde_delta_h_coll_init(t_mde_delta_h_coll* dhc, const t_inputrec& inputrec)
     snew(dhc->subblock_i, c_subblockINumPreEntries + dhc->n_lambda_vec);
 
     /* now decide which data to write out */
-    dhc->nlambda     = 0;
-    dhc->ndhdl       = 0;
-    dhc->dh_expanded = nullptr;
-    dhc->dh_energy   = nullptr;
-    dhc->dh_pv       = nullptr;
+    dhc->nlambda           = 0;
+    dhc->ndhdl             = 0;
+    dhc->dh_expanded_index = -1;
+    dhc->dh_energy_index   = -1;
+    dhc->dh_pv_index       = -1;
 
     /* total number of raw data point collections in the sample */
     dhc->ndh = 0;
@@ -506,7 +506,7 @@ void mde_delta_h_coll_init(t_mde_delta_h_coll* dhc, const t_inputrec& inputrec)
         n = 0;
         if (bExpanded)
         {
-            dhc->dh_expanded = dhc->dh + n;
+            dhc->dh_expanded_index = n;
             mde_delta_h_init(dhc->dh + n,
                              inputrec.fepvals->dh_hist_size,
                              inputrec.fepvals->dh_hist_spacing,
@@ -519,7 +519,7 @@ void mde_delta_h_coll_init(t_mde_delta_h_coll* dhc, const t_inputrec& inputrec)
         }
         if (bEnergy)
         {
-            dhc->dh_energy = dhc->dh + n;
+            dhc->dh_energy_index = n;
             mde_delta_h_init(dhc->dh + n,
                              inputrec.fepvals->dh_hist_size,
                              inputrec.fepvals->dh_hist_spacing,
@@ -534,7 +534,7 @@ void mde_delta_h_coll_init(t_mde_delta_h_coll* dhc, const t_inputrec& inputrec)
         n_lambda_components = 0;
         if (fep->dhdl_derivatives == DhDlDerivativeCalculation::Yes)
         {
-            dhc->dh_dhdl = dhc->dh + n;
+            dhc->dh_dhdl_index = n;
             for (auto i : keysOf(fep->separate_dvdl))
             {
                 if (inputrec.fepvals->separate_dvdl[i])
@@ -564,7 +564,7 @@ void mde_delta_h_coll_init(t_mde_delta_h_coll* dhc, const t_inputrec& inputrec)
             }
         }
         /* add the lambdas */
-        dhc->dh_du = dhc->dh + n;
+        dhc->dh_du_index = n;
         snew(lambda_vec, n_lambda_components);
         for (i = inputrec.fepvals->lambda_start_n; i < inputrec.fepvals->lambda_stop_n; i++)
         {
@@ -591,7 +591,7 @@ void mde_delta_h_coll_init(t_mde_delta_h_coll* dhc, const t_inputrec& inputrec)
         sfree(lambda_vec);
         if (bPV)
         {
-            dhc->dh_pv = dhc->dh + n;
+            dhc->dh_pv_index = n;
             mde_delta_h_init(dhc->dh + n,
                              inputrec.fepvals->dh_hist_size,
                              inputrec.fepvals->dh_hist_spacing,
@@ -642,23 +642,23 @@ void mde_delta_h_coll_add_dh(t_mde_delta_h_coll*   dhc,
 
     for (i = 0; i < dhc->ndhdl; i++)
     {
-        mde_delta_h_add_dh(dhc->dh_dhdl + i, dhdl[i]);
+        mde_delta_h_add_dh(&dhc->dh[dhc->dh_dhdl_index + i], dhdl[i]);
     }
     for (i = 0; i < dhc->nlambda; i++)
     {
-        mde_delta_h_add_dh(dhc->dh_du + i, foreign_dU[i]);
+        mde_delta_h_add_dh(&dhc->dh[dhc->dh_du_index + i], foreign_dU[i]);
     }
-    if (dhc->dh_pv != nullptr)
+    if (dhc->dh_pv_index >= 0)
     {
-        mde_delta_h_add_dh(dhc->dh_pv, pV);
+        mde_delta_h_add_dh(&dhc->dh[dhc->dh_pv_index], pV);
     }
-    if (dhc->dh_energy != nullptr)
+    if (dhc->dh_energy_index >= 0)
     {
-        mde_delta_h_add_dh(dhc->dh_energy, energy);
+        mde_delta_h_add_dh(&dhc->dh[dhc->dh_energy_index], energy);
     }
-    if (dhc->dh_expanded != nullptr)
+    if (dhc->dh_expanded_index >= 0)
     {
-        mde_delta_h_add_dh(dhc->dh_expanded, fep_state);
+        mde_delta_h_add_dh(&dhc->dh[dhc->dh_expanded_index], fep_state);
     }
 }
 
index ba56c53922d29f834c9bf3b59f366068660696ff..5964c2fe62f3eeed8759767afb92890642941d49 100644 (file)
@@ -94,16 +94,15 @@ struct t_mde_delta_h_coll
     t_mde_delta_h* dh;  /* the delta h data */
     int            ndh; /* the number of delta_h structures */
 
-    int            nlambda; /* number of bar dU delta_h structures */
-    t_mde_delta_h* dh_du;   /* the delta h data (pointer into dh) */
+    int nlambda;     /* number of bar dU delta_h structures */
+    int dh_du_index; /* the delta h data (index into dh) */
 
-    int            ndhdl;   /* number of bar dU delta_h structures */
-    t_mde_delta_h* dh_dhdl; /* the dhdl data (pointer into dh) */
+    int ndhdl;         /* number of bar dU delta_h structures */
+    int dh_dhdl_index; /* the dhdl data (index into dh) */
 
-    t_mde_delta_h* dh_energy;   /* energy output block (pointer into dh) */
-    t_mde_delta_h* dh_pv;       /* pV output block (pointer into dh) */
-    t_mde_delta_h* dh_expanded; /* expanded ensemble output block (pointer
-                                   into dh) */
+    int dh_energy_index;   /* energy output block (index into dh) */
+    int dh_pv_index;       /* pV output block (index into dh) */
+    int dh_expanded_index; /* expanded ensemble output block (index into dh) */
 
     double start_time;     /* start time of the current dh collection */
     double delta_time;     /* time difference between samples */
@@ -122,10 +121,6 @@ struct t_mde_delta_h_coll
 
     double* subblock_d; /* for writing a metadata mdebin subblock for I/O */
     int*    subblock_i; /* for writing a metadata mdebin subblock for I/O */
-
-    double* lambda_vec_subblock; /* native lambda vector data subblock for
-                                    I/O */
-    int* lambda_index_subblock;  /* lambda vector index data subblock for I/O */
 };