2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * Implements function to compute many autocorrelation functions
39 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40 * \ingroup module_correlationfunctions
44 #include "manyautocorrelation.h"
53 #include "gromacs/fft/fft.h"
54 #include "gromacs/legacyheaders/macros.h"
55 #include "gromacs/utility/gmxomp.h"
56 #include "gromacs/utility/smalloc.h"
58 int many_auto_correl(int nfunc, int ndata, int nfft, real **c)
62 typedef real complex[2];
67 int nthreads, thread_id;
69 nthreads = gmx_omp_get_max_threads();
70 thread_id = gmx_omp_get_thread_num();
73 // fprintf(stderr, "There are %d threads for correlation functions\n", nthreads);
75 i0 = thread_id*nfunc/nthreads;
76 i1 = std::min(nfunc, (thread_id+1)*nfunc/nthreads);
78 gmx_fft_init_1d(&fft1, nfft, GMX_FFT_FLAG_CONSERVATIVE);
79 /* Allocate temporary arrays */
82 for (i = i0; (i < i1); i++)
84 for (j = 0; j < ndata; j++)
89 for (; (j < nfft); j++)
91 in[j][0] = in[j][1] = 0;
94 gmx_fft_1d(fft1, GMX_FFT_BACKWARD, (void *)in, (void *)out);
95 for (j = 0; j < nfft; j++)
97 in[j][0] = (out[j][0]*out[j][0] + out[j][1]*out[j][1])/nfft;
100 for (; (j < nfft); j++)
102 in[j][0] = in[j][1] = 0;
105 gmx_fft_1d(fft1, GMX_FFT_FORWARD, (void *)in, (void *)out);
106 for (j = 0; (j < nfft); j++)
108 c[i][j] = out[j][0]/ndata;
111 /* Free the memory */
112 gmx_fft_destroy(fft1);
116 // gmx_fft_cleanup();