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