2 * This file is part of the GROMACS molecular simulation package.
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 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
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.
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.
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.
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.
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.
39 * Implement routines for integrating a data set
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
45 #include "integrate.h"
50 #include "gromacs/math/functions.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/utility/fatalerror.h"
54 /*! \brief Integrate a function and printe the integral value. */
55 real print_and_integrate(FILE* fp, int n, real dt, const real c[], const real* fit, int nskip)
60 /* Use trapezoidal rule for calculating integral */
62 for (j = 0; (j < n); j++)
65 if (fp && (nskip == 0 || j % nskip == 0))
67 fprintf(fp, "%10.3f %10.5f\n", j * dt, c0);
71 sum += dt * (c0 + c[j - 1]);
79 for (j = 0; (j < n); j++)
81 if (nskip == 0 || j % nskip == 0)
83 fprintf(fp, "%10.3f %10.5f\n", j * dt, fit[j]);
92 /*! \brief Compute and return the integral of a function. */
93 real evaluate_integral(int n, const real x[], const real y[], const real dy[], real aver_start, real* stddev)
95 double sum, sum_var, w;
96 double sum_tail = 0, sum2_tail = 0;
99 /* Use trapezoidal rule for calculating integral */
102 gmx_fatal(FARGS, "Evaluating integral: n = %d (file %s, line %d)", n, __FILE__, __LINE__);
107 for (j = 0; (j < n); j++)
112 w += 0.5 * (x[j] - x[j - 1]);
116 w += 0.5 * (x[j + 1] - x[j]);
121 /* Assume all errors are uncorrelated */
122 sum_var += gmx::square(w * dy[j]);
125 if ((aver_start > 0) && (x[j] >= aver_start))
128 sum2_tail += std::sqrt(sum_var);
135 sum = sum_tail / nsum_tail;
136 /* This is a worst case estimate, assuming all stddev's are correlated. */
137 *stddev = sum2_tail / nsum_tail;
141 *stddev = std::sqrt(sum_var);