Remove unnecessary includes of arrayref.h
[alexxy/gromacs.git] / src / gromacs / gmxana / powerspect.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2013,2014,2015 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 #include "gmxpre.h"
38
39 #include "powerspect.h"
40
41 #include "gromacs/fft/fft.h"
42 #include "gromacs/gmxana/interf.h"
43 #include "gromacs/utility/arrayref.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/utility/futil.h"
46 #include "gromacs/utility/real.h"
47 #include "gromacs/utility/smalloc.h"
48
49 static void addtoavgenergy(t_complex* list, real* result, int size, int tsteps)
50 {
51     int i;
52     for (i = 0; i < size; i++)
53     {
54         result[i] += cabs2(list[i]) / tsteps;
55     }
56 }
57
58
59 void powerspectavg(real*** intftab, int tsteps, int xbins, int ybins, gmx::ArrayRef<const std::string> outfiles)
60 {
61     /*Fourier plans and output;*/
62     gmx_fft_t  fftp;
63     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 */
64     t_complex* ftspect2;
65     real*      pspectavg1; /*power -spectrum 1st interface*/
66     real*      pspectavg2; /* -------------- 2nd interface*/
67     real*      temp;
68     FILE *     datfile1, *datfile2; /*data-files with interface data*/
69     int        n;                   /*time index*/
70     int        fy  = ybins / 2 + 1; /* number of (symmetric) fourier y elements; */
71     int        rfl = xbins * fy;    /*length of real - DFT == Symmetric 2D matrix*/
72
73     /*Prepare data structures for FFT, with time averaging of power spectrum*/
74     if (gmx_fft_init_2d_real(&fftp, xbins, ybins, GMX_FFT_FLAG_NONE) != 0)
75     {
76         gmx_fatal(FARGS, "Error allocating FFT");
77     }
78
79     /*Initialize arrays*/
80     snew(ftspect1, rfl);
81     snew(ftspect2, rfl);
82     snew(temp, rfl);
83     snew(pspectavg1, rfl);
84     snew(pspectavg2, rfl);
85
86     /*Fouriertransform directly (no normalization or anything)*/
87     /*NB! Check carefully indexes here*/
88
89     for (n = 0; n < tsteps; n++)
90     {
91         gmx_fft_2d_real(fftp, GMX_FFT_REAL_TO_COMPLEX, intftab[0][n], ftspect1);
92         gmx_fft_2d_real(fftp, GMX_FFT_REAL_TO_COMPLEX, intftab[1][n], ftspect2);
93
94         /*Add to average for interface 1 here*/
95         addtoavgenergy(ftspect1, pspectavg1, rfl, tsteps);
96         /*Add to average for interface 2 here*/
97         addtoavgenergy(ftspect2, pspectavg2, rfl, tsteps);
98     }
99     /*Print out average energy-spectrum to outfiles[0] and outfiles[1];*/
100
101     datfile1 = gmx_ffopen(outfiles[0], "w");
102     datfile2 = gmx_ffopen(outfiles[1], "w");
103
104     /*Filling dat files with spectral data*/
105     fprintf(datfile1, "%s\n", "kx\t ky\t\tPower(kx,ky)");
106     fprintf(datfile2, "%s\n", "kx\t ky\t\tPower(kx,ky)");
107     for (n = 0; n < rfl; n++)
108     {
109         fprintf(datfile1, "%d\t%d\t %8.6f\n", (n / fy), (n % fy), pspectavg1[n]);
110         fprintf(datfile2, "%d\t%d\t %8.6f\n", (n / fy), (n % fy), pspectavg2[n]);
111     }
112     gmx_ffclose(datfile1);
113     gmx_ffclose(datfile2);
114
115     sfree(ftspect1);
116     sfree(ftspect2);
117 }
118
119 void powerspectavg_intf(t_interf*** if1, t_interf*** if2, int t, int xb, int yb, gmx::ArrayRef<const std::string> outfiles)
120 {
121     real*** surf;
122
123     int xy = xb * yb;
124     int i, n;
125
126     snew(surf, 2);
127     snew(surf[0], t);
128     snew(surf[1], t);
129     for (n = 0; n < t; n++)
130     {
131         snew(surf[0][n], xy);
132         snew(surf[1][n], xy);
133         for (i = 0; i < xy; i++)
134         {
135             surf[0][n][i] = if1[n][i]->Z;
136             surf[1][n][i] = if2[n][i]->Z;
137         }
138     }
139     powerspectavg(surf, t, xb, yb, outfiles);
140 }