5f3dbf996047fda1f1e334dead77b633689b2909
[alexxy/gromacs.git] / src / gromacs / ewald / calculate_spline_moduli.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,2016,2017 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "calculate_spline_moduli.h"
41
42 #include <cmath>
43
44 #include <algorithm>
45
46 #include "gromacs/math/utilities.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/utility/exceptions.h"
49 #include "gromacs/utility/gmxassert.h"
50 #include "gromacs/utility/smalloc.h"
51
52 #include "pme_internal.h"
53
54 static void make_dft_mod(real* mod, const double* data, int splineOrder, int ndata)
55 {
56     for (int i = 0; i < ndata; i++)
57     {
58         /* We use double precision, since this is only called once per grid.
59          * But for single precision bsp_mod, single precision also seems
60          * to give full accuracy.
61          */
62         double sc = 0;
63         double ss = 0;
64         for (int j = 0; j < splineOrder; j++)
65         {
66             double arg = (2.0 * M_PI * i * (j + 1)) / ndata;
67             sc += data[j] * cos(arg);
68             ss += data[j] * sin(arg);
69         }
70         mod[i] = sc * sc + ss * ss;
71     }
72     if (splineOrder % 2 == 0 && ndata % 2 == 0)
73     {
74         /* Note that pme_order = splineOrder + 1 */
75         GMX_RELEASE_ASSERT(mod[ndata / 2] < GMX_DOUBLE_EPS,
76                            "With even spline order and even grid size (ndata), dft_mod[ndata/2] "
77                            "should first come out as zero");
78         /* This factor causes a division by zero. But since this occurs in
79          * the tail of the distribution, the term with this factor can
80          * be ignored (see Essmann et al. JCP 103, 8577).
81          * Using the average of the neighbors probably originates from
82          * Tom Darden's original PME code. It seems to give slighlty better
83          * accuracy than using a large value.
84          */
85         mod[ndata / 2] = (mod[ndata / 2 - 1] + mod[ndata / 2 + 1]) * 0.5;
86     }
87 }
88
89 void make_bspline_moduli(splinevec bsp_mod, int nx, int ny, int nz, int pme_order)
90 {
91     /* We use double precision, since this is only called once per grid.
92      * But for single precision bsp_mod, single precision also seems
93      * to give full accuracy.
94      */
95     double* data;
96
97     /* In GROMACS we, confusingly, defined pme-order as the order
98      * of the cardinal B-spline + 1. This probably happened because
99      * the smooth PME paper only talks about "n" which is the number
100      * of points we spread to and that was chosen to be pme-order.
101      */
102     const int splineOrder = pme_order - 1;
103
104     snew(data, splineOrder);
105
106     data[0] = 1;
107     for (int k = 1; k < splineOrder; k++)
108     {
109         data[k] = 0;
110     }
111
112     for (int k = 2; k <= splineOrder; k++)
113     {
114         double div = 1.0 / k;
115         for (int m = k - 1; m > 0; m--)
116         {
117             data[m] = div * ((k - m) * data[m - 1] + (m + 1) * data[m]);
118         }
119         data[0] = div * data[0];
120     }
121
122     make_dft_mod(bsp_mod[XX], data, splineOrder, nx);
123     make_dft_mod(bsp_mod[YY], data, splineOrder, ny);
124     make_dft_mod(bsp_mod[ZZ], data, splineOrder, nz);
125
126     sfree(data);
127 }
128
129 /* Return the P3M optimal influence function */
130 static double do_p3m_influence(double z, int order)
131 {
132     double z2, z4;
133
134     z2 = z * z;
135     z4 = z2 * z2;
136
137     /* The formula and most constants can be found in:
138      * Ballenegger et al., JCTC 8, 936 (2012)
139      */
140     switch (order)
141     {
142         case 2: return 1.0 - 2.0 * z2 / 3.0;
143         case 3: return 1.0 - z2 + 2.0 * z4 / 15.0;
144         case 4: return 1.0 - 4.0 * z2 / 3.0 + 2.0 * z4 / 5.0 + 4.0 * z2 * z4 / 315.0;
145         case 5:
146             return 1.0 - 5.0 * z2 / 3.0 + 7.0 * z4 / 9.0 - 17.0 * z2 * z4 / 189.0 + 2.0 * z4 * z4 / 2835.0;
147         case 6:
148             return 1.0 - 2.0 * z2 + 19.0 * z4 / 15.0 - 256.0 * z2 * z4 / 945.0
149                    + 62.0 * z4 * z4 / 4725.0 + 4.0 * z2 * z4 * z4 / 155925.0;
150         case 7:
151             return 1.0 - 7.0 * z2 / 3.0 + 28.0 * z4 / 15.0 - 16.0 * z2 * z4 / 27.0 + 26.0 * z4 * z4 / 405.0
152                    - 2.0 * z2 * z4 * z4 / 1485.0 + 4.0 * z4 * z4 * z4 / 6081075.0;
153         case 8:
154             return 1.0 - 8.0 * z2 / 3.0 + 116.0 * z4 / 45.0 - 344.0 * z2 * z4 / 315.0
155                    + 914.0 * z4 * z4 / 4725.0 - 248.0 * z4 * z4 * z2 / 22275.0
156                    + 21844.0 * z4 * z4 * z4 / 212837625.0 - 8.0 * z4 * z4 * z4 * z2 / 638512875.0;
157     }
158
159     return 0.0;
160 }
161
162 /* Calculate the P3M B-spline moduli for one dimension */
163 static void make_p3m_bspline_moduli_dim(real* bsp_mod, int n, int order)
164 {
165     double zarg, zai, sinzai, infl;
166     int    maxk, i;
167
168     if (order > 8)
169     {
170         GMX_THROW(gmx::InconsistentInputError("The current P3M code only supports orders up to 8"));
171     }
172
173     zarg = M_PI / n;
174
175     maxk = (n + 1) / 2;
176
177     for (i = -maxk; i < 0; i++)
178     {
179         zai            = zarg * i;
180         sinzai         = sin(zai);
181         infl           = do_p3m_influence(sinzai, order);
182         bsp_mod[n + i] = infl * infl * pow(sinzai / zai, -2.0 * order);
183     }
184     bsp_mod[0] = 1.0;
185     for (i = 1; i < maxk; i++)
186     {
187         zai        = zarg * i;
188         sinzai     = sin(zai);
189         infl       = do_p3m_influence(sinzai, order);
190         bsp_mod[i] = infl * infl * pow(sinzai / zai, -2.0 * order);
191     }
192 }
193
194 /* Calculate the P3M B-spline moduli */
195 void make_p3m_bspline_moduli(splinevec bsp_mod, int nx, int ny, int nz, int order)
196 {
197     make_p3m_bspline_moduli_dim(bsp_mod[XX], nx, order);
198     make_p3m_bspline_moduli_dim(bsp_mod[YY], ny, order);
199     make_p3m_bspline_moduli_dim(bsp_mod[ZZ], nz, order);
200 }