Clarify memory management for selection values
authorTeemu Murtola <teemu.murtola@gmail.com>
Sat, 1 Nov 2014 10:57:00 +0000 (12:57 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Fri, 28 Nov 2014 22:32:38 +0000 (23:32 +0100)
Add asserts for cases where memory would be leaked when assigning the
storage pointer for gmx_ana_selvalue_t.  Add functions for cases where
the assignment is combined with some other operation that takes care of
the memory deallocation or ownership transfer to make the asserts not
give false positives.

The above clarifies the memory management and makes it less prone to
errors, but does not fix the memory leak that was revealed by "z of ..."
implementation.  To actually fix that, additionally split the logic for
freeing memory for gmx_ana_selparam_t, and call that also from the place
that earlier set the value pointer to NULL.

Move Doxygen documentation for some of the affected methods and their
neighbors to the header per newer guideline.

Add some extra debug output.

Change-Id: I4e5dfa1248a20cab8a242be7209b5b0779204e64

src/gromacs/selection/evaluate.cpp
src/gromacs/selection/selelem.cpp
src/gromacs/selection/selelem.h
src/gromacs/selection/selvalue.cpp
src/gromacs/selection/selvalue.h

index cd98dd5d9d32fdaa3e685ad66ba0691df32074d4..db2fd2cc666c7aeca092e7c65be588211f1e2673 100644 (file)
@@ -240,8 +240,7 @@ class SelelemTemporaryValueAssigner
                                "Can only assign one element with one instance");
             GMX_RELEASE_ASSERT(sel->v.type == vsource.v.type,
                                "Mismatching selection value types");
-            old_ptr_    = sel->v.u.ptr;
-            old_nalloc_ = sel->v.nalloc;
+            _gmx_selvalue_getstore_and_release(&sel->v, &old_ptr_, &old_nalloc_);
             _gmx_selvalue_setstore(&sel->v, vsource.v.u.ptr);
             sel_ = sel;
         }
index 76f65be97163fdfe2d6db303ef3eda0297ed12fc..b3b65409799ca6b4d6dacebee73aa82149bb868a 100644 (file)
@@ -195,21 +195,13 @@ void SelectionTreeElement::freeValues()
                 break;
         }
     }
-    if (v.nalloc > 0)
+    _gmx_selvalue_free(&v);
+    if (type == SEL_SUBEXPRREF && u.param != NULL)
     {
-        if (v.type == POS_VALUE)
-        {
-            delete [] v.u.p;
-        }
-        else
-        {
-            sfree(v.u.ptr);
-        }
-    }
-    _gmx_selvalue_setstore(&v, NULL);
-    if (type == SEL_SUBEXPRREF && u.param)
-    {
-        u.param->val.u.ptr = NULL;
+        // TODO: This is now called from two different locations.
+        // It is likely that one of them is unnecessary, but that requires
+        // extra analysis to clarify.
+        _gmx_selelem_free_param(u.param);
     }
 }
 
@@ -472,10 +464,22 @@ _gmx_selelem_set_vtype(const gmx::SelectionTreeElementPointer &sel,
     }
 }
 
-/*!
- * \param[in] method Method to free.
- * \param[in] mdata  Method data to free.
- */
+void
+_gmx_selelem_free_param(gmx_ana_selparam_t *param)
+{
+    if (param->val.u.ptr != NULL)
+    {
+        if (param->val.type == GROUP_VALUE)
+        {
+            for (int i = 0; i < param->val.nr; ++i)
+            {
+                gmx_ana_index_deinit(&param->val.u.g[i]);
+            }
+        }
+        _gmx_selvalue_free(&param->val);
+    }
+}
+
 void
 _gmx_selelem_free_method(gmx_ana_selmethod_t *method, void *mdata)
 {
@@ -493,37 +497,11 @@ _gmx_selelem_free_method(gmx_ana_selmethod_t *method, void *mdata)
      * access them. */
     if (method)
     {
-        int  i, j;
-
         /* Free the memory allocated for the parameters that are not managed
          * by the selection method itself. */
-        for (i = 0; i < method->nparams; ++i)
+        for (int i = 0; i < method->nparams; ++i)
         {
-            gmx_ana_selparam_t *param = &method->param[i];
-
-            if (param->val.u.ptr)
-            {
-                if (param->val.type == GROUP_VALUE)
-                {
-                    for (j = 0; j < param->val.nr; ++j)
-                    {
-                        gmx_ana_index_deinit(&param->val.u.g[j]);
-                    }
-                }
-                else if (param->val.type == POS_VALUE)
-                {
-                    if (param->val.nalloc > 0)
-                    {
-                        delete[] param->val.u.p;
-                        _gmx_selvalue_setstore(&param->val, NULL);
-                    }
-                }
-
-                if (param->val.nalloc > 0)
-                {
-                    sfree(param->val.u.ptr);
-                }
-            }
+            _gmx_selelem_free_param(&method->param[i]);
         }
         sfree(method->param);
         sfree(method);
@@ -632,7 +610,7 @@ _gmx_selelem_print_tree(FILE *fp, const gmx::SelectionTreeElement &sel,
         fprintf(fp, " eval=");
         _gmx_sel_print_evalfunc_name(fp, sel.evaluate);
     }
-    if (!(sel.flags & SEL_ALLOCVAL))
+    if (sel.v.nalloc < 0)
     {
         fprintf(fp, " (ext)");
     }
@@ -674,6 +652,23 @@ _gmx_selelem_print_tree(FILE *fp, const gmx::SelectionTreeElement &sel,
             fprintf(fp, "\n");
         }
     }
+    else if (sel.type == SEL_SUBEXPRREF && sel.u.param != NULL)
+    {
+        fprintf(fp, "%*c param", level*2+1, ' ');
+        if (sel.u.param->name != NULL)
+        {
+            fprintf(fp, " \"%s\"", sel.u.param->name);
+        }
+        if (sel.u.param->val.nalloc < 0)
+        {
+            fprintf(fp, " (ext)");
+        }
+        else
+        {
+            fprintf(fp, " nalloc: %d", sel.u.param->val.nalloc);
+        }
+        fprintf(fp, "\n");
+    }
 
     if (sel.cdata)
     {
index db1882f5a756ef0ffc26a03f8300341a2173d28f..5fa419cab9abdd9aabc78b2401b853d253a0fa35 100644 (file)
@@ -465,7 +465,19 @@ void
 _gmx_selelem_set_vtype(const gmx::SelectionTreeElementPointer &sel,
                        e_selvalue_t                            vtype);
 
-/** Frees the memory allocated for a selection method. */
+/*! \brief
+ * Frees the memory allocated for a selection method parameter.
+ *
+ * \param[in] param Parameter to free.
+ */
+void
+_gmx_selelem_free_param(struct gmx_ana_selparam_t *param);
+/*! \brief
+ * Frees the memory allocated for a selection method.
+ *
+ * \param[in] method Method to free.
+ * \param[in] mdata  Method data to free.
+ */
 void
 _gmx_selelem_free_method(struct gmx_ana_selmethod_t *method, void *mdata);
 
index 86edfba070fa72c3cc207e17c3c6094c53aed62c..80e1cfea57b57615367f5e93d218ee69c2ed5594 100644 (file)
 #include "gromacs/utility/gmxassert.h"
 #include "gromacs/utility/smalloc.h"
 
-/*!
- * \param[out] val  Output structure
- *
- * The type of \p val is not touched.
- * Any contents of \p val are discarded without freeing.
- */
 void
 _gmx_selvalue_clear(gmx_ana_selvalue_t *val)
 {
@@ -62,28 +56,33 @@ _gmx_selvalue_clear(gmx_ana_selvalue_t *val)
     val->nalloc = 0;
 }
 
-/*!
- * \param[in,out] val  Value structure to allocate.
- * \param[in]     n    Maximum number of values needed.
- * \returns       Zero on success.
- *
- * Reserves memory for the values within \p val to store at least \p n values,
- * of the type specified in the \p val structure.
- *
- * If the type is \ref POS_VALUE or \ref GROUP_VALUE, memory is reserved for
- * the data structures, but no memory is reserved inside these newly allocated
- * data structures.
- * Similarly, for \ref STR_VALUE values, the pointers are set to NULL.
- * For other values, the memory is uninitialized.
- */
-int
+void
+_gmx_selvalue_free(gmx_ana_selvalue_t *val)
+{
+    if (val->nalloc > 0)
+    {
+        if (val->type == POS_VALUE)
+        {
+            delete[] val->u.p;
+        }
+        else
+        {
+            sfree(val->u.ptr);
+        }
+    }
+    // TODO: It causes a memory leak somewhere if val->nr is assigned zero here...
+    val->u.ptr  = NULL;
+    val->nalloc = 0;
+}
+
+void
 _gmx_selvalue_reserve(gmx_ana_selvalue_t *val, int n)
 {
     int  i;
 
     if (val->nalloc == -1)
     {
-        return 0;
+        return;
     }
 
     if (!val->u.ptr || val->nalloc < n)
@@ -115,34 +114,31 @@ _gmx_selvalue_reserve(gmx_ana_selvalue_t *val, int n)
         }
         val->nalloc = n;
     }
-    return 0;
 }
 
-/*!
- * \param[in,out] val    Value structure to allocate.
- * \param[in]     ptr    Pointer where the values should be stored.
- * \returns       Zero on success.
- *
- * Automatic memory management is disabled for \p ptr, unless \p ptr is NULL.
- */
-int
+void
+_gmx_selvalue_getstore_and_release(gmx_ana_selvalue_t *val, void **ptr, int *nalloc)
+{
+    *ptr        = val->u.ptr;
+    *nalloc     = val->nalloc;
+    val->u.ptr  = NULL;
+    val->nalloc = 0;
+}
+
+void
 _gmx_selvalue_setstore(gmx_ana_selvalue_t *val, void *ptr)
 {
+    GMX_ASSERT(val->nalloc <= 0,
+               "Memory leak from discarding an existing value");
     val->u.ptr  = ptr;
     val->nalloc = (ptr ? -1 : 0);
-    return 0;
 }
 
-/*!
- * \param[in,out] val    Value structure to allocate.
- * \param[in]     ptr    Pointer where the values should be stored.
- * \param[in]     nalloc Number of values allocated for \p ptr.
- * \returns       Zero on success.
- */
-int
+void
 _gmx_selvalue_setstore_alloc(gmx_ana_selvalue_t *val, void *ptr, int nalloc)
 {
+    GMX_ASSERT(val->nalloc <= 0 || (ptr == val->u.ptr && nalloc == val->nalloc),
+               "Memory leak from discarding an existing value");
     val->u.ptr  = ptr;
     val->nalloc = nalloc;
-    return 0;
 }
index 63f260b2b417c296027b52723eb8e5f5ba448c42..d4a1e13318cc6a59d39572877e20239b4c94b74a 100644 (file)
@@ -105,17 +105,79 @@ typedef struct gmx_ana_selvalue_t
     int                         nalloc;
 } gmx_ana_selvalue_t;
 
-/** Initializes an empty selection value structure. */
+/*! \brief
+ * Initializes an empty selection value structure.
+ *
+ * \param[out] val  Output structure
+ *
+ * The type of \p val is not touched.
+ * Any contents of \p val are discarded without freeing.
+ */
 void
 _gmx_selvalue_clear(gmx_ana_selvalue_t *val);
-/** Reserve memory for storing selection values. */
-int
+/*! \brief
+ * Frees memory allocated for a selection value structure.
+ *
+ * \param[in,out] val  Values to free.
+ *
+ * The type of \p val is not touched.
+ * If memory is not allocated, the value pointers are simply cleared without
+ * freeing.
+ */
+void
+_gmx_selvalue_free(gmx_ana_selvalue_t *val);
+/*! \brief
+ * Reserve memory for storing selection values.
+ *
+ * \param[in,out] val  Value structure to allocate.
+ * \param[in]     n    Maximum number of values needed.
+ *
+ * Reserves memory for the values within \p val to store at least \p n values,
+ * of the type specified in the \p val structure.
+ *
+ * If the type is \ref POS_VALUE or \ref GROUP_VALUE, memory is reserved for
+ * the data structures, but no memory is reserved inside these newly allocated
+ * data structures.
+ * Similarly, for \ref STR_VALUE values, the pointers are set to NULL.
+ * For other values, the memory is uninitialized.
+ */
+void
 _gmx_selvalue_reserve(gmx_ana_selvalue_t *val, int n);
-/** Sets the memory for storing selection values. */
-int
+/*! \brief
+ * Gets and releases the memory pointer for storing selection values.
+ *
+ * \param[in,out] val    Value structure to release.
+ * \param[out]    ptr    Pointer where the values are stored.
+ * \param[out]    nalloc Pointer where the values are stored.
+ *
+ * Returns the pointer where values of \p val were stored in \p ptr and
+ * \p nalloc, and clears the memory in \p val.
+ */
+void
+_gmx_selvalue_getstore_and_release(gmx_ana_selvalue_t *val, void **ptr, int *nalloc);
+/*! \brief
+ * Sets the memory for storing selection values.
+ *
+ * \param[in,out] val    Value structure to set storage for.
+ * \param[in]     ptr    Pointer where the values should be stored.
+ *
+ * Automatic memory management is disabled for \p ptr.
+ * Asserts if \p val had a previous storage that it owned, as that would result
+ * in a memory leak.
+ */
+void
 _gmx_selvalue_setstore(gmx_ana_selvalue_t *val, void *ptr);
-/** Sets the memory for storing selection values and marks it for automatic freeing. */
-int
+/*! \brief
+ * Sets the memory for storing selection values and marks it for automatic freeing.
+ *
+ * \param[in,out] val    Value structure to set storage for.
+ * \param[in]     ptr    Pointer where the values should be stored.
+ * \param[in]     nalloc Number of values allocated for \p ptr.
+ *
+ * Asserts if \p val had a previous storage that it owned, as that would result
+ * in a memory leak.
+ */
+void
 _gmx_selvalue_setstore_alloc(gmx_ana_selvalue_t *val, void *ptr, int nalloc);
 
 #endif