Move symtab.* to topology/
[gromacs.git] / src / gromacs / gmxana / correl.c
blobb4a2c9b47516863d2f2ec01e81730e26b9335c99
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2008, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
42 #include <stdio.h>
43 #include <stdlib.h>
44 #include <math.h>
46 #include "gromacs/fft/fft.h"
47 #include "gromacs/utility/smalloc.h"
48 #include "correl.h"
50 #define SWAP(a, b) tempr = (a); (a) = (b); (b) = tempr
52 void four1(real data[], int nn, int isign)
54 int n, mmax, m, j, istep, i;
55 double wtemp, wr, wpr, wpi, wi, theta;
56 real tempr, tempi;
58 n = nn << 1;
59 j = 1;
60 for (i = 1; i < n; i += 2)
62 if (j > i)
64 SWAP(data[j], data[i]);
65 SWAP(data[j+1], data[i+1]);
67 m = n >> 1;
68 while (m >= 2 && j > m)
70 j -= m;
71 m >>= 1;
73 j += m;
75 mmax = 2;
76 while (n > mmax)
78 istep = 2*mmax;
79 theta = 6.28318530717959/(isign*mmax);
80 wtemp = sin(0.5*theta);
81 wpr = -2.0*wtemp*wtemp;
82 wpi = sin(theta);
83 wr = 1.0;
84 wi = 0.0;
85 for (m = 1; m < mmax; m += 2)
87 for (i = m; i <= n; i += istep)
89 j = i+mmax;
90 tempr = wr*data[j]-wi*data[j+1];
91 tempi = wr*data[j+1]+wi*data[j];
92 data[j] = data[i]-tempr;
93 data[j+1] = data[i+1]-tempi;
94 data[i] += tempr;
95 data[i+1] += tempi;
97 wr = (wtemp = wr)*wpr-wi*wpi+wr;
98 wi = wi*wpr+wtemp*wpi+wi;
100 mmax = istep;
104 #undef SWAP
106 static void realft(real data[], int n, int isign)
108 int i, i1, i2, i3, i4, n2p3;
109 real c1 = 0.5, c2, h1r, h1i, h2r, h2i;
110 double wr, wi, wpr, wpi, wtemp, theta;
112 theta = 3.141592653589793/(double) n;
113 if (isign == 1)
115 c2 = -0.5;
116 four1(data, n, 1);
118 else
120 c2 = 0.5;
121 theta = -theta;
123 wtemp = sin(0.5*theta);
124 wpr = -2.0*wtemp*wtemp;
125 wpi = sin(theta);
126 wr = 1.0+wpr;
127 wi = wpi;
128 n2p3 = 2*n+3;
129 for (i = 2; i <= n/2; i++)
131 i4 = 1+(i3 = n2p3-(i2 = 1+(i1 = i+i-1)));
132 h1r = c1*(data[i1]+data[i3]);
133 h1i = c1*(data[i2]-data[i4]);
134 h2r = -c2*(data[i2]+data[i4]);
135 h2i = c2*(data[i1]-data[i3]);
136 data[i1] = h1r+wr*h2r-wi*h2i;
137 data[i2] = h1i+wr*h2i+wi*h2r;
138 data[i3] = h1r-wr*h2r+wi*h2i;
139 data[i4] = -h1i+wr*h2i+wi*h2r;
140 wr = (wtemp = wr)*wpr-wi*wpi+wr;
141 wi = wi*wpr+wtemp*wpi+wi;
143 if (isign == 1)
145 data[1] = (h1r = data[1])+data[2];
146 data[2] = h1r-data[2];
148 else
150 data[1] = c1*((h1r = data[1])+data[2]);
151 data[2] = c1*(h1r-data[2]);
152 four1(data, n, -1);
156 static void twofft(real data1[], real data2[], real fft1[], real fft2[], int n)
158 int nn3, nn2, jj, j;
159 real rep, rem, aip, aim;
161 nn3 = 1+(nn2 = 2+n+n);
162 for (j = 1, jj = 2; j <= n; j++, jj += 2)
164 fft1[jj-1] = data1[j];
165 fft1[jj] = data2[j];
167 four1(fft1, n, 1);
168 fft2[1] = fft1[2];
169 fft1[2] = fft2[2] = 0.0;
170 for (j = 3; j <= n+1; j += 2)
172 rep = 0.5*(fft1[j]+fft1[nn2-j]);
173 rem = 0.5*(fft1[j]-fft1[nn2-j]);
174 aip = 0.5*(fft1[j+1]+fft1[nn3-j]);
175 aim = 0.5*(fft1[j+1]-fft1[nn3-j]);
176 fft1[j] = rep;
177 fft1[j+1] = aim;
178 fft1[nn2-j] = rep;
179 fft1[nn3-j] = -aim;
180 fft2[j] = aip;
181 fft2[j+1] = -rem;
182 fft2[nn2-j] = aip;
183 fft2[nn3-j] = rem;
187 void correl(real data1[], real data2[], int n, real ans[])
189 int no2, i;
190 real dum, *fft;
192 snew(fft, 2*n+1);
193 twofft(data1, data2, fft, ans, n);
194 no2 = n/2;
195 for (i = 2; i <= n+2; i += 2)
197 dum = ans[i-1];
198 ans[i-1] = (fft[i-1]*dum+fft[i]*ans[i])/no2;
199 ans[i] = (fft[i]*dum-fft[i-1]*ans[i])/no2;
201 ans[2] = ans[n+1];
202 realft(ans, no2, -1);
203 sfree(fft);