Move general functions out of typedefs.h
authorTeemu Murtola <teemu.murtola@gmail.com>
Sun, 25 May 2014 18:44:50 +0000 (21:44 +0300)
committerTeemu Murtola <teemu.murtola@gmail.com>
Sun, 25 May 2014 19:38:16 +0000 (22:38 +0300)
- Move overallocation routines to smalloc.h.  Perhaps the DD routines
  could go higher up the dependency stack, but then again, they weren't
  particularly high earlier.
- Move gmx_step_str() to cstringutil.h.  This could probably go away
  completely in favor of GMX_PRId64, but here it will force less
  dependencies on code that still uses it.

These changes allow files that use one of these functions to still get
rid of all unnecessary dependencies that typedefs.h brings.

Change-Id: Ic7973799c14f1748867965fda5bf8e33fae4fd0c

src/gromacs/gmxlib/typedefs.c
src/gromacs/legacyheaders/typedefs.h
src/gromacs/mdlib/nlistheuristics.c
src/gromacs/utility/cstringutil.c
src/gromacs/utility/cstringutil.h
src/gromacs/utility/smalloc.c
src/gromacs/utility/smalloc.h

index 390863ba91366f87c8ce9c708fb32f93da657a2a..cc8d1bd9bdd2de2ab6cb8adacaaa47ee0ff2e602 100644 (file)
@@ -43,8 +43,6 @@
 
 #include <string.h>
 
-#include "thread_mpi/threads.h"
-
 #include "macros.h"
 #include "gromacs/math/vec.h"
 #include "gromacs/pbcutil/pbc.h"
 /* The source code in this file should be thread-safe.
       Please keep it that way. */
 
-
-
-static gmx_bool            bOverAllocDD     = FALSE;
-static tMPI_Thread_mutex_t over_alloc_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
-
-
-void set_over_alloc_dd(gmx_bool set)
-{
-    tMPI_Thread_mutex_lock(&over_alloc_mutex);
-    /* we just make sure that we don't set this at the same time.
-       We don't worry too much about reading this rarely-set variable */
-    bOverAllocDD = set;
-    tMPI_Thread_mutex_unlock(&over_alloc_mutex);
-}
-
-int over_alloc_dd(int n)
-{
-    if (bOverAllocDD)
-    {
-        return OVER_ALLOC_FAC*n + 100;
-    }
-    else
-    {
-        return n;
-    }
-}
-
 int gmx_int64_to_int(gmx_int64_t step, const char *warn)
 {
     int i;
@@ -98,13 +69,6 @@ int gmx_int64_to_int(gmx_int64_t step, const char *warn)
     return i;
 }
 
-char *gmx_step_str(gmx_int64_t i, char *buf)
-{
-    sprintf(buf, "%"GMX_PRId64, i);
-
-    return buf;
-}
-
 void init_inputrec(t_inputrec *ir)
 {
     memset(ir, 0, (size_t)sizeof(*ir));
index c47f5c1fbefd071e9ce2a57e0e141e874e58ddd6..242435d4e94f7059007924d6bcf903c12b819e0f 100644 (file)
 extern "C" {
 #endif
 
-/*
- * Memory (re)allocation can be VERY slow, especially with some
- * MPI libraries that replace the standard malloc and realloc calls.
- * To avoid slow memory allocation we use over_alloc to set the memory
- * allocation size for large data blocks. Since this scales the size
- * with a factor, we use log(n) realloc calls instead of n.
- * This can reduce allocation times from minutes to seconds.
- */
-/* This factor leads to 4 realloc calls to double the array size */
-#define OVER_ALLOC_FAC 1.19
-
-void set_over_alloc_dd(gmx_bool set);
-/* Turns over allocation for variable size atoms/cg/top arrays on or off,
- * default is off.
- */
-
-int over_alloc_dd(int n);
-/* Returns n when domain decomposition over allocation is off.
- * Returns OVER_ALLOC_FAC*n + 100 when over allocation in on.
- * This is to avoid frequent reallocation
- * during domain decomposition in mdrun.
- */
-
-/* Over allocation for small data types: int, real etc. */
-#define over_alloc_small(n) (int)(OVER_ALLOC_FAC*(n) + 8000)
-
-/* Over allocation for large data types: complex structs */
-#define over_alloc_large(n) (int)(OVER_ALLOC_FAC*(n) + 1000)
-
 int gmx_int64_to_int(gmx_int64_t step, const char *warn);
 /* Convert a gmx_int64_t value to int.
  * If warn!=NULL a warning message will be written
@@ -106,15 +77,6 @@ int gmx_int64_to_int(gmx_int64_t step, const char *warn);
  * "WARNING during %s:", where warn is printed in %s.
  */
 
-#define STEPSTRSIZE 22
-
-char *gmx_step_str(gmx_int64_t i, char *buf);
-/* Prints a gmx_int64_t value in buf and returns the pointer to buf.
- * buf should be large enough to contain i: STEPSTRSIZE (22) chars.
- * When multiple gmx_int64_t values are printed in the same printf call,
- * be sure to call gmx_step_str with different buffers.
- */
-
 /* Functions to initiate and delete structures *
  * These functions are defined in gmxlib/typedefs.c
  */
index 42d81ece8400c521a9db58d50fe0b87dbefcc402..2830818d57b88a1feafa2726733131d2869f52d3 100644 (file)
 
 #include "typedefs.h"
 #include "types/nlistheuristics.h"
-#include "gromacs/utility/fatalerror.h"
+
 #include "gromacs/math/vec.h"
+#include "gromacs/utility/cstringutil.h"
+#include "gromacs/utility/fatalerror.h"
 
 void reset_nlistheuristics(gmx_nlheur_t *nlh, gmx_int64_t step)
 {
index fc280a5eb48615a347b8a34358ccb818ce7f125f..3f7fd492b110c4038f0f41bd634f0331e12bca8a 100644 (file)
@@ -586,6 +586,12 @@ str_to_int64_t(const char *str, char **endptr)
 #endif
 }
 
+char *gmx_step_str(gmx_int64_t i, char *buf)
+{
+    sprintf(buf, "%"GMX_PRId64, i);
+    return buf;
+}
+
 void parse_digits_from_plain_string(const char *digitstring, int *ndigits, int **digitlist)
 {
     int i;
index a02d2b7ec171654d27975fca4420e4b89108cbde..bfd1d6f0adbabd7d08dfbbfea369a301f0e67a4c 100644 (file)
@@ -203,6 +203,18 @@ char *wrap_lines(const char *buf, int line_width, int indent,
  */
 gmx_int64_t str_to_int64_t(const char *str, char **endptr);
 
+/** Minimum size of buffer to pass to gmx_step_str(). */
+#define STEPSTRSIZE 22
+
+/*! \brief
+ * Prints a gmx_int64_t value in buf and returns the pointer to buf.
+ *
+ * buf should be large enough to contain i: STEPSTRSIZE (22) chars.
+ * When multiple gmx_int64_t values are printed in the same printf call,
+ * be sure to call gmx_step_str with different buffers.
+ */
+char *gmx_step_str(gmx_int64_t i, char *buf);
+
 #ifdef GMX_NATIVE_WINDOWS
 #define snprintf _snprintf
 #endif
index 8e08a341072d9d848d2d44848eca7eed504f616e..41bfc63153553785980894516f6559b0ce128f4e 100644 (file)
 #include <dmalloc.h>
 #endif
 
+#include "thread_mpi/threads.h"
+
 #include "gromacs/utility/fatalerror.h"
 #ifdef PRINT_ALLOC_KB
 #include "gromacs/utility/gmxmpi.h"
 #endif
 
-#ifdef DEBUG
-#include "thread_mpi/threads.h"
+static gmx_bool            g_bOverAllocDD     = FALSE;
+static tMPI_Thread_mutex_t g_over_alloc_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
 
+#ifdef DEBUG
 static void log_action(int bMal, const char *what, const char *file, int line,
                        int nelem, int size, void *ptr)
 {
@@ -358,3 +361,24 @@ void save_free_aligned(const char *name, const char *file, int line, void *ptr)
 #endif
     }
 }
+
+void set_over_alloc_dd(gmx_bool set)
+{
+    tMPI_Thread_mutex_lock(&g_over_alloc_mutex);
+    /* we just make sure that we don't set this at the same time.
+       We don't worry too much about reading this rarely-set variable */
+    g_bOverAllocDD = set;
+    tMPI_Thread_mutex_unlock(&g_over_alloc_mutex);
+}
+
+int over_alloc_dd(int n)
+{
+    if (g_bOverAllocDD)
+    {
+        return OVER_ALLOC_FAC*n + 100;
+    }
+    else
+    {
+        return n;
+    }
+}
index 0d1aa3aaed042c91c9cd19a63b34a6d88dbbceaa..ea91947cf738910e83a4f5aa200f514272156a5c 100644 (file)
@@ -73,6 +73,8 @@
 
 #include <stddef.h>
 
+#include "basedefinitions.h"
+
 #ifdef __cplusplus
 extern "C" {
 #endif
@@ -347,4 +349,43 @@ void gmx_snew_aligned_impl(const char *name, const char *file, int line,
  */
 #define sfree_aligned(ptr) save_free_aligned(#ptr, __FILE__, __LINE__, (ptr))
 
+/*! \brief
+ * Over allocation factor for memory allocations.
+ *
+ * Memory (re)allocation can be VERY slow, especially with some
+ * MPI libraries that replace the standard malloc and realloc calls.
+ * To avoid slow memory allocation we use over_alloc to set the memory
+ * allocation size for large data blocks. Since this scales the size
+ * with a factor, we use log(n) realloc calls instead of n.
+ * This can reduce allocation times from minutes to seconds.
+ *
+ * This factor leads to 4 realloc calls to double the array size.
+ */
+#define OVER_ALLOC_FAC 1.19
+
+/*! \brief
+ * Turns over allocation for variable size atoms/cg/top arrays on or off,
+ * default is off.
+ *
+ * \todo
+ * This is mdrun-specific, so it might be better to put this and
+ * over_alloc_dd() much higher up.
+ */
+void set_over_alloc_dd(gmx_bool set);
+
+/*! \brief
+ * Returns new allocation count for domain decomposition allocations.
+ *
+ * Returns n when domain decomposition over allocation is off.
+ * Returns OVER_ALLOC_FAC*n + 100 when over allocation in on.
+ * This is to avoid frequent reallocation during domain decomposition in mdrun.
+ */
+int over_alloc_dd(int n);
+
+/** Over allocation for small data types: int, real etc. */
+#define over_alloc_small(n) (int)(OVER_ALLOC_FAC*(n) + 8000)
+
+/** Over allocation for large data types: complex structs */
+#define over_alloc_large(n) (int)(OVER_ALLOC_FAC*(n) + 1000)
+
 #endif