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, 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.
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.
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.
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.
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.
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.
39 #include "gromacs/legacyheaders/calcgrid.h"
43 #include "gromacs/legacyheaders/typedefs.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/utility/smalloc.h"
47 /* The grid sizes below are based on timing of a 3D cubic grid in fftw
48 * compiled with SSE using 4 threads in fft5d.c.
49 * A grid size is removed when a larger grid is faster.
52 /* Small grid size array */
54 const int grid_init
[g_initNR
] = { 6, 8, 10, 12, 14, 16, 20, 24, 25, 28, 32, 36, 40, 42, 44 };
56 /* For larger grid sizes, a prefactor with any power of 2 can be added.
57 * Only sizes divisible by 4 should be used, 90 is allowed, 140 not.
60 const int grid_base
[g_baseNR
] = { 45, 48, 50, 52, 54, 56, 60, 64, 70, 72, 75, 80, 81, 84 };
62 real
calc_grid(FILE *fp
, matrix box
, real gr_sp
,
63 int *nx
, int *ny
, int *nz
)
72 if ((*nx
<= 0 || *ny
<= 0 || *nz
<= 0) && gr_sp
<= 0)
74 gmx_fatal(FARGS
, "invalid fourier grid spacing: %g", gr_sp
);
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
] = sqrt(box_size
[d
]);
108 if ((*nx
<= 0) || (*ny
<= 0) || (*nz
<= 0))
112 fprintf(fp
, "Calculating fourier grid dimensions for%s%s%s\n",
113 *nx
> 0 ? "" : " X", *ny
> 0 ? "" : " Y", *nz
> 0 ? "" : " Z");
118 for (d
= 0; d
< DIM
; d
++)
122 nmin
= (int)(box_size
[d
]/gr_sp
+ 0.999);
125 if (grid_init
[i
] >= nmin
)
127 /* Take the smallest possible grid in the list */
128 while (i
> 0 && grid_init
[i
-1] >= nmin
)
136 /* Determine how many pre-factors of 2 we need */
139 while (fac2
*grid_base
[i
] < nmin
)
143 /* Find the smallest grid that is >= nmin */
146 try = fac2
*grid_base
[i
];
147 /* We demand a factor of 4, avoid 140, allow 90 */
148 if (((try % 4 == 0 && try != 140) || try == 90) &&
159 spacing
[d
] = box_size
[d
]/n
[d
];
160 if (spacing
[d
] > max_spacing
)
162 max_spacing
= spacing
[d
];
170 fprintf(fp
, "Using a fourier grid of %dx%dx%d, spacing %.3f %.3f %.3f\n",
171 *nx
, *ny
, *nz
, spacing
[XX
], spacing
[YY
], spacing
[ZZ
]);