c23d4c694e45069b8fea0de97f160108a75b082c
[alexxy/gromacs.git] / src / gromacs / correlationfunctions / manyautocorrelation.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2018,2019, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements function to compute many autocorrelation functions
38  *
39  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40  * \ingroup module_correlationfunctions
41  */
42 #include "gmxpre.h"
43
44 #include "manyautocorrelation.h"
45
46 #include <algorithm>
47
48 #include "gromacs/fft/fft.h"
49 #include "gromacs/utility/exceptions.h"
50 #include "gromacs/utility/gmxomp.h"
51
52 int many_auto_correl(std::vector<std::vector<real>>* c)
53 {
54     size_t nfunc = (*c).size();
55     if (nfunc == 0)
56     {
57         GMX_THROW(gmx::InconsistentInputError("Empty array of vectors supplied"));
58     }
59     size_t ndata = (*c)[0].size();
60     if (ndata == 0)
61     {
62         GMX_THROW(gmx::InconsistentInputError("Empty vector supplied"));
63     }
64 #ifndef NDEBUG
65     for (size_t i = 1; i < nfunc; i++)
66     {
67         if ((*c)[i].size() != ndata)
68         {
69             char buf[256];
70             snprintf(buf, sizeof(buf), "Vectors of different lengths supplied (%d %d)",
71                      static_cast<int>((*c)[i].size()), static_cast<int>(ndata));
72             GMX_THROW(gmx::InconsistentInputError(buf));
73         }
74     }
75 #endif
76     // Add buffer size to the arrays.
77     size_t nfft = (3 * ndata / 2) + 1;
78     // Pad arrays with zeros
79     for (auto& i : *c)
80     {
81         i.resize(nfft, 0);
82     }
83 #pragma omp parallel
84     {
85         try
86         {
87             gmx_fft_t         fft1;
88             std::vector<real> in, out;
89
90             int nthreads  = gmx_omp_get_max_threads();
91             int thread_id = gmx_omp_get_thread_num();
92             int i0        = (thread_id * nfunc) / nthreads;
93             int i1        = std::min(nfunc, ((thread_id + 1) * nfunc) / nthreads);
94
95             gmx_fft_init_1d(&fft1, nfft, GMX_FFT_FLAG_CONSERVATIVE);
96             /* Allocate temporary arrays */
97             in.resize(2 * nfft, 0);
98             out.resize(2 * nfft, 0);
99             for (int i = i0; (i < i1); i++)
100             {
101                 for (size_t j = 0; j < ndata; j++)
102                 {
103                     in[2 * j + 0] = (*c)[i][j];
104                     in[2 * j + 1] = 0;
105                 }
106                 gmx_fft_1d(fft1, GMX_FFT_BACKWARD, in.data(), out.data());
107                 for (size_t j = 0; j < nfft; j++)
108                 {
109                     in[2 * j + 0] =
110                             (out[2 * j + 0] * out[2 * j + 0] + out[2 * j + 1] * out[2 * j + 1]) / nfft;
111                     in[2 * j + 1] = 0;
112                 }
113                 gmx_fft_1d(fft1, GMX_FFT_FORWARD, in.data(), out.data());
114                 for (size_t j = 0; (j < nfft); j++)
115                 {
116                     (*c)[i][j] = out[2 * j + 0];
117                 }
118             }
119             /* Free the memory */
120             gmx_fft_destroy(fft1);
121         }
122         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
123     }
124     for (auto& i : *c)
125     {
126         i.resize(ndata);
127     }
128
129     return 0;
130 }