840f8166754de073b1a1439b01ac348f6dc79c4f
[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,2018,2019, 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 <cmath>
47 #include <cstdio>
48
49 #include "gromacs/math/functions.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/utility/fatalerror.h"
52
53 /*! \brief Integrate a function and printe the integral value. */
54 real print_and_integrate(FILE* fp, int n, real dt, const real c[], const real* fit, int nskip)
55 {
56     real c0, sum;
57     int  j;
58
59     /* Use trapezoidal rule for calculating integral */
60     sum = 0.0;
61     for (j = 0; (j < n); j++)
62     {
63         c0 = c[j];
64         if (fp && (nskip == 0 || j % nskip == 0))
65         {
66             fprintf(fp, "%10.3f  %10.5f\n", j * dt, c0);
67         }
68         if (j > 0)
69         {
70             sum += dt * (c0 + c[j - 1]);
71         }
72     }
73     if (fp)
74     {
75         fprintf(fp, "&\n");
76         if (fit)
77         {
78             for (j = 0; (j < n); j++)
79             {
80                 if (nskip == 0 || j % nskip == 0)
81                 {
82                     fprintf(fp, "%10.3f  %10.5f\n", j * dt, fit[j]);
83                 }
84             }
85             fprintf(fp, "&\n");
86         }
87     }
88     return sum * 0.5;
89 }
90
91 /*! \brief Compute and return the integral of a function. */
92 real evaluate_integral(int n, const real x[], const real y[], const real dy[], real aver_start, real* stddev)
93 {
94     double sum, sum_var, w;
95     double sum_tail = 0, sum2_tail = 0;
96     int    j, nsum_tail            = 0;
97
98     /* Use trapezoidal rule for calculating integral */
99     if (n <= 0)
100     {
101         gmx_fatal(FARGS, "Evaluating integral: n = %d (file %s, line %d)", n, __FILE__, __LINE__);
102     }
103
104     sum     = 0;
105     sum_var = 0;
106     for (j = 0; (j < n); j++)
107     {
108         w = 0;
109         if (j > 0)
110         {
111             w += 0.5 * (x[j] - x[j - 1]);
112         }
113         if (j < n - 1)
114         {
115             w += 0.5 * (x[j + 1] - x[j]);
116         }
117         sum += w * y[j];
118         if (dy)
119         {
120             /* Assume all errors are uncorrelated */
121             sum_var += gmx::square(w * dy[j]);
122         }
123
124         if ((aver_start > 0) && (x[j] >= aver_start))
125         {
126             sum_tail += sum;
127             sum2_tail += std::sqrt(sum_var);
128             nsum_tail += 1;
129         }
130     }
131
132     if (nsum_tail > 0)
133     {
134         sum = sum_tail / nsum_tail;
135         /* This is a worst case estimate, assuming all stddev's are correlated. */
136         *stddev = sum2_tail / nsum_tail;
137     }
138     else
139     {
140         *stddev = std::sqrt(sum_var);
141     }
142
143     return sum;
144 }