3e2f6b2602a0b6239173f40aa7a43a07f4ca76f0
[alexxy/gromacs.git] / src / gromacs / gmxana / dens_filter.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2011,2013,2014,2015, 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 <cmath>
45
46 #include "gromacs/legacyheaders/typedefs.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/utility/smalloc.h"
49
50 gmx_bool convolution(int dataSize, real *x, int kernelSize, real* kernel)
51 {
52     int   i, j, k;
53     real *out;
54     snew(out, dataSize);
55     /* check validity of params */
56     if (!x || !kernel)
57     {
58         return FALSE;
59     }
60     if (dataSize <= 0 || kernelSize <= 0)
61     {
62         return FALSE;
63     }
64
65     /* start convolution from out[kernelSize-1] to out[dataSize-1] (last) */
66     for (i = kernelSize-1; i < dataSize; ++i)
67     {
68         for (j = i, k = 0; k < kernelSize; --j, ++k)
69         {
70             out[i] += x[j] * kernel[k];
71         }
72     }
73
74     /* convolution from out[0] to out[kernelSize-2] */
75     for (i = 0; i < kernelSize - 1; ++i)
76     {
77         for (j = i, k = 0; j >= 0; --j, ++k)
78         {
79             out[i] += x[j] * kernel[k];
80         }
81     }
82
83     for (i = 0; i < dataSize; i++)
84     {
85         x[i] = out[i];
86     }
87     sfree(out);
88     return TRUE;
89 }
90
91 /* Assuming kernel is shorter than x */
92
93 gmx_bool periodic_convolution(int datasize, real *x, int kernelsize,
94                               real *kernel)
95 {
96     int   i, j, idx;
97     real *filtered;
98
99     if (!x || !kernel)
100     {
101         return FALSE;
102     }
103     if (kernelsize <= 0 || datasize <= 0 || kernelsize > datasize)
104     {
105         return FALSE;
106     }
107
108     snew(filtered, datasize);
109
110     for (i = 0; (i < datasize); i++)
111     {
112         for (j = 0; (j < kernelsize); j++)
113         {
114             // add datasize in case i-j is <0
115             idx          = i-j + datasize;
116             filtered[i] += kernel[j]*x[idx % datasize];
117         }
118     }
119     for (i = 0; i < datasize; i++)
120     {
121         x[i] = filtered[i];
122     }
123     sfree(filtered);
124
125     return TRUE;
126 }
127
128
129 /* 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
130  * NO checks are performed
131  */
132 void gausskernel(real *out, int n, real var)
133 {
134     int  i, j = 0, k;
135     real arg, tot = 0;
136     k = n/2;
137
138     for (i = -k; i <= k; i++)
139     {
140         arg  = (i*i)/(2*var);
141         tot += out[j++] = std::exp(-arg);
142     }
143     for (i = 0; i < j; i++)
144     {
145         out[i] /= tot;
146     }
147 }