2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2013,2014,2015,2017,2018, 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.
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.
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.
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.
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.
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.
38 #include "powerspect.h"
40 #include "gromacs/fft/fft.h"
41 #include "gromacs/gmxana/interf.h"
42 #include "gromacs/utility/fatalerror.h"
43 #include "gromacs/utility/futil.h"
44 #include "gromacs/utility/real.h"
45 #include "gromacs/utility/smalloc.h"
47 static void addtoavgenergy(t_complex
*list
, real
*result
, int size
, int tsteps
)
50 for (i
= 0; i
< size
; i
++)
52 result
[i
] += cabs2(list
[i
])/tsteps
;
58 void powerspectavg(real
***intftab
, int tsteps
, int xbins
, int ybins
,
59 gmx::ArrayRef
<const std::string
> outfiles
)
61 /*Fourier plans and output;*/
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 */
65 real
*pspectavg1
; /*power -spectrum 1st interface*/
66 real
*pspectavg2
; /* -------------- 2nd interface*/
68 FILE *datfile1
, *datfile2
; /*data-files with interface data*/
70 int fy
= ybins
/2+1; /* number of (symmetric) fourier y elements; */
71 int rfl
= xbins
*fy
; /*length of real - DFT == Symmetric 2D matrix*/
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)
76 gmx_fatal(FARGS
, "Error allocating FFT");
83 snew(pspectavg1
, rfl
);
84 snew(pspectavg2
, rfl
);
86 /*Fouriertransform directly (no normalization or anything)*/
87 /*NB! Check carefully indexes here*/
89 for (n
= 0; n
< tsteps
; n
++)
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
);
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
);
99 /*Print out average energy-spectrum to outfiles[0] and outfiles[1];*/
101 datfile1
= gmx_ffopen(outfiles
[0], "w");
102 datfile2
= gmx_ffopen(outfiles
[1], "w");
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
++)
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
]);
112 gmx_ffclose(datfile1
);
113 gmx_ffclose(datfile2
);
120 void powerspectavg_intf(t_interf
***if1
, t_interf
***if2
, int t
, int xb
, int yb
,
121 gmx::ArrayRef
<const std::string
> outfiles
)
131 for (n
= 0; n
< t
; n
++)
133 snew(surf
[0][n
], xy
);
134 snew(surf
[1][n
], xy
);
135 for (i
= 0; i
< xy
; i
++)
137 surf
[0][n
][i
] = if1
[n
][i
]->Z
;
138 surf
[1][n
][i
] = if2
[n
][i
]->Z
;
141 powerspectavg(surf
, t
, xb
, yb
, outfiles
);