Added CMake support for Extrae.
authorRossen Apostolov <rossen@kth.se>
Wed, 5 Mar 2014 16:12:05 +0000 (17:12 +0100)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Mon, 31 Mar 2014 17:33:47 +0000 (19:33 +0200)
Website of the project is http://www.bsc.es/computer-sciences/extrae

Change-Id: Ifc71637ae7fbe855f3879b19f0968e929a92daf0

CMakeLists.txt
cmake/FindEXTRAE.cmake [new file with mode: 0644]
src/config.h.cmakein
src/gromacs/CMakeLists.txt
src/gromacs/gmxlib/copyrite.cpp

index 63e61763c3dc068c497fc8070f4a779dc5aa7681..ff512eeab0c0e07cfe6abccfccce5ec26617fe73 100644 (file)
@@ -440,6 +440,19 @@ if (GMX_GSL)
   endif()
 endif()
 
+option(GMX_EXTRAE "Add support for tracing using EXTRAE" OFF)
+mark_as_advanced(GMX_EXTRAE)
+
+if (GMX_EXTRAE)
+  find_package(EXTRAE)
+  if(EXTRAE_FOUND)
+    include_directories(${EXTRAE_INCLUDE_DIR})
+    set(HAVE_EXTRAE 1)
+  else()
+    message(FATAL_ERROR "EXTRAE library was not found. Please add the correct path to CMAKE_PREFIX_PATH")
+  endif()
+endif()
+
 option(GMX_X11 "Use X window system" OFF)
 if (GMX_X11)
     find_package(X11)
diff --git a/cmake/FindEXTRAE.cmake b/cmake/FindEXTRAE.cmake
new file mode 100644 (file)
index 0000000..02042f7
--- /dev/null
@@ -0,0 +1,106 @@
+#
+# 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.
+
+# - Try to find EXTRAE
+# Once done this will define
+#  EXTRAE_FOUND - System has EXTRAE
+#  EXTRAE_INCLUDE_DIRS - The EXTRAE include directories
+#  EXTRAE_LIBRARIES - The libraries needed to use EXTRAE
+
+find_path(EXTRAE_INCLUDE_DIR extrae_user_events.h)
+
+# EXTRAE libraries have different names depending on the supported features,
+# and in extrae-2.5.0 the following combinations are supported:
+
+# seqtrace:        seq code
+# mpitrace:        MPI
+# omptrace:        OpenMP
+# ompitrace:       MPI + OpenMP
+# pttrace:         pthreads
+# ptmpitrace       pthreads + MPI (unsupported combination in Gromacs)
+# cudatrace:       CUDA
+# cudampitrace:    CUDA + MPI
+# cudaompitrace:   CUDA+OPENMP+MPI (in the dev version)
+
+# TODO: Add support for the following combinations when available in a future release:
+
+# cudaomptrace:    CUDA + OPENMP
+# cudapttrace:        CUDA + pthreads
+# cudaptmpitrace:    CUDA + pthreads + MPI (unsupported combination in Gromacs)
+
+set (extraelib "trace")
+
+# libs with MPI support
+if (GMX_MPI)
+  if (GMX_OPENMP)
+    set (extraelib "ompi${extraelib}")
+  else()
+    set (extraelib "mpi${extraelib}")
+  endif()
+  if (GMX_GPU)
+    set (extraelib "cuda${extraelib}")
+  endif()
+
+# other libs with OpenMP support
+elseif (GMX_OPENMP)
+  set (extraelib "omp${extraelib}")
+    if (GMX_GPU)
+      set (extraelib "cuda${extraelib}")
+    endif()
+
+# library with CUDA only support
+elseif (GMX_GPU)
+    set (extraelib "cuda${extraelib}")
+
+# library with PThreads support
+elseif (GMX_THREAD_MPI)
+    set (extraelib "pt${extraelib}")
+
+else()
+  set (extraelib "seq${extraelib}")
+endif()
+
+find_library(EXTRAE_LIBRARY NAMES  ${extraelib})
+
+set(EXTRAE_LIBRARIES ${EXTRAE_LIBRARY} )
+set(EXTRAE_INCLUDE_DIRS ${EXTRAE_INCLUDE_DIR} )
+
+# handle the QUIETLY and REQUIRED arguments and set EXTRAE_FOUND to TRUE
+# if all listed variables are TRUE
+include(FindPackageHandleStandardArgs)
+find_package_handle_standard_args(EXTRAE  DEFAULT_MSG
+                                  EXTRAE_LIBRARY EXTRAE_INCLUDE_DIR)
+
+mark_as_advanced(EXTRAE_INCLUDE_DIR EXTRAE_LIBRARY )
+
index 98c19e39468cdb115ee11a65afc05ef902b53e76..0a75d07ac35c9137dfe4ba4588cb9fbe672d776c 100644 (file)
 /* Compile to use TNG library */
 #cmakedefine GMX_USE_TNG
 
+/* Add support for tracing using Extrae */
+#cmakedefine HAVE_EXTRAE
+
 /* Use MPI (with mpicc) for parallelization */
 #cmakedefine GMX_LIB_MPI
 
index ec3b87ec71c577c8defce380cc9e954e58b700ce..21ba8166a63e2a0b9e5fc3159bd896fbb75e55ea 100644 (file)
@@ -156,6 +156,7 @@ set_source_files_properties(selection/scanner.cpp PROPERTIES COMPILE_FLAGS "${_s
 target_link_libraries(libgromacs ${GMX_GPU_LIBRARIES}
                       ${GMX_EXTRA_LIBRARIES}
                       ${GMX_TNG_LIBRARIES}
+                      ${EXTRAE_LIBRARIES}
                       ${FFT_LIBRARIES} ${LINEAR_ALGEBRA_LIBRARIES}
                       ${XML_LIBRARIES} ${GSL_LIBRARIES}
                       ${THREAD_LIB} ${GMX_SHARED_LINKER_FLAGS})
index d3d67f9b8284de188b01e51a309c458793b7931f..e3f3033ca1ecafba3ef08e2565730ec645dfd37e 100644 (file)
 #include <mkl.h>
 #endif
 
+#ifdef HAVE_EXTRAE
+ #include "extrae_user_events.h"
+#endif
+
 #include <boost/version.hpp>
 
 /* This file is completely threadsafe - keep it that way! */
@@ -693,6 +697,14 @@ static void gmx_print_version_info(FILE *fp)
 #else
     fprintf(fp, "TNG support:        disabled\n");
 #endif
+#ifdef HAVE_EXTRAE
+    unsigned major, minor, revision;
+    Extrae_get_version(&major, &minor, &revision);
+    fprintf(fp, "Tracing support:    enabled. Using Extrae-%d.%d.%d\n", major, minor, revision);
+#else
+    fprintf(fp, "Tracing support:    disabled\n");
+#endif
+
 
     fprintf(fp, "Built on:           %s\n", BUILD_TIME);
     fprintf(fp, "Built by:           %s\n", BUILD_USER);