Merge "Fix bug in selection subexpression handling."
[alexxy/gromacs.git] / src / gromacs / gmxana / 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 {
59     if (x < 0)
60     {
61         x = x+y;
62     }
63     return (mod(x, y));
64 }
65
66
67 gmx_bool convolution(int dataSize, real *x, int kernelSize, real* kernel)
68 {
69     int   i, j, k;
70     real *out;
71     snew(out, dataSize);
72     /* check validity of params */
73     if (!x || !kernel)
74     {
75         return FALSE;
76     }
77     if (dataSize <= 0 || kernelSize <= 0)
78     {
79         return FALSE;
80     }
81
82     /* start convolution from out[kernelSize-1] to out[dataSize-1] (last) */
83     for (i = kernelSize-1; i < dataSize; ++i)
84     {
85         for (j = i, k = 0; k < kernelSize; --j, ++k)
86         {
87             out[i] += x[j] * kernel[k];
88         }
89     }
90
91     /* convolution from out[0] to out[kernelSize-2] */
92     for (i = 0; i < kernelSize - 1; ++i)
93     {
94         for (j = i, k = 0; j >= 0; --j, ++k)
95         {
96             out[i] += x[j] * kernel[k];
97         }
98     }
99
100     for (i = 0; i < dataSize; i++)
101     {
102         x[i] = out[i];
103     }
104     sfree(out);
105     return TRUE;
106 }
107
108 /* Assuming kernel is shorter than x */
109
110 gmx_bool periodic_convolution(int datasize, real *x, int kernelsize,
111                               real *kernel)
112 {
113     int   i, j, k, nj;
114     real *filtered;
115
116     if (!x || !kernel)
117     {
118         return FALSE;
119     }
120     if (kernelsize <= 0 || datasize <= 0 || kernelsize > datasize)
121     {
122         return FALSE;
123     }
124
125     snew(filtered, datasize);
126
127     for (i = 0; (i < datasize); i++)
128     {
129         for (j = 0; (j < kernelsize); j++)
130         {
131             filtered[i] += kernel[j]*x[MOD((i-j), datasize)];
132         }
133     }
134     for (i = 0; i < datasize; i++)
135     {
136         x[i] = filtered[i];
137     }
138     sfree(filtered);
139
140     return TRUE;
141 }
142
143
144 /* 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
145  * NO checks are performed
146  */
147 void gausskernel(real *out, int n, real var)
148 {
149     int  i, j = 0, k;
150     real arg, tot = 0;
151     k = n/2;
152
153     for (i = -k; i <= k; i++)
154     {
155         arg  = (i*i)/(2*var);
156         tot += out[j++] = EXP(-arg);
157     }
158     for (i = 0; i < j; i++)
159     {
160         out[i] /= tot;
161     }
162 }