Make using C++ type with snew() a compiler error
authorTeemu Murtola <teemu.murtola@gmail.com>
Sun, 4 Oct 2015 06:33:01 +0000 (09:33 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Mon, 5 Oct 2015 12:32:32 +0000 (14:32 +0200)
Make it impossible to use snew() and friends on a C++ type (now gives a
compiler error).  This is possible with C++11, and makes it much less
error-prone to start introducing C++ code into existing C structs.

As a side effect, this enforces also proper encapsulation of memory
management for opaque structs.  Fixed a few instances where this wasn't
the case.

CUDA code skips these assertions.

Change-Id: I470dc077bb82c6ca7da5fb9a2daa64211b6d632c

src/gromacs/gmxana/gmx_anadock.cpp
src/gromacs/gmxana/gmx_dipoles.cpp
src/gromacs/gmxana/gmx_msd.cpp
src/gromacs/statistics/statistics.cpp
src/gromacs/statistics/statistics.h
src/gromacs/utility/basedefinitions.h
src/gromacs/utility/smalloc.h
src/programs/mdrun/membed.cpp
src/programs/mdrun/membed.h
src/programs/mdrun/runner.cpp

index 9a20466442b12e8a070047e0b6b919d033b6c495..f67f55a03ed96d86ec1f1073ae1f6e3d8f10d010 100644 (file)
@@ -212,10 +212,8 @@ static void clust_stat(FILE *fp, int start, int end, t_pdbfile *pdbf[])
     fprintf(fp, "  <%12s> = %8.3f (+/- %6.3f)\n", etitles[FALSE], aver, sigma);
     gmx_stats_get_ase(ef, &aver, &sigma, NULL);
     fprintf(fp, "  <%12s> = %8.3f (+/- %6.3f)\n", etitles[TRUE], aver, sigma);
-    gmx_stats_done(ed);
-    gmx_stats_done(ef);
-    sfree(ed);
-    sfree(ef);
+    gmx_stats_free(ed);
+    gmx_stats_free(ef);
 }
 
 static real rmsd_dist(t_pdbfile *pa, t_pdbfile *pb, gmx_bool bRMSD)
index 53c23a8418d8e1d1b88c082e24a7fb5b69bf9abf..b4313049bc304bdf75997079ed1f0c7903248de3 100644 (file)
@@ -1316,7 +1316,7 @@ static void do_dip(t_topology *top, int ePBC, real volume,
                 fprintf(outeps, "%10g  %12.8e\n", t, epsilon);
             }
         }
-        gmx_stats_done(muframelsq);
+        gmx_stats_free(muframelsq);
 
         if (bMU)
         {
index 361c52b113aa13cab6ff66d8b1e06cf8f1e73ac8..4edcb51b26f2c840fb8af8ff98a83487ff31d072 100644 (file)
@@ -583,8 +583,7 @@ void printmol(t_corr *curr, const char *fn,
             }
         }
         gmx_stats_get_ab(lsq1, elsqWEIGHT_NONE, &a, &b, NULL, NULL, NULL, NULL);
-        gmx_stats_done(lsq1);
-        sfree(lsq1);
+        gmx_stats_free(lsq1);
         D     = a*FACTOR/curr->dim_factor;
         if (D < 0)
         {
index 623a47a05a02501603e11138f3d6af70389cffc5..434d55fb6d047cb01ac676b7e714a731f2709a33 100644 (file)
@@ -76,20 +76,15 @@ int gmx_stats_get_npoints(gmx_stats_t gstats, int *N)
     return estatsOK;
 }
 
-int gmx_stats_done(gmx_stats_t gstats)
+void gmx_stats_free(gmx_stats_t gstats)
 {
     gmx_stats *stats = (gmx_stats *) gstats;
 
     sfree(stats->x);
-    stats->x = NULL;
     sfree(stats->y);
-    stats->y = NULL;
     sfree(stats->dx);
-    stats->dx = NULL;
     sfree(stats->dy);
-    stats->dy = NULL;
-
-    return estatsOK;
+    sfree(stats);
 }
 
 int gmx_stats_add_point(gmx_stats_t gstats, double x, double y,
@@ -752,11 +747,7 @@ int lsq_y_ax_b_error(int n, real x[], real y[], real dy[],
     {
         return ok;
     }
-    if ((ok = gmx_stats_done(lsq)) != estatsOK)
-    {
-        return ok;
-    }
-    sfree(lsq);
+    gmx_stats_free(lsq);
 
     return estatsOK;
 }
index c09881398a1d43c7e074fb198e4ae043515c46fe..d50bdf14e89116ac2daf5fbcd26855ca7ffc1855 100644 (file)
@@ -77,9 +77,8 @@ gmx_stats_t gmx_stats_init();
 /*! \brief
  * Destroy a data structure
  * \param stats The data structure
- * \return error code
  */
-int gmx_stats_done(gmx_stats_t stats);
+void gmx_stats_free(gmx_stats_t stats);
 
 /*! \brief
  * Remove outliers from a straight line, where level in units of
index 740bbdc234e241f6ae32c2f393c82a47e7663a5d..1d2528e273bfb61493cbbd7c8084326bc23896f9 100644 (file)
@@ -171,6 +171,28 @@ typedef uint64_t gmx_uint64_t;
 #define gmx_cxx_const
 #endif
 
+/*! \def GMX_CXX11_COMPILATION
+ * \brief
+ * Defined to 1 when compiling as C++11.
+ *
+ * While \Gromacs only supports C++11 compilation, there are some parts of the
+ * code that are compiled with other tools than the actual C++ compiler, and
+ * these may not support C++11.  Most notable such case is all of CUDA code
+ * (with CUDA versions older than 6.5), but other types of kernels might also
+ * have similar limitations in the future.
+ *
+ * The define is intended for conditional compilation in low-level headers that
+ * need to support inclusion from such non-C++11 files, but get significant
+ * benefit (e.g., for correctness checking or more convenient use) from C++11.
+ * It should only be used for features that do not influence the ABI of the
+ * header; e.g., static_asserts or additional helper methods.
+ */
+#if defined __cplusplus && __cplusplus >= 201103L
+#    define GMX_CXX11_COMPILATION 1
+#else
+#    define GMX_CXX11_COMPILATION 0
+#endif
+
 /*! \def gmx_unused
  * \brief
  * Attribute to suppress compiler warnings about unused function parameters.
index 74d71a558642b28d2863554eec754883e0e0d723..465c1a4130e52cac85078f002b2226feb7d8c1ba 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,2014, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015, 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.
@@ -195,6 +195,11 @@ void save_free_aligned(const char *name, const char *file, int line, void *ptr);
 #endif
 
 #ifdef __cplusplus
+
+#if GMX_CXX11_COMPILATION
+#include <type_traits>
+#endif
+
 /*! \cond internal */
 /*! \name Implementation templates for C++ memory allocation macros
  *
@@ -222,6 +227,9 @@ template <typename T> static inline
 void gmx_snew_impl(const char *name, const char *file, int line,
                    T * &ptr, size_t nelem)
 {
+#if GMX_CXX11_COMPILATION
+    static_assert(std::is_pod<T>::value, "snew() called on C++ type");
+#endif
     ptr = (T *)save_calloc(name, file, line, nelem, sizeof(T));
 }
 /** C++ helper for srenew(). */
@@ -229,6 +237,9 @@ template <typename T> static inline
 void gmx_srenew_impl(const char *name, const char *file, int line,
                      T * &ptr, size_t nelem)
 {
+#if GMX_CXX11_COMPILATION
+    static_assert(std::is_pod<T>::value, "srenew() called on C++ type");
+#endif
     ptr = (T *)save_realloc(name, file, line, ptr, nelem, sizeof(T));
 }
 /** C++ helper for smalloc(). */
@@ -236,6 +247,9 @@ template <typename T> static inline
 void gmx_smalloc_impl(const char *name, const char *file, int line,
                       T * &ptr, size_t size)
 {
+#if GMX_CXX11_COMPILATION
+    static_assert(std::is_pod<T>::value, "smalloc() called on C++ type");
+#endif
     ptr = (T *)save_malloc(name, file, line, size);
 }
 /** C++ helper for snew_aligned(). */
@@ -243,9 +257,32 @@ template <typename T> static inline
 void gmx_snew_aligned_impl(const char *name, const char *file, int line,
                            T * &ptr, size_t nelem, size_t alignment)
 {
+#if GMX_CXX11_COMPILATION
+    static_assert(std::is_pod<T>::value, "snew_aligned() called on C++ type");
+#endif
     ptr = (T *)save_calloc_aligned(name, file, line, nelem, sizeof(T), alignment);
 }
-/*! \] */
+/** C++ helper for sfree(). */
+template <typename T> static inline
+void gmx_sfree_impl(const char *name, const char *file, int line, T *ptr)
+{
+#if GMX_CXX11_COMPILATION
+    static_assert(std::is_pod<T>::value || std::is_void<T>::value,
+                  "sfree() called on C++ type");
+#endif
+    save_free(name, file, line, ptr);
+}
+/** C++ helper for sfree_aligned(). */
+template <typename T> static inline
+void gmx_sfree_aligned_impl(const char *name, const char *file, int line, T *ptr)
+{
+#if GMX_CXX11_COMPILATION
+    static_assert(std::is_pod<T>::value || std::is_void<T>::value,
+                  "sfree_aligned() called on C++ type");
+#endif
+    save_free_aligned(name, file, line, ptr);
+}
+/*! \} */
 /*! \endcond */
 #endif /* __cplusplus */
 
@@ -304,6 +341,23 @@ void gmx_snew_aligned_impl(const char *name, const char *file, int line,
  *
  * \hideinitializer
  */
+/*! \def sfree
+ * \brief
+ * Frees memory referenced by \p ptr.
+ *
+ * \p ptr is allowed to be NULL, in which case nothing is done.
+ *
+ * \hideinitializer
+ */
+/*! \def sfree_aligned
+ * \brief
+ * Frees aligned memory referenced by \p ptr.
+ *
+ * This must only be called with a pointer obtained through snew_aligned().
+ * \p ptr is allowed to be NULL, in which case nothing is done.
+ *
+ * \hideinitializer
+ */
 #ifdef __cplusplus
 
 /* C++ implementation */
@@ -315,6 +369,10 @@ void gmx_snew_aligned_impl(const char *name, const char *file, int line,
     gmx_smalloc_impl(#ptr, __FILE__, __LINE__, (ptr), (size))
 #define snew_aligned(ptr, nelem, alignment) \
     gmx_snew_aligned_impl(#ptr, __FILE__, __LINE__, (ptr), (nelem), alignment)
+#define sfree(ptr) \
+    gmx_sfree_impl(#ptr, __FILE__, __LINE__, (ptr))
+#define sfree_aligned(ptr) \
+    gmx_sfree_aligned_impl(#ptr, __FILE__, __LINE__, (ptr))
 
 #else
 
@@ -327,6 +385,8 @@ void gmx_snew_aligned_impl(const char *name, const char *file, int line,
     (ptr) = save_malloc(#ptr, __FILE__, __LINE__, size)
 #define snew_aligned(ptr, nelem, alignment) \
     (ptr) = save_calloc_aligned(#ptr, __FILE__, __LINE__, (nelem), sizeof(*(ptr)), alignment)
+#define sfree(ptr) save_free(#ptr, __FILE__, __LINE__, (ptr))
+#define sfree_aligned(ptr) save_free_aligned(#ptr, __FILE__, __LINE__, (ptr))
 
 #endif
 
@@ -334,25 +394,6 @@ void gmx_snew_aligned_impl(const char *name, const char *file, int line,
 extern "C" {
 #endif
 
-/*! \brief
- * Frees memory referenced by \p ptr.
- *
- * \p ptr is allowed to be NULL, in which case nothing is done.
- *
- * \hideinitializer
- */
-#define sfree(ptr) save_free(#ptr, __FILE__, __LINE__, (ptr))
-
-/*! \brief
- * Frees aligned memory referenced by \p ptr.
- *
- * This must only be called with a pointer obtained through snew_aligned().
- * \p ptr is allowed to be NULL, in which case nothing is done.
- *
- * \hideinitializer
- */
-#define sfree_aligned(ptr) save_free_aligned(#ptr, __FILE__, __LINE__, (ptr))
-
 /*! \brief
  * Over allocation factor for memory allocations.
  *
index 947581a57bc582a1365d56ceec7b05abc2d00b2a..63b8bbaded9a41402024e2d7a580f1495bc011ee 100644 (file)
@@ -1336,3 +1336,8 @@ gmx_membed_t *init_membed(FILE *fplog, int nfile, const t_filenm fnm[], gmx_mtop
 
     return membed;
 }
+
+void free_membed(gmx_membed_t *membed)
+{
+    sfree(membed);
+}
index 8b59520367e789a47a2e6dc9d93015add996b411..4a3e1e431f79d0b520a50c4ee6a4c46cf563d081 100644 (file)
@@ -49,4 +49,6 @@ gmx_membed_t *init_membed(FILE *fplog, int nfile, const t_filenm fnm[], gmx_mtop
 /* rescaling the coordinates voor de membed code */
 void rescale_membed(int step_rel, gmx_membed_t *membed, rvec *x);
 
+void free_membed(gmx_membed_t *membed);
+
 #endif
index 23d087fb3be33df7c2f64754be10a45efc48521a..bbc19ca24b9a9cb21b3df194ae30eea871e6a27b 100644 (file)
@@ -1381,9 +1381,9 @@ int mdrunner(gmx_hw_opt_t *hw_opt,
     /* Free GPU memory and context */
     free_gpu_resources(fr, cr, &hwinfo->gpu_info, fr ? fr->gpu_opt : NULL);
 
-    if (opt2bSet("-membed", nfile, fnm))
+    if (membed != nullptr)
     {
-        sfree(membed);
+        free_membed(membed);
     }
 
     gmx_hardware_info_free(hwinfo);