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