Converted correlationfuntcions to C++.
authorRossen Apostolov <rossen@kth.se>
Wed, 3 Jun 2015 10:13:01 +0000 (12:13 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Fri, 24 Jul 2015 14:34:40 +0000 (16:34 +0200)
Change-Id: If5a478eb05fc295670cfa982525838c5e1e255b7

src/gromacs/correlationfunctions/CMakeLists.txt
src/gromacs/correlationfunctions/autocorr.cpp [moved from src/gromacs/correlationfunctions/autocorr.c with 97% similarity]
src/gromacs/correlationfunctions/crosscorr.cpp [moved from src/gromacs/correlationfunctions/crosscorr.c with 98% similarity]
src/gromacs/correlationfunctions/integrate.cpp [moved from src/gromacs/correlationfunctions/integrate.c with 98% similarity]
src/gromacs/correlationfunctions/manyautocorrelation.cpp [moved from src/gromacs/correlationfunctions/manyautocorrelation.c with 89% similarity]
src/gromacs/correlationfunctions/polynomials.cpp [moved from src/gromacs/correlationfunctions/polynomials.c with 97% similarity]

index 6c6718af5d6b8c522d1604bc3c1624408645bda7..88d3bd88b34f1705d6599700052bdebb7f3ac058 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2014, by the GROMACS development team, led by
+# Copyright (c) 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.
@@ -32,7 +32,7 @@
 # To help us fund GROMACS development, we humbly ask that you cite
 # the research papers on the package. Check out http://www.gromacs.org.
 
-file(GLOB GMXCORRFUNC_SOURCES *.c *.cpp)
+file(GLOB GMXCORRFUNC_SOURCES *.cpp)
 file(GLOB LMFIT_SOURCES ${CMAKE_SOURCE_DIR}/src/external/lmfit/*.c)
 
 set(LIBGROMACS_SOURCES ${LIBGROMACS_SOURCES} ${GMXCORRFUNC_SOURCES}  ${LMFIT_SOURCES} PARENT_SCOPE)
similarity index 97%
rename from src/gromacs/correlationfunctions/autocorr.c
rename to src/gromacs/correlationfunctions/autocorr.cpp
index 72966978ce38d0a90f5fef87c4e35d1f6ffbcf61..7cdf48d57796f5b184fa53249d04a480e270e6ea 100644 (file)
 
 #include "autocorr.h"
 
-#include <math.h>
 #include <stdio.h>
 #include <string.h>
 
+#include <cmath>
+
+#include <algorithm>
+
 #include "gromacs/correlationfunctions/expfit.h"
 #include "gromacs/correlationfunctions/integrate.h"
 #include "gromacs/correlationfunctions/manyautocorrelation.h"
@@ -87,16 +90,12 @@ static void low_do_four_core(int nfour, int nframes, real c1[], real cfour[],
                              int nCos)
 {
     int  i = 0;
-    int  fftcode;
-    real aver, *ans;
 
-    aver = 0.0;
     switch (nCos)
     {
         case enNorm:
             for (i = 0; (i < nframes); i++)
             {
-                aver    += c1[i];
                 cfour[i] = c1[i];
             }
             break;
@@ -116,7 +115,7 @@ static void low_do_four_core(int nfour, int nframes, real c1[], real cfour[],
             gmx_fatal(FARGS, "nCos = %d, %s %d", nCos, __FILE__, __LINE__);
     }
 
-    fftcode = many_auto_correl(1, nframes, nfour, &cfour);
+    many_auto_correl(1, nframes, nfour, &cfour);
 }
 
 /*! \brief Routine to comput ACF without FFT. */
@@ -577,8 +576,8 @@ void low_do_autocorr(const char *fn, const output_env_t oenv, const char *title,
     }
     if (bFour)
     {
-        c0 = log((double)nframes)/log(2.0);
-        k  = c0;
+        c0 = std::log(static_cast<double>(nframes))/std::log(2.0);
+        k  = static_cast<int>(c0);
         if (k < c0)
         {
             k++;
@@ -606,10 +605,9 @@ void low_do_autocorr(const char *fn, const output_env_t oenv, const char *title,
      * In this loop the actual correlation functions are computed, but without
      * normalizing them.
      */
-    k = max(1, pow(10, (int)(log(nitem)/log(100))));
-    for (i = 0; i < nitem; i++)
+    for (int i = 0; i < nitem; i++)
     {
-        if (bVerbose && ((i%k == 0 || i == nitem-1)))
+        if (bVerbose && (((i % 100) == 0) || (i == nitem-1)))
         {
             fprintf(stderr, "\rThingie %d", i+1);
         }
@@ -660,10 +658,10 @@ void low_do_autocorr(const char *fn, const output_env_t oenv, const char *title,
         else
         {
             sum = print_and_integrate(fp, nout, dt, c1[0], NULL, 1);
-            if (bVerbose)
-            {
-                printf("Correlation time (integral over corrfn): %g (ps)\n", sum);
-            }
+        }
+        if (bVerbose)
+        {
+            printf("Correlation time (integral over corrfn): %g (ps)\n", sum);
         }
     }
     else
similarity index 98%
rename from src/gromacs/correlationfunctions/crosscorr.c
rename to src/gromacs/correlationfunctions/crosscorr.cpp
index 1911b5b41e176ab0e7e6f67db618efc1e6227c9a..b588bd8f86ad5cec0dd657293f236a1bf8ed0f2f 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2014, by the GROMACS development team, led by
+ * Copyright (c) 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.
similarity index 98%
rename from src/gromacs/correlationfunctions/integrate.c
rename to src/gromacs/correlationfunctions/integrate.cpp
index dc0a7ae0834ee8f834576fe8459b2e447c19b076..d58c5add9d8ba5d2f6d322242450493fedce526d 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.
similarity index 89%
rename from src/gromacs/correlationfunctions/manyautocorrelation.c
rename to src/gromacs/correlationfunctions/manyautocorrelation.cpp
index c4fb7db465c5716f7d5a32ac7e21759cea165f85..8474c867a35ace234da337d420bee75d393805c8 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2014, by the GROMACS development team, led by
+ * Copyright (c) 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.
 
 #include "manyautocorrelation.h"
 
-#include <math.h>
 #include <stdio.h>
 #include <stdlib.h>
 
+#include <cmath>
+
+#include <algorithm>
+
 #include "gromacs/fft/fft.h"
 #include "gromacs/legacyheaders/macros.h"
 #include "gromacs/utility/gmxomp.h"
@@ -57,7 +60,7 @@ int many_auto_correl(int nfunc, int ndata, int nfft, real **c)
     #pragma omp parallel
     {
         typedef real complex[2];
-        int          i, t, j, fftcode;
+        int          i, j;
         gmx_fft_t    fft1;
         complex     *in, *out;
         int          i0, i1;
@@ -70,9 +73,9 @@ int many_auto_correl(int nfunc, int ndata, int nfft, real **c)
             // fprintf(stderr, "There are %d threads for correlation functions\n", nthreads);
         }
         i0 = thread_id*nfunc/nthreads;
-        i1 = min(nfunc, (thread_id+1)*nfunc/nthreads);
+        i1 = std::min(nfunc, (thread_id+1)*nfunc/nthreads);
 
-        fftcode = gmx_fft_init_1d(&fft1, nfft, GMX_FFT_FLAG_CONSERVATIVE);
+        gmx_fft_init_1d(&fft1, nfft, GMX_FFT_FLAG_CONSERVATIVE);
         /* Allocate temporary arrays */
         snew(in, nfft);
         snew(out, nfft);
@@ -88,7 +91,7 @@ int many_auto_correl(int nfunc, int ndata, int nfft, real **c)
                 in[j][0] = in[j][1] = 0;
             }
 
-            fftcode = gmx_fft_1d(fft1, GMX_FFT_BACKWARD, (void *)in, (void *)out);
+            gmx_fft_1d(fft1, GMX_FFT_BACKWARD, (void *)in, (void *)out);
             for (j = 0; j < nfft; j++)
             {
                 in[j][0] = (out[j][0]*out[j][0] + out[j][1]*out[j][1])/nfft;
@@ -99,7 +102,7 @@ int many_auto_correl(int nfunc, int ndata, int nfft, real **c)
                 in[j][0] = in[j][1] = 0;
             }
 
-            fftcode = gmx_fft_1d(fft1, GMX_FFT_FORWARD, (void *)in, (void *)out);
+            gmx_fft_1d(fft1, GMX_FFT_FORWARD, (void *)in, (void *)out);
             for (j = 0; (j < nfft); j++)
             {
                 c[i][j] = out[j][0]/ndata;
similarity index 97%
rename from src/gromacs/correlationfunctions/polynomials.c
rename to src/gromacs/correlationfunctions/polynomials.cpp
index c5f499a17daf50f92fbd1b2a13d3db484427c02c..18765ab1ff7f8c37c3daee3c9a685cf892314b09 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.