0711b8c3a1a762736cbf1c45efc753af67c8c1d0
[alexxy/gromacs.git] / src / tools / dens_filter.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2001
5  * BIOSON Research Institute, Dept. of Biophysical Chemistry
6  * University of Groningen, The Netherlands
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42 /* dens_filter.c
43  * Routines for Filters and convolutions
44  */
45
46 #include <math.h>
47 #include "typedefs.h"
48 #include "dens_filter.h"
49 #include "smalloc.h"
50 #include "vec.h"
51
52 #ifdef GMX_DOUBLE
53 #define EXP(x) (exp(x))
54 #else
55 #define EXP(x) (expf(x))
56 #endif
57
58 /* Modulo function assuming y>0 but with arbitrary INTEGER x */
59 static int MOD(int x, int y){
60     if (x<0) x=x+y;
61     return (mod(x,y));
62 }
63
64
65 gmx_bool convolution(int dataSize,real *x, int kernelSize, real* kernel)
66 {
67     int i, j, k;
68     real *out;
69     snew(out,dataSize);
70     /* check validity of params */
71     if(!x || !kernel) return FALSE;
72     if(dataSize <=0 || kernelSize <= 0) return FALSE;
73
74     /* start convolution from out[kernelSize-1] to out[dataSize-1] (last) */
75     for(i = kernelSize-1; i < dataSize; ++i)
76     {
77         for(j = i, k = 0; k < kernelSize; --j, ++k)
78             out[i] += x[j] * kernel[k];
79     }
80
81     /* convolution from out[0] to out[kernelSize-2] */
82     for(i = 0; i < kernelSize - 1; ++i)
83     {
84         for(j = i, k = 0; j >= 0; --j, ++k)
85             out[i] += x[j] * kernel[k];
86     }
87
88     for(i=0;i<dataSize;i++) x[i]=out[i];
89     sfree(out); 
90     return TRUE;
91 }
92
93 /* Assuming kernel is shorter than x */
94
95 gmx_bool periodic_convolution(int datasize, real *x, int kernelsize, 
96                               real *kernel)
97 {
98     int i,j,k,nj;
99     real *filtered;
100
101     if (!x || !kernel) return FALSE;
102     if (kernelsize<=0|| datasize<=0|| kernelsize > datasize) return FALSE;
103
104     snew(filtered,datasize);
105     
106     for(i=0;(i<datasize); i++){
107         for(j=0; (j<kernelsize); j++){
108             filtered[i] += kernel[j]*x[MOD((i-j),datasize)];
109         }
110     }
111     for(i=0; i<datasize; i++) x[i]=filtered[i];
112     sfree(filtered);
113
114     return TRUE;
115 }
116
117
118 /* 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
119  * NO checks are performed 
120  */
121 void gausskernel(real *out, int n, real var){
122         int i,j=0,k;
123         real arg,tot=0;
124         k=n/2;
125         
126         for(i=-k; i<=k; i++){
127                 arg=(i*i)/(2*var);
128                 tot+=out[j++]=EXP(-arg);
129         }
130         for(i=0;i<j;i++){
131                 out[i]/=tot;
132         }
133 }
134
135