Add gmx_getpid()
authorTeemu Murtola <teemu.murtola@gmail.com>
Thu, 11 Sep 2014 02:59:11 +0000 (05:59 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Thu, 11 Sep 2014 14:19:48 +0000 (16:19 +0200)
This allows getting rid of all conditional includes and some conditional
compilation in two different source files.

Change-Id: I3095843a2155b087730362a4f3ef409ccf24018c

src/gromacs/gmxlib/main.cpp
src/gromacs/random/random.c
src/gromacs/utility.h
src/gromacs/utility/sysinfo.cpp [new file with mode: 0644]
src/gromacs/utility/sysinfo.h [new file with mode: 0644]

index 9a3078310dfad84821c6d5c76fc25737b80e5b87..deea126674bb394a24a63eb457d5f125ac782925 100644 (file)
 #include <string.h>
 #include <time.h>
 
-#ifdef GMX_NATIVE_WINDOWS
-#include <process.h>
-#endif
-#ifdef HAVE_SYS_TIME_H
-#include <sys/time.h>
-#endif
-#ifdef HAVE_UNISTD_H
-#include <unistd.h>
-#endif
-
 #include "gromacs/fileio/filenm.h"
 #include "gromacs/fileio/gmxfio.h"
 #include "gromacs/legacyheaders/copyrite.h"
@@ -70,6 +60,7 @@
 #include "gromacs/utility/gmxmpi.h"
 #include "gromacs/utility/programcontext.h"
 #include "gromacs/utility/smalloc.h"
+#include "gromacs/utility/sysinfo.h"
 
 /* The source code in this file should be thread-safe.
          Please keep it that way. */
@@ -243,15 +234,7 @@ void gmx_log_open(const char *lognm, const t_commrec *cr,
 
     time(&t);
 
-#ifndef NO_GETPID
-#   ifdef GMX_NATIVE_WINDOWS
-    pid = _getpid();
-#   else
-    pid = getpid();
-#   endif
-#else
-    pid = 0;
-#endif
+    pid = gmx_getpid();
 
     if (bAppendFiles)
     {
index b8295a73d72b2c1a387c7baa6a51b4ea5a6d1474..1e8390abe24ffcb323a0c2c261709d1afec06e03 100644 (file)
 #include <stdlib.h>
 #include <time.h>
 
-#ifdef HAVE_UNISTD_H
-#include <unistd.h>
-#endif
-#ifdef GMX_NATIVE_WINDOWS
-#include <process.h>
-#endif
-
 #include "external/Random123-1.08/include/Random123/threefry.h"
 
 #include "gromacs/math/utilities.h"
 #include "gromacs/random/random_gausstable.h"
+#include "gromacs/utility/sysinfo.h"
 
 #define RNG_N 624
 #define RNG_M 397
@@ -237,12 +231,8 @@ gmx_rng_make_seed(void)
     else
     {
         /* No random device available, use time-of-day and process id */
-#ifdef GMX_NATIVE_WINDOWS
-        my_pid = (long)_getpid();
-#else
-        my_pid = (long)getpid();
-#endif
-        data = (unsigned int)(((long)time(NULL)+my_pid) % (long)1000000);
+        my_pid = gmx_getpid();
+        data   = (unsigned int)(((long)time(NULL)+my_pid) % (long)1000000);
     }
     return data;
 }
index a25f70cb10d6db246b409ffa13615e27c5d9b000..c0bdd6a91cdda968c6449ae462f2023ccdb23ac4 100644 (file)
  * <H3>Other Functionality</H3>
  *
  * The header init.h declares gmx::init() and gmx::finalize() for initializing
- * and deinitializing the Gromacs library.
+ * and deinitializing the \Gromacs library.
  *
  * The header arrayref.h implements a gmx::ConstArrayRef class for exposing a
  * C array or part of a std::vector (basically, any continuous stretch of
  * The header messagestringcollector.h declares a gmx::MessageStringCollector
  * class for composing messages with context information.
  *
+ * The header sysinfo.h declares gmx_getpid() for getting the current process
+ * id.
+ *
  * The header programcontext.h declares a gmx::ProgramContextInterface that is
  * used to
  * initialize and access information about the running program, such as the
diff --git a/src/gromacs/utility/sysinfo.cpp b/src/gromacs/utility/sysinfo.cpp
new file mode 100644 (file)
index 0000000..64f9ff2
--- /dev/null
@@ -0,0 +1,62 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2014, 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.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \internal \file
+ * \brief
+ * Implements functions from sysinfo.h.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \ingroup module_utility
+ */
+#include "gmxpre.h"
+
+#include "sysinfo.h"
+
+#include "config.h"
+
+#ifdef GMX_NATIVE_WINDOWS
+#include <process.h>
+#endif
+#ifdef HAVE_UNISTD_H
+#include <unistd.h>
+#endif
+
+int gmx_getpid()
+{
+#ifdef GMX_NATIVE_WINDOWS
+    return _getpid();
+#else
+    return getpid();
+#endif
+}
diff --git a/src/gromacs/utility/sysinfo.h b/src/gromacs/utility/sysinfo.h
new file mode 100644 (file)
index 0000000..919c7b6
--- /dev/null
@@ -0,0 +1,64 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2014, 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.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \libinternal \file
+ * \brief
+ * Declares functions for obtaining information about the operating environment
+ * and the current process.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \inlibraryapi
+ * \ingroup module_utility
+ */
+#ifndef GMX_UTILITY_SYSINFO_H
+#define GMX_UTILITY_SYSINFO_H
+
+#ifdef __cplusplus
+extern "C" {
+#endif
+
+/*! \brief
+ * Returns the process ID of the current process.
+ *
+ * Does not throw.
+ *
+ * \ingroup module_utility
+ */
+int gmx_getpid();
+
+#ifdef __cplusplus
+}
+#endif
+
+#endif