Add conserved quantity for Berendsen P-couple
[gromacs.git] / src / gromacs / math / calculate-ewald-splitting-coefficient.cpp
blobd653d9ffea553f1bb394a74795ab952b9e9046fe
1 /*
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, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "calculate-ewald-splitting-coefficient.h"
41 #include <cmath>
43 #include "gromacs/math/utilities.h"
44 #include "gromacs/utility/real.h"
46 real calc_ewaldcoeff_q(real rc, real rtol)
48 real beta = 5, low, high;
49 int n, i = 0;
53 i++;
54 beta *= 2;
56 while (std::erfc(beta*rc) > rtol);
58 /* Do a binary search with tolerance 2^-60 */
59 n = i+60;
60 low = 0;
61 high = beta;
62 for (i = 0; i < n; i++)
64 beta = (low+high)/2;
65 if (std::erfc(beta*rc) > rtol)
67 low = beta;
69 else
71 high = beta;
74 return beta;
77 static real compute_lj_function(real beta, real rc)
79 real xrc, xrc2, xrc4, result;
80 xrc = beta*rc;
81 xrc2 = xrc*xrc;
82 xrc4 = xrc2*xrc2;
83 result = std::exp(-xrc2)*(1 + xrc2 + xrc4/2.0);
85 return result;
88 real calc_ewaldcoeff_lj(real rc, real rtol)
90 real beta = 5, low, high;
91 int n, i = 0;
95 i++;
96 beta *= 2.0;
98 while (compute_lj_function(beta, rc) > rtol);
100 /* Do a binary search with tolerance 2^-60 */
101 n = i + 60;
102 low = 0;
103 high = beta;
104 for (i = 0; i < n; ++i)
106 beta = (low + high) / 2.0;
107 if (compute_lj_function(beta, rc) > rtol)
109 low = beta;
111 else
113 high = beta;
116 return beta;