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