2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
38 * Routines for Filters and convolutions
43 #include "dens_filter.h"
44 #include "gromacs/utility/smalloc.h"
45 #include "gromacs/math/vec.h"
48 #define EXP(x) (exp(x))
50 #define EXP(x) (expf(x))
53 /* Modulo function assuming y>0 but with arbitrary INTEGER x */
54 static int MOD(int x, int y)
64 gmx_bool convolution(int dataSize, real *x, int kernelSize, real* kernel)
69 /* check validity of params */
74 if (dataSize <= 0 || kernelSize <= 0)
79 /* start convolution from out[kernelSize-1] to out[dataSize-1] (last) */
80 for (i = kernelSize-1; i < dataSize; ++i)
82 for (j = i, k = 0; k < kernelSize; --j, ++k)
84 out[i] += x[j] * kernel[k];
88 /* convolution from out[0] to out[kernelSize-2] */
89 for (i = 0; i < kernelSize - 1; ++i)
91 for (j = i, k = 0; j >= 0; --j, ++k)
93 out[i] += x[j] * kernel[k];
97 for (i = 0; i < dataSize; i++)
105 /* Assuming kernel is shorter than x */
107 gmx_bool periodic_convolution(int datasize, real *x, int kernelsize,
117 if (kernelsize <= 0 || datasize <= 0 || kernelsize > datasize)
122 snew(filtered, datasize);
124 for (i = 0; (i < datasize); i++)
126 for (j = 0; (j < kernelsize); j++)
128 filtered[i] += kernel[j]*x[MOD((i-j), datasize)];
131 for (i = 0; i < datasize; i++)
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
144 void gausskernel(real *out, int n, real var)
150 for (i = -k; i <= k; i++)
153 tot += out[j++] = EXP(-arg);
155 for (i = 0; i < j; i++)