2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2011,2013,2014,2015,2018 by the GROMACS development team.
5 * Copyright (c) 2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
40 * Routines for Filters and convolutions
43 #include "dens_filter.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/utility/smalloc.h"
50 gmx_bool
convolution(int dataSize
, real
* x
, int kernelSize
, const real
* kernel
)
55 /* check validity of params */
60 if (dataSize
<= 0 || kernelSize
<= 0)
65 /* start convolution from out[kernelSize-1] to out[dataSize-1] (last) */
66 for (i
= kernelSize
- 1; i
< dataSize
; ++i
)
68 for (j
= i
, k
= 0; k
< kernelSize
; --j
, ++k
)
70 out
[i
] += x
[j
] * kernel
[k
];
74 /* convolution from out[0] to out[kernelSize-2] */
75 for (i
= 0; i
< kernelSize
- 1; ++i
)
77 for (j
= i
, k
= 0; j
>= 0; --j
, ++k
)
79 out
[i
] += x
[j
] * kernel
[k
];
83 for (i
= 0; i
< dataSize
; i
++)
91 /* Assuming kernel is shorter than x */
93 gmx_bool
periodic_convolution(int datasize
, real
* x
, int kernelsize
, const real
* kernel
)
102 if (kernelsize
<= 0 || datasize
<= 0 || kernelsize
> datasize
)
107 snew(filtered
, datasize
);
109 for (i
= 0; (i
< datasize
); i
++)
111 for (j
= 0; (j
< kernelsize
); j
++)
113 // add datasize in case i-j is <0
114 idx
= i
- j
+ datasize
;
115 filtered
[i
] += kernel
[j
] * x
[idx
% datasize
];
118 for (i
= 0; i
< datasize
; i
++)
128 /* 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
129 * order NO checks are performed
131 void gausskernel(real
* out
, int n
, real var
)
137 for (i
= -k
; i
<= k
; i
++)
139 arg
= (i
* i
) / (2 * var
);
140 tot
+= out
[j
++] = std::exp(-arg
);
142 for (i
= 0; i
< j
; i
++)