Remove unnecessary config.h includes
[alexxy/gromacs.git] / src / gromacs / gmxana / dens_filter.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 #include "gmxpre.h"
37
38 /* dens_filter.c
39  * Routines for Filters and convolutions
40  */
41
42 #include <math.h>
43 #include "gromacs/legacyheaders/typedefs.h"
44 #include "dens_filter.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/math/vec.h"
47
48 #ifdef GMX_DOUBLE
49 #define EXP(x) (exp(x))
50 #else
51 #define EXP(x) (expf(x))
52 #endif
53
54 /* Modulo function assuming y>0 but with arbitrary INTEGER x */
55 static int MOD(int x, int y)
56 {
57     if (x < 0)
58     {
59         x = x+y;
60     }
61     return (mod(x, y));
62 }
63
64
65 gmx_bool convolution(int dataSize, real *x, int kernelSize, real* kernel)
66 {
67     int   i, j, k;
68     real *out;
69     snew(out, dataSize);
70     /* check validity of params */
71     if (!x || !kernel)
72     {
73         return FALSE;
74     }
75     if (dataSize <= 0 || kernelSize <= 0)
76     {
77         return FALSE;
78     }
79
80     /* start convolution from out[kernelSize-1] to out[dataSize-1] (last) */
81     for (i = kernelSize-1; i < dataSize; ++i)
82     {
83         for (j = i, k = 0; k < kernelSize; --j, ++k)
84         {
85             out[i] += x[j] * kernel[k];
86         }
87     }
88
89     /* convolution from out[0] to out[kernelSize-2] */
90     for (i = 0; i < kernelSize - 1; ++i)
91     {
92         for (j = i, k = 0; j >= 0; --j, ++k)
93         {
94             out[i] += x[j] * kernel[k];
95         }
96     }
97
98     for (i = 0; i < dataSize; i++)
99     {
100         x[i] = out[i];
101     }
102     sfree(out);
103     return TRUE;
104 }
105
106 /* Assuming kernel is shorter than x */
107
108 gmx_bool periodic_convolution(int datasize, real *x, int kernelsize,
109                               real *kernel)
110 {
111     int   i, j, k, nj;
112     real *filtered;
113
114     if (!x || !kernel)
115     {
116         return FALSE;
117     }
118     if (kernelsize <= 0 || datasize <= 0 || kernelsize > datasize)
119     {
120         return FALSE;
121     }
122
123     snew(filtered, datasize);
124
125     for (i = 0; (i < datasize); i++)
126     {
127         for (j = 0; (j < kernelsize); j++)
128         {
129             filtered[i] += kernel[j]*x[MOD((i-j), datasize)];
130         }
131     }
132     for (i = 0; i < datasize; i++)
133     {
134         x[i] = filtered[i];
135     }
136     sfree(filtered);
137
138     return TRUE;
139 }
140
141
142 /* returns discrete gaussian kernel of size n in *out, where n=2k+1=3,5,7,9,11 and k=1,2,3 is the order
143  * NO checks are performed
144  */
145 void gausskernel(real *out, int n, real var)
146 {
147     int  i, j = 0, k;
148     real arg, tot = 0;
149     k = n/2;
150
151     for (i = -k; i <= k; i++)
152     {
153         arg  = (i*i)/(2*var);
154         tot += out[j++] = EXP(-arg);
155     }
156     for (i = 0; i < j; i++)
157     {
158         out[i] /= tot;
159     }
160 }