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