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) 2012,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.
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.
46 #include "gromacs/utility/fatalerror.h"
48 /* The grid sizes below are based on timing of a 3D cubic grid in fftw
49 * compiled with SSE using 4 threads in fft5d.c.
50 * A grid size is removed when a larger grid is faster.
53 /* Small grid size array */
55 const int grid_init[g_initNR] = { 6, 8, 10, 12, 14, 16, 20, 24, 25, 28, 32, 36, 40, 42, 44 };
57 /* For larger grid sizes, a prefactor with any power of 2 can be added.
58 * Only sizes divisible by 4 should be used, 90 is allowed, 140 not.
61 const int grid_base[g_baseNR] = { 45, 48, 50, 52, 54, 56, 60, 64, 70, 72, 75, 80, 81, 84 };
63 real calcFftGrid(FILE* fp, const matrix box, real gridSpacing, int minGridPointsPerDim, int* nx, int* ny, int* nz)
68 int nmin, fac2, attempt;
72 if ((*nx <= 0 || *ny <= 0 || *nz <= 0) && gridSpacing <= 0)
74 gmx_fatal(FARGS, "invalid fourier grid spacing: %g", gridSpacing);
77 if (grid_base[g_baseNR - 1] % 4 != 0)
79 gmx_incons("the last entry in grid_base is not a multiple of 4");
82 /* New grid calculation setup:
84 * To maintain similar accuracy for triclinic PME grids as for rectangular
85 * ones, the max grid spacing should set along the box vectors rather than
86 * cartesian X/Y/Z directions. This will lead to slightly larger grids, but
87 * it is much better than having to go to pme_order=6.
89 * Thus, instead of just extracting the diagonal elements to box_size[d], we
90 * now calculate the cartesian length of the vectors.
92 * /Erik Lindahl, 20060402.
94 for (d = 0; d < DIM; d++)
97 for (i = 0; i < DIM; i++)
99 box_size[d] += box[d][i] * box[d][i];
101 box_size[d] = std::sqrt(box_size[d]);
108 if ((*nx <= 0) || (*ny <= 0) || (*nz <= 0))
112 fprintf(fp, "Calculating fourier grid dimensions for%s%s%s\n", *nx > 0 ? "" : " X",
113 *ny > 0 ? "" : " Y", *nz > 0 ? "" : " Z");
118 for (d = 0; d < DIM; d++)
122 nmin = static_cast<int>(box_size[d] / gridSpacing + 0.999);
123 nmin = std::max(nmin, minGridPointsPerDim);
126 if (grid_init[i] >= nmin)
128 /* Take the smallest possible grid in the list */
129 while (i > 0 && grid_init[i - 1] >= nmin)
137 /* Determine how many pre-factors of 2 we need */
140 while (fac2 * grid_base[i] < nmin)
144 /* Find the smallest grid that is >= nmin */
147 attempt = fac2 * grid_base[i];
148 /* We demand a factor of 4, avoid 140, allow 90 */
149 if (((attempt % 4 == 0 && attempt != 140) || attempt == 90) && attempt >= nmin)
158 spacing[d] = box_size[d] / n[d];
159 max_spacing = std::max(max_spacing, spacing[d]);
166 fprintf(fp, "Using a fourier grid of %dx%dx%d, spacing %.3f %.3f %.3f\n", *nx, *ny, *nz,
167 spacing[XX], spacing[YY], spacing[ZZ]);