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