5bb256f6868e794d6b58dc2d19fb5406bb9e5ab5
[alexxy/gromacs.git] / src / gromacs / correlationfunctions / integrate.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2017, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /*! \brief
38  * Implement routines for integrating a data set
39  *
40  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41  */
42 #include "gmxpre.h"
43
44 #include "integrate.h"
45
46 #include <stdio.h>
47
48 #include <cmath>
49
50 #include "gromacs/math/functions.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/utility/fatalerror.h"
53
54 /*! \brief Integrate a function and printe the integral value. */
55 real print_and_integrate(FILE *fp, int n, real dt, const real c[],
56                          const real *fit, int nskip)
57 {
58     real c0, sum;
59     int  j;
60
61     /* Use trapezoidal rule for calculating integral */
62     sum = 0.0;
63     for (j = 0; (j < n); j++)
64     {
65         c0 = c[j];
66         if (fp && (nskip == 0 || j % nskip == 0))
67         {
68             fprintf(fp, "%10.3f  %10.5f\n", j*dt, c0);
69         }
70         if (j > 0)
71         {
72             sum += dt*(c0+c[j-1]);
73         }
74     }
75     if (fp)
76     {
77         fprintf(fp, "&\n");
78         if (fit)
79         {
80             for (j = 0; (j < n); j++)
81             {
82                 if (nskip == 0 || j % nskip == 0)
83                 {
84                     fprintf(fp, "%10.3f  %10.5f\n", j*dt, fit[j]);
85                 }
86             }
87             fprintf(fp, "&\n");
88         }
89     }
90     return sum*0.5;
91 }
92
93 /*! \brief Compute and return the integral of a function. */
94 real evaluate_integral(int n, const real x[], const real y[],
95                        const real dy[], real aver_start,
96                        real *stddev)
97 {
98     double sum, sum_var, w;
99     double sum_tail = 0, sum2_tail = 0;
100     int    j, nsum_tail = 0;
101
102     /* Use trapezoidal rule for calculating integral */
103     if (n <= 0)
104     {
105         gmx_fatal(FARGS, "Evaluating integral: n = %d (file %s, line %d)",
106                   n, __FILE__, __LINE__);
107     }
108
109     sum     = 0;
110     sum_var = 0;
111     for (j = 0; (j < n); j++)
112     {
113         w = 0;
114         if (j > 0)
115         {
116             w += 0.5*(x[j] - x[j-1]);
117         }
118         if (j < n-1)
119         {
120             w += 0.5*(x[j+1] - x[j]);
121         }
122         sum += w*y[j];
123         if (dy)
124         {
125             /* Assume all errors are uncorrelated */
126             sum_var += gmx::square(w*dy[j]);
127         }
128
129         if ((aver_start > 0) && (x[j] >= aver_start))
130         {
131             sum_tail  += sum;
132             sum2_tail += std::sqrt(sum_var);
133             nsum_tail += 1;
134         }
135     }
136
137     if (nsum_tail > 0)
138     {
139         sum = sum_tail/nsum_tail;
140         /* This is a worst case estimate, assuming all stddev's are correlated. */
141         *stddev = sum2_tail/nsum_tail;
142     }
143     else
144     {
145         *stddev = std::sqrt(sum_var);
146     }
147
148     return sum;
149 }