Merge branch 'master' into pygromacs
[alexxy/gromacs.git] / src / gromacs / gmxlib / calcgrid.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) 2012,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.
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 #include "gmxpre.h"
38
39 #include "gromacs/legacyheaders/calcgrid.h"
40
41 #include <cmath>
42
43 #include <algorithm>
44
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/smalloc.h"
48
49 /* The grid sizes below are based on timing of a 3D cubic grid in fftw
50  * compiled with SSE using 4 threads in fft5d.c.
51  * A grid size is removed when a larger grid is faster.
52  */
53
54 /* Small grid size array */
55 #define g_initNR 15
56 const int grid_init[g_initNR] = { 6, 8, 10, 12, 14, 16, 20, 24, 25, 28, 32, 36, 40, 42, 44 };
57
58 /* For larger grid sizes, a prefactor with any power of 2 can be added.
59  * Only sizes divisible by 4 should be used, 90 is allowed, 140 not.
60  */
61 #define g_baseNR 14
62 const int grid_base[g_baseNR] = { 45, 48, 50, 52, 54, 56, 60, 64, 70, 72, 75, 80, 81, 84 };
63
64 real calc_grid(FILE *fp, matrix box, real gr_sp,
65                int *nx, int *ny, int *nz)
66 {
67     int  d, n[DIM];
68     int  i;
69     rvec box_size;
70     int  nmin, fac2, attempt;
71     rvec spacing;
72     real max_spacing;
73
74     if ((*nx <= 0 || *ny <= 0 || *nz <= 0) && gr_sp <= 0)
75     {
76         gmx_fatal(FARGS, "invalid fourier grid spacing: %g", gr_sp);
77     }
78
79     if (grid_base[g_baseNR-1] % 4 != 0)
80     {
81         gmx_incons("the last entry in grid_base is not a multiple of 4");
82     }
83
84     /* New grid calculation setup:
85      *
86      * To maintain similar accuracy for triclinic PME grids as for rectangular
87      * ones, the max grid spacing should set along the box vectors rather than
88      * cartesian X/Y/Z directions. This will lead to slightly larger grids, but
89      * it is much better than having to go to pme_order=6.
90      *
91      * Thus, instead of just extracting the diagonal elements to box_size[d], we
92      * now calculate the cartesian length of the vectors.
93      *
94      * /Erik Lindahl, 20060402.
95      */
96     for (d = 0; d < DIM; d++)
97     {
98         box_size[d] = 0;
99         for (i = 0; i < DIM; i++)
100         {
101             box_size[d] += box[d][i]*box[d][i];
102         }
103         box_size[d] = std::sqrt(box_size[d]);
104     }
105
106     n[XX] = *nx;
107     n[YY] = *ny;
108     n[ZZ] = *nz;
109
110     if ((*nx <= 0) || (*ny <= 0) || (*nz <= 0))
111     {
112         if (NULL != fp)
113         {
114             fprintf(fp, "Calculating fourier grid dimensions for%s%s%s\n",
115                     *nx > 0 ? "" : " X", *ny > 0 ? "" : " Y", *nz > 0 ? "" : " Z");
116         }
117     }
118
119     max_spacing = 0;
120     for (d = 0; d < DIM; d++)
121     {
122         if (n[d] <= 0)
123         {
124             nmin = static_cast<int>(box_size[d]/gr_sp + 0.999);
125
126             i = g_initNR - 1;
127             if (grid_init[i] >= nmin)
128             {
129                 /* Take the smallest possible grid in the list */
130                 while (i > 0 && grid_init[i-1] >= nmin)
131                 {
132                     i--;
133                 }
134                 n[d] = grid_init[i];
135             }
136             else
137             {
138                 /* Determine how many pre-factors of 2 we need */
139                 fac2 = 1;
140                 i    = g_baseNR - 1;
141                 while (fac2*grid_base[i] < nmin)
142                 {
143                     fac2 *= 2;
144                 }
145                 /* Find the smallest grid that is >= nmin */
146                 do
147                 {
148                     attempt = fac2*grid_base[i];
149                     /* We demand a factor of 4, avoid 140, allow 90 */
150                     if (((attempt % 4 == 0 && attempt != 140) || attempt == 90) &&
151                         attempt >= nmin)
152                     {
153                         n[d] = attempt;
154                     }
155                     i--;
156                 }
157                 while (i > 0);
158             }
159         }
160
161         spacing[d]  = box_size[d]/n[d];
162         max_spacing = std::max(max_spacing, spacing[d]);
163     }
164     *nx = n[XX];
165     *ny = n[YY];
166     *nz = n[ZZ];
167     if (NULL != fp)
168     {
169         fprintf(fp, "Using a fourier grid of %dx%dx%d, spacing %.3f %.3f %.3f\n",
170                 *nx, *ny, *nz, spacing[XX], spacing[YY], spacing[ZZ]);
171     }
172
173     return max_spacing;
174 }