Merge release-4-6 into release-5-0
[alexxy/gromacs.git] / src / gromacs / gmxana / powerspect.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2013,2014, 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
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include "gromacs/utility/smalloc.h"
41 #include "gromacs/fft/fft.h"
42 #include "gmx_fatal.h"
43 #include "gromacs/fileio/futil.h"
44 #include "interf.h"
45 #include "powerspect.h"
46
47 void addtoavgenergy(t_complex *list, real *result, int size, int tsteps)
48 {
49     int i;
50     for (i = 0; i < size; i++)
51     {
52         result[i] += cabs2(list[i])/tsteps;
53     }
54
55 }
56
57
58 void powerspectavg(real ***intftab, int tsteps, int xbins, int ybins, char **outfiles)
59 {
60     /*Fourier plans and output;*/
61     gmx_fft_t  fftp;
62     t_complex *ftspect1;            /* Spatial FFT of interface for each time frame and interface ftint[time,xycoord][0], ftintf[time,xycoord][1] for interface 1 and 2                 respectively */
63     t_complex *ftspect2;
64     real      *pspectavg1;          /*power -spectrum 1st interface*/
65     real      *pspectavg2;          /* -------------- 2nd interface*/
66     real      *temp;
67     FILE      *datfile1, *datfile2; /*data-files with interface data*/
68     int        n;                   /*time index*/
69     int        fy  = ybins/2+1;     /* number of (symmetric) fourier y elements; */
70     int        rfl = xbins*fy;      /*length of real - DFT == Symmetric 2D matrix*/
71     int        status;
72
73 /*Prepare data structures for FFT, with time averaging of power spectrum*/
74     if ( (status = gmx_fft_init_2d_real(&fftp, xbins, ybins, GMX_FFT_FLAG_NONE) ) != 0)
75     {
76         free(fftp);
77         gmx_fatal(status, __FILE__, __LINE__, "Error allocating FFT");
78     }
79
80 /*Initialize arrays*/
81     snew(ftspect1, rfl);
82     snew(ftspect2, rfl);
83     snew(temp, rfl);
84     snew(pspectavg1, rfl);
85     snew(pspectavg2, rfl);
86
87 /*Fouriertransform directly (no normalization or anything)*/
88 /*NB! Check carefully indexes here*/
89
90     for (n = 0; n < tsteps; n++)
91     {
92         gmx_fft_2d_real(fftp, GMX_FFT_REAL_TO_COMPLEX, intftab[0][n], ftspect1);
93         gmx_fft_2d_real(fftp, GMX_FFT_REAL_TO_COMPLEX, intftab[1][n], ftspect2);
94
95         /*Add to average for interface 1 here*/
96         addtoavgenergy(ftspect1, pspectavg1, rfl, tsteps);
97         /*Add to average for interface 2 here*/
98         addtoavgenergy(ftspect2, pspectavg2, rfl, tsteps);
99     }
100     /*Print out average energy-spectrum to outfiles[0] and outfiles[1];*/
101
102     datfile1 = gmx_ffopen(outfiles[0], "w");
103     datfile2 = gmx_ffopen(outfiles[1], "w");
104
105 /*Filling dat files with spectral data*/
106     fprintf(datfile1, "%s\n", "kx\t ky\t\tPower(kx,ky)");
107     fprintf(datfile2, "%s\n", "kx\t ky\t\tPower(kx,ky)");
108     for (n = 0; n < rfl; n++)
109     {
110         fprintf(datfile1, "%d\t%d\t %8.6f\n", (n / fy), (n % fy), pspectavg1[n]);
111         fprintf(datfile2, "%d\t%d\t %8.6f\n", (n /fy), (n % fy), pspectavg2[n]);
112     }
113     gmx_ffclose(datfile1);
114     gmx_ffclose(datfile2);
115
116     free(ftspect1);
117     free(ftspect2);
118
119 }
120
121 void powerspectavg_intf(t_interf ***if1, t_interf ***if2, int t, int xb, int yb, char **fnms)
122 {
123     real ***surf;
124
125     int     xy = xb*yb;
126     int     i, n;
127
128     snew(surf, 2);
129     snew(surf[0], t);
130     snew(surf[1], t);
131     for (n = 0; n < t; n++)
132     {
133         snew(surf[0][n], xy);
134         snew(surf[1][n], xy);
135         for (i = 0; i < xy; i++)
136         {
137             surf[0][n][i] = if1[n][i]->Z;
138             surf[1][n][i] = if2[n][i]->Z;
139         }
140     }
141     powerspectavg(surf, t, xb, yb, fnms);
142 }