Merge branch release-5-1
[alexxy/gromacs.git] / src / gromacs / correlationfunctions / crosscorr.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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
36  * \file
37  * \brief
38  * Implements routine for computing a cross correlation between two data sets
39  *
40  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41  * \ingroup module_correlationfunctions
42  */
43 #include "gmxpre.h"
44
45 #include "crosscorr.h"
46
47 #include "gromacs/fft/fft.h"
48 #include "gromacs/legacyheaders/macros.h"
49 #include "gromacs/utility/smalloc.h"
50
51 /*! \brief
52  * return the size witch the array should be after zero padding
53  *
54  * \param[in] n Factor to multiply by 2
55  * \return zeroPaddingSize
56  */
57 static int zeroPaddingSize(int n)
58 {
59     return 2*n;
60 }
61
62 /*! \brief
63  * Compute complex conjugate. Output in the first input variable.
64  *
65  * \param[in] in1 first complex number
66  * \param[in] in2 second complex number
67  */
68 static void complexConjugatMult(double in1[], double in2[])
69 {
70     double res[2];
71     res[0]  = in1[0] * in2[0] + in1[1] * in2[1];
72     res[1]  = in1[0] * -in2[1] + in1[1] * in2[0];
73     in1[0]  = res[0];
74     in1[1]  = res[1];
75 }
76
77 /*! \brief
78  * Compute one cross correlation corr = f x g using FFT.
79  *
80  * \param[in] n number of data point
81  * \param[in] f first function
82  * \param[in] g second function
83  * \param[out] corr output correlation
84  * \param[in] fft FFT data structure
85  */
86 static void cross_corr_low(int n, real f[], real g[], real corr[], gmx_fft_t fft)
87 {
88     int             i;
89     const int       size = zeroPaddingSize(n);
90     double    **    in1, ** in2;
91     snew(in1, size);
92     snew(in2, size);
93
94     for (i = 0; i < size; i++)
95     {
96         snew(in1[i], 2);
97         snew(in2[i], 2);
98     }
99
100     for (i = 0; i < n; i++)
101     {
102         in1[i][0]  = (double)f[i];
103         in1[i][1]  = 0;
104         in2[i][0]  = (double)g[i];
105         in2[i][1]  = 0;
106     }
107     for (; i < size; i++)
108     {
109         in1[i][0]  = 0;
110         in1[i][1]  = 0;
111         in2[i][0]  = 0;
112         in2[i][1]  = 0;
113     }
114
115
116     gmx_fft_1d(fft, GMX_FFT_FORWARD, in1, in1);
117     gmx_fft_1d(fft, GMX_FFT_FORWARD, in2, in2);
118
119     for (i = 0; i < size; i++)
120     {
121         complexConjugatMult(in1[i], in2[i]);
122         in1[i][0] /= size;
123     }
124     gmx_fft_1d(fft, GMX_FFT_BACKWARD, in1, in1);
125
126     for (i = 0; i < n; i++)
127     {
128         corr[i] = (real)(in1[i][0]);
129     }
130
131     for (i = 0; i < size; i++)
132     {
133         sfree(in1[i]);
134         sfree(in2[i]);
135     }
136
137     sfree(in1);
138     sfree(in2);
139
140 }
141
142 void cross_corr(int n, real f[], real g[], real corr[])
143 {
144     gmx_fft_t fft;
145     gmx_fft_init_1d(&fft, zeroPaddingSize(n), GMX_FFT_FLAG_CONSERVATIVE);
146     cross_corr_low( n,  f,  g, corr, fft);
147     gmx_fft_destroy(fft);
148     gmx_fft_cleanup();
149 }
150
151 void many_cross_corr(int nFunc, int * nData, real ** f, real ** g, real ** corr)
152 {
153 #pragma omp parallel
154     //gmx_fft_t is not thread safe, so structure are allocated per thread.
155     {
156         int i;
157
158 #pragma omp for
159         for (i = 0; i < nFunc; i++)
160         {
161             gmx_fft_t fft;
162             gmx_fft_init_1d(&fft, zeroPaddingSize(nData[i]), GMX_FFT_FLAG_CONSERVATIVE);
163             cross_corr_low( nData[i],  f[i],  g[i], corr[i], fft);
164             gmx_fft_destroy(fft);
165         }
166     }
167     gmx_fft_cleanup();
168
169 }