introduce gmx_omp wrapper for the OpenMP API
authorSzilard Pall <pszilard@cbr.su.se>
Fri, 15 Jun 2012 15:03:41 +0000 (17:03 +0200)
committerSzilard Pall <pszilard@cbr.su.se>
Mon, 25 Jun 2012 17:03:00 +0000 (19:03 +0200)
The gmx_omp_* OpenMP API wrappers should be used instead of directly
calling OpenMP functions. Therefore, gmx_omp.h should be included
instead of omp.h so OpenMP API functions don't need #ifdef-ing.

Additionally, moved GMX_OPENMP from command-line define to config.h.

Change-Id: If30f433946fef908c26a29f1bff3671580e90629

CMakeLists.txt
include/gmx_omp.h [new file with mode: 0644]
src/config.h.cmakein
src/gmxlib/gmx_omp.c [new file with mode: 0644]
src/kernel/runner.c
src/mdlib/fft5d.c
src/mdlib/pme.c
src/tools/geminate.c
src/tools/gmx_hbond.c
src/tools/gmx_sans.c
src/tools/nsfactor.c

index f1ed2488e2fd9d39f8e21b9b0cc11563161d798f..9760d3a83a0a3755d67ab44192c23e76460e86eb 100644 (file)
@@ -185,7 +185,6 @@ if(GMX_OPENMP)
     if(OPENMP_FOUND)
         set(GROMACS_C_FLAGS "${OpenMP_C_FLAGS} ${GROMACS_C_FLAGS}")
         set(GROMACS_CXX_FLAGS "${OpenMP_CXX_FLAGS} ${GROMACS_CXX_FLAGS}")
-        add_definitions(-DGMX_OPENMP)
     else(OPENMP_FOUND)
         message(WARNING
                 "Compiler not supporting OpenMP. This might hurt your performance a lot, "
diff --git a/include/gmx_omp.h b/include/gmx_omp.h
new file mode 100644 (file)
index 0000000..3fe53b0
--- /dev/null
@@ -0,0 +1,48 @@
+/* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
+ *
+ *
+ *                This source code is part of
+ *
+ *                 G   R   O   M   A   C   S
+ *
+ *          GROningen MAchine for Chemical Simulations
+ *
+ * Written by the Gromacs development team under coordination of
+ * David van der Spoel, Berk Hess, and Erik Lindahl.
+ *
+ * This library 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
+ * of the License, or (at your option) any later version.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org
+ *
+ * And Hey:
+ * GROup of MAchos and Cynical Suckers
+ */
+
+#ifndef GMX_OMP_H
+#define GMX_OMP_H
+
+/* This module defines wrappers for OpenMP API functions and enables compiling
+ * code even when OpenMP is turned off in the build system.
+ * Therfore, OpenMP API functions should always be used through these wrappers
+ * and omp.h should never be directly included. Instead, this header should be
+ * used whnever OpenMP API functions are needed.
+ */
+
+/*! Sets the number of threads in subsequent parallel regions, unless overridden
+ *  by a num_threads clause. Acts as a wrapper for omp_get_max_threads(void). */
+int  gmx_omp_get_max_threads(void);
+
+/*! Returns the thread number of the thread executing within its thread team.
+ *  Acts as a warpper for omp_get_thread_num(void). */
+int  gmx_omp_get_thread_num(void);
+
+/*! Returns an integer that is equal to or greater than the number of threads
+ * that would be available if a parallel region without num_threads were
+ * defined at that point in the code. Acts as a wapepr for omp_set_num_threads(void). */
+void gmx_omp_set_num_threads(int num_threads);
+
+#endif /* GMX_OMP_H */
index 0472b8412d7dd4db6456d8fbbed85508e5d45c31..81fd750a999abbe08747995941e122cd6bd78bd3 100644 (file)
    (MPI or thread_mpi) */
 #cmakedefine GMX_MPI
 
-/* Use threads for parallelization */
+/* Use threads_mpi for parallelization */
 #cmakedefine GMX_THREAD_MPI
 
+/* Use OpenMP multithreading */
+#cmakedefine GMX_OPENMP
+
 /* Use old threading (domain decomp force calc) code */
 #cmakedefine GMX_THREAD_SHM_FDECOMP 
 
diff --git a/src/gmxlib/gmx_omp.c b/src/gmxlib/gmx_omp.c
new file mode 100644 (file)
index 0000000..09b2a19
--- /dev/null
@@ -0,0 +1,61 @@
+/* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
+ *
+ *
+ *                This source code is part of
+ *
+ *                 G   R   O   M   A   C   S
+ *
+ *          GROningen MAchine for Chemical Simulations
+ *
+ * Written by the Gromacs development team under coordination of
+ * David van der Spoel, Berk Hess, and Erik Lindahl.
+ *
+ * This library 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
+ * of the License, or (at your option) any later version.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org
+ *
+ * And Hey:
+ * Gnomes, ROck Monsters And Chili Sauce
+ */
+
+#ifdef HAVE_CONFIG
+#include "config.h"
+#endif
+
+#ifdef GMX_OPENMP
+#include <omp.h>
+#endif
+
+#include "gmx_omp.h"
+
+int gmx_omp_get_max_threads(void)
+{
+#ifdef GMX_OPENMP
+    return omp_get_max_threads();
+#else
+    return 1;
+#endif
+}
+
+
+int gmx_omp_get_thread_num(void)
+{
+#ifdef GMX_OPENMP
+    return omp_get_thread_num();
+#else
+    return 0;
+#endif
+}
+
+void gmx_omp_set_num_threads(int num_threads)
+{
+#ifdef GMX_OPENMP
+    omp_set_num_threads(num_threads);
+#else
+    return;
+#endif
+}
index 1cb9a4bbcad085a6a0597a26a15b4d5486e89ff0..87e2568e1ed5ac7f5f3268dd81a354b5e1239830 100644 (file)
@@ -73,6 +73,8 @@
 #include "membed.h"
 #include "md_openmm.h"
 
+#include "gmx_omp.h"
+
 #ifdef GMX_LIB_MPI
 #include <mpi.h>
 #endif
 #include "md_openmm.h"
 #endif
 
-#ifdef GMX_OPENMP
-#include <omp.h>
-#endif
-
 
 typedef struct { 
     gmx_integrator_t *func;
@@ -829,7 +827,7 @@ int mdrunner(int nthreads_requested, FILE *fplog,t_commrec *cr,int nfile,
             {
                 cpu_set_t mask;
                 CPU_ZERO(&mask);
-                core+=omp_get_thread_num();
+                core+=gmx_omp_get_thread_num();
                 CPU_SET(core,&mask);
                 sched_setaffinity((pid_t) syscall (SYS_gettid),sizeof(cpu_set_t),&mask);
             }
index 1337c1119574ebc42b87e6f09da3e6d3b12e8276..c606b6ad6724eeda32a3205f972666c1b05250d2 100644 (file)
 #endif
 
 #ifdef GMX_OPENMP
+/* TODO: Do we still need this? Are we still planning ot use fftw + OpenMP? */
 #define FFT5D_THREADS
 #endif
 #ifdef FFT5D_THREADS
-#include <omp.h>
+#include "gmx_omp.h"
 /* requires fftw compiled with openmp */
 /* #define FFT5D_FFTW_THREADS (now set by cmake) */
 #endif
index 55e4270397a674e049c155ea72c6005cf124d76a..fac7dbb962688394da9fff3e434122731f1d6a31 100644 (file)
@@ -67,9 +67,7 @@
 #include "tmpi.h"
 #endif
 
-#ifdef GMX_OPENMP
-#include <omp.h>
-#endif
+#include "gmx_omp.h"
 
 #include <stdio.h>
 #include <string.h>
index 0423d7028b00144fce0833eae4adc6f0855845a2..672422932cf3d79dd33ee19d06659ae4deaecf49 100644 (file)
 #include "smalloc.h"
 #include "vec.h"
 #include "geminate.h"
+#include "gmx_omp.h"
 
-#ifdef DOUSEOPENMP
-#define HAVE_OPENMP
-#endif
-#ifdef HAVE_OPENMP
-#include <omp.h>
-#endif
 
 /* The first few sections of this file contain functions that were adopted,
  * and to some extent modified, by Erik Marklund (erikm[aT]xray.bmc.uu.se,
@@ -473,13 +468,11 @@ static double eq10v2(double theoryCt[], double time[], int manytimes,
   part3 = gem_cxmul(gamma, gem_cxmul(gem_cxadd(alpha, beta) , gem_cxsub(alpha, beta)));  /* 3(1+2)(1-2) */
   part4 = gem_cxmul(gem_cxsub(gamma, alpha), gem_cxmul(gem_cxsub(alpha, beta), gem_cxsub(beta, gamma))); /* (3-1)(1-2)(2-3) */
 
-#ifdef HAVE_OPENMP
 #pragma omp parallel for                               \
   private(i, tsqrt, oma, omb, omc, c1, c2, c3, c4),    \
   reduction(+:sumimaginary),                           \
   default(shared),                                     \
   schedule(guided)
-#endif
   for (i=0; i<manytimes; i++){
     tsqrt = sqrt(time[i]);
     oma   = gem_comega(gem_cxrmul(alpha, tsqrt));
@@ -594,13 +587,11 @@ static double gemFunc_residual2(const gsl_vector *p, void *data)
   fixGemACF(GD->ctTheory, nFitPoints);
 
   /* Removing a bunch of points from the log-part. */
-#ifdef HAVE_OPENMP
 #pragma omp parallel for schedule(dynamic)     \
   firstprivate(nData, ctTheory, y, nFitPoints) \
   private (i, iLog, r)                 \
   reduction(+:residual2)                       \
   default(shared)
-#endif
   for(i=0; i<nFitPoints; i++)
     {
       iLog = GD->logtime[i];
@@ -667,8 +658,8 @@ extern real fitGemRecomb(double *ct, double *time, double **ctFit,
 
 #ifdef HAVE_LIBGSL
 #ifdef HAVE_OPENMP
-  nThreads = omp_get_num_procs();
-  omp_set_num_threads(nThreads);
+  nThreads = gmx_omp_get_max_threads();
+  gmx_omp_set_num_threads(nThreads);
   fprintf(stdout, "We will be using %i threads.\n", nThreads);
 #endif
 
index 3193e67588548f1ed55b7e1c0ce8ac66ce492ad5..5981ffe7afdee809182fdf1e02e7e0fac595a408 100644 (file)
@@ -2332,10 +2332,10 @@ static void do_hbac(const char *fn,t_hbdata *hb,
 /* #if (_OPENMP >= 200805) /\* =====================\ *\/ */
 /*         nThreads = min((nThreads <= 0) ? INT_MAX : nThreads, omp_get_thread_limit()); */
 /* #else */
-        nThreads = min((nThreads <= 0) ? INT_MAX : nThreads, omp_get_num_procs());
+        nThreads = min((nThreads <= 0) ? INT_MAX : nThreads, gmx_omp_get_num_procs());
 /* #endif /\* _OPENMP >= 200805 ====================/ *\/ */
 
-        omp_set_num_threads(nThreads);
+        gmx_omp_set_num_threads(nThreads);
         snew(dondata, nThreads);
         for (i=0; i<nThreads; i++)
             dondata[i] = -1;
@@ -2369,7 +2369,7 @@ static void do_hbac(const char *fn,t_hbdata *hb,
     default(shared)
         {
 #pragma omp barrier
-            thisThread = omp_get_thread_num();
+            thisThread = gmx_omp_get_thread_num();
             rhbex = NULL;
 #endif /* ==============================================/ */
 
@@ -2491,7 +2491,7 @@ static void do_hbac(const char *fn,t_hbdata *hb,
         { /* ##########  THE START OF THE ENORMOUS PARALLELIZED BLOCK!  ########## */
             h = NULL;
             g = NULL;
-            thisThread = omp_get_thread_num();
+            thisThread = gmx_omp_get_thread_num();
             snew(h,hb->maxhydro);
             snew(g,hb->maxhydro);
             mMax = INT_MIN;
@@ -3521,9 +3521,9 @@ int gmx_hbond(int argc,char *argv[])
 /* #if (_OPENMP > 200805) */
 /*         actual_nThreads = min((nThreads <= 0) ? INT_MAX : nThreads, omp_get_thread_limit()); */
 /* #else */
-        actual_nThreads = min((nThreads <= 0) ? INT_MAX : nThreads, omp_get_num_procs());
+        actual_nThreads = min((nThreads <= 0) ? INT_MAX : nThreads, gmx_omp_get_num_procs());
 /* #endif */
-        omp_set_num_threads(actual_nThreads);
+        gmx_omp_set_num_threads(actual_nThreads);
         printf("Frame loop parallelized with OpenMP using %i threads.\n", actual_nThreads);
         fflush(stdout);
     }
@@ -3591,7 +3591,7 @@ int gmx_hbond(int argc,char *argv[])
            bGem, oenv, fnm, fpnhb, trrStatus, natoms,   \
            status, nabin, nrbin, adist, rdist, debug)
     {    /* Start of parallel region */
-        threadNr = omp_get_thread_num();
+        threadNr = gmx_omp_get_thread_num();
 #endif /* HAVE_OPENMP ================================================= */
         do
         {
index 2cae68baa5c3d37628331e7b544e81541d72e812..fcf99c8e752bd9292b7f5080ce4b787f704894a6 100644 (file)
 #include "matio.h"
 #include "gmx_ana.h"
 #include "nsfactor.h"
-
-#ifdef GMX_OPENMP
-#include <omp.h>
-#endif
-
+#include "gmx_omp.h"
 
 int gmx_sans(int argc,char *argv[])
 {
@@ -153,9 +149,7 @@ int gmx_sans(int argc,char *argv[])
       { efXVG, "-pr",         "pr",   ffWRITE }
   };
 
-#ifdef GMX_OPENMP
-    nthreads = omp_get_max_threads();
-#endif
+  nthreads = gmx_omp_get_max_threads();
 
   CopyRight(stderr,argv[0]);
   parse_common_args(&argc,argv,PCA_BE_NICE,
@@ -166,9 +160,8 @@ int gmx_sans(int argc,char *argv[])
   check_mcover(mcover);
 
   /* setting number of omp threads globaly */
-#ifdef GMX_OPENMP
-  omp_set_num_threads(nthreads);
-#endif
+  gmx_omp_set_num_threads(nthreads);
+
   /* Now try to parse opts for modes */
   switch(emethod[0][0]) {
   case 'd':
index e7879e7ca0bf613b5e3b4d7a4fb90c75f9b0e235..863b658e6e3e687685a3e19b5c19ed98d7c7b16a 100644 (file)
 #include "strdb.h"
 #include "vec.h"
 #include "nsfactor.h"
-
-#ifdef GMX_OPENMP
-#include <omp.h>
-#endif
+#include "gmx_omp.h"
 
 void check_binwidth(real binwidth) {
     real smallest_bin=0.1;
@@ -202,7 +199,7 @@ gmx_radial_distribution_histogram_t *calc_radial_distribution_histogram (
         }
         rng=gmx_rng_init(seed);
 #ifdef GMX_OPENMP
-        nthreads = omp_get_max_threads();
+        nthreads = gmx_omp_get_max_threads();
         snew(tgr,nthreads);
         snew(trng,nthreads);
         for(i=0;i<nthreads;i++){
@@ -211,7 +208,7 @@ gmx_radial_distribution_histogram_t *calc_radial_distribution_histogram (
         }
         #pragma omp parallel shared(tgr,trng,mc) private(tid,i,j)
         {
-            tid = omp_get_thread_num();
+            tid = gmx_omp_get_thread_num();
             /* now starting parallel threads */
             #pragma omp for
             for(mc=0;mc<max;mc++) {
@@ -246,7 +243,7 @@ gmx_radial_distribution_histogram_t *calc_radial_distribution_histogram (
         gmx_rng_destroy(rng);
     } else {
 #ifdef GMX_OPENMP
-        nthreads = omp_get_max_threads();
+        nthreads = gmx_omp_get_max_threads();
         /* Allocating memory for tgr arrays */
         snew(tgr,nthreads);
         for(i=0;i<nthreads;i++) {
@@ -254,7 +251,7 @@ gmx_radial_distribution_histogram_t *calc_radial_distribution_histogram (
         }
         #pragma omp parallel shared(tgr) private(tid,i,j)
         {
-            tid = omp_get_thread_num();
+            tid = gmx_omp_get_thread_num();
             /* starting parallel threads */
             #pragma omp for
             for(i=0;i<isize;i++) {