d6cb073d6897e15df21a6680c5699ca9b64ae221
[alexxy/gromacs.git] / src / tools / dens_filter.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  * $Id: densorder.c,v 0.9
3  * 
4  *                This source code is part of
5  * 
6  *                 G   R   O   M   A   C   S
7  * 
8  *          GROningen MAchine for Chemical Simulations
9  * 
10  *                        VERSION 3.0
11  * 
12  * Copyright (c) 1991-2001
13  * BIOSON Research Institute, Dept. of Biophysical Chemistry
14  * University of Groningen, The Netherlands
15  * 
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  * 
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  * 
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  * 
31  * Do check out http://www.gromacs.org , or mail us at gromacs@gromacs.org .
32  * 
33  * And Hey:
34  * Gyas ROwers Mature At Cryogenic Speed
35  */
36
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40 /* dens_filter.c
41  * Routines for Filters and convolutions
42  */
43
44 #include <math.h>
45 #include "typedefs.h"
46 #include "dens_filter.h"
47 #include "smalloc.h"
48 #include "vec.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     if (x<0) x=x+y;
59     return (mod(x,y));
60 }
61
62
63 gmx_bool convolution(int dataSize,real *x, int kernelSize, real* kernel)
64 {
65     int i, j, k;
66     real *out;
67     snew(out,dataSize);
68     /* check validity of params */
69     if(!x || !kernel) return FALSE;
70     if(dataSize <=0 || kernelSize <= 0) return FALSE;
71
72     /* start convolution from out[kernelSize-1] to out[dataSize-1] (last) */
73     for(i = kernelSize-1; i < dataSize; ++i)
74     {
75         for(j = i, k = 0; k < kernelSize; --j, ++k)
76             out[i] += x[j] * kernel[k];
77     }
78
79     /* convolution from out[0] to out[kernelSize-2] */
80     for(i = 0; i < kernelSize - 1; ++i)
81     {
82         for(j = i, k = 0; j >= 0; --j, ++k)
83             out[i] += x[j] * kernel[k];
84     }
85
86     for(i=0;i<dataSize;i++) x[i]=out[i];
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,k,nj;
97     real *filtered;
98
99     if (!x || !kernel) return FALSE;
100     if (kernelsize<=0|| datasize<=0|| kernelsize > datasize) return FALSE;
101
102     snew(filtered,datasize);
103     
104     for(i=0;(i<datasize); i++){
105         for(j=0; (j<kernelsize); j++){
106             filtered[i] += kernel[j]*x[MOD((i-j),datasize)];
107         }
108     }
109     for(i=0; i<datasize; i++) x[i]=filtered[i];
110     sfree(filtered);
111
112     return TRUE;
113 }
114
115
116 /* 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
117  * NO checks are performed 
118  */
119 void gausskernel(real *out, int n, real var){
120         int i,j=0,k;
121         real arg,tot=0;
122         k=n/2;
123         
124         for(i=-k; i<=k; i++){
125                 arg=(i*i)/(2*var);
126                 tot+=out[j++]=EXP(-arg);
127         }
128         for(i=0;i<j;i++){
129                 out[i]/=tot;
130         }
131 }
132
133