3 # This file is part of the GROMACS molecular simulation package.
5 # Copyright (c) 2013,2014, by the GROMACS development team, led by
6 # Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 # and including many others, as listed in the AUTHORS file in the
8 # top-level source directory and at http://www.gromacs.org.
10 # GROMACS is free software; you can redistribute it and/or
11 # modify it under the terms of the GNU Lesser General Public License
12 # as published by the Free Software Foundation; either version 2.1
13 # of the License, or (at your option) any later version.
15 # GROMACS is distributed in the hope that it will be useful,
16 # but WITHOUT ANY WARRANTY; without even the implied warranty of
17 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 # Lesser General Public License for more details.
20 # You should have received a copy of the GNU Lesser General Public
21 # License along with GROMACS; if not, see
22 # http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 # Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 # If you want to redistribute modifications to GROMACS, please
26 # consider that scientific software is very special. Version
27 # control is crucial - bugs must be traceable. We will be happy to
28 # consider code for inclusion in the official distribution, but
29 # derived work must not be called official GROMACS. Details are found
30 # in the README & COPYING files - if they are missing, get the
31 # official version at http://www.gromacs.org.
33 # To help us fund GROMACS development, we humbly ask that you cite
34 # the research papers on the package. Check out http://www.gromacs.org.
36 # This script is used by the GROMACS developers to build most of the
37 # files from which the nbnxn kernels are compiled. It is not called at
38 # CMake time, and users should never need to use it. It currently
39 # works for nbnxn kernel structure types 2xnn and 4xn. The generated
40 # files are versions of the *.pre files in this directory, customized
41 # for the kernel structure type and/or the detailed kernel type. These
44 # A single header file that declares all the kernel functions for
45 # this nbnxn kernel structure type, including the function that does
46 # the dispatch via the function pointer table.
48 # A single C kernel dispatcher file that defines the function that
49 # decides at run time which kernel to call.
51 # Many C kernel files, each defining a single kernel function. These
52 # functions can take a noticeable time to compile, and should tend
53 # to be in seperate files to take advantage of make-time
56 # This script should be run from the directory in which it is
57 # located. The generated files are located in ../simd_<type>. There
58 # are three other files in those locations that are not generated. These
61 # setup logic peculiar to the kernel structure type but common to
62 # all the kernels within that type, and
64 # the logic for the outer and inner loops of the kernels, as
65 # customized by numerous preprocessor defines to suit the hardware
68 # Note that while functions for both nbnxn kernel structures are
69 # compiled and built into an mdrun executable, because that executable
70 # is not portable, only the functions for the useful nbnxn kernel
71 # structure for the hardware selected at CMake time contain real
72 # kernel logic. A run-time error occurs if an inappropriate kernel
73 # dispatcher function is called (but that is normally impossible).
78 import collections # Requires Python 2.7
79 sys.path.append('../../../../../admin')
80 from copyright import create_copyright_header
82 FileHeader = create_copyright_header('2012,2013,2014')
84 * Note: this file was generated by the Verlet kernel generator for
90 def read_kernel_template(filename):
91 with open(filename, "r") as TemplateFile:
92 TemplateText = TemplateFile.read()
93 copyright_re = r'/\*\n \* This file is part of the GROMACS molecular simulation package\.\n( \*.*\n)* \*/\n'
94 match = re.match(copyright_re, TemplateText)
96 TemplateText = TemplateText[match.end():]
99 # The dict order must match the order of an enumeration in
100 # nbnxn_kernel_simd_template.c.pre
101 ElectrostaticsDict = collections.OrderedDict()
102 ElectrostaticsDict['ElecRF'] = { 'define' : '#define CALC_COUL_RF' }
103 ElectrostaticsDict['ElecQSTab'] = { 'define' : '#define CALC_COUL_TAB' }
104 ElectrostaticsDict['ElecQSTabTwinCut'] = { 'define' : '#define CALC_COUL_TAB\n#define VDW_CUTOFF_CHECK /* Use twin-range cut-off */' }
105 ElectrostaticsDict['ElecEw'] = { 'define' : '#define CALC_COUL_EWALD' }
106 ElectrostaticsDict['ElecEwTwinCut'] = { 'define' : '#define CALC_COUL_EWALD\n#define VDW_CUTOFF_CHECK /* Use twin-range cut-off */' }
108 # The dict order must match the order of a C enumeration.
109 VdwTreatmentDict = collections.OrderedDict()
110 VdwTreatmentDict['VdwLJCombGeom'] = { 'define' : '#define LJ_CUT\n#define LJ_COMB_GEOM' }
111 VdwTreatmentDict['VdwLJCombLB'] = { 'define' : '#define LJ_CUT\n#define LJ_COMB_LB' }
112 VdwTreatmentDict['VdwLJ'] = { 'define' : '#define LJ_CUT\n/* Use full LJ combination matrix */' }
113 VdwTreatmentDict['VdwLJFSw'] = { 'define' : '#define LJ_FORCE_SWITCH\n/* Use full LJ combination matrix */' }
114 VdwTreatmentDict['VdwLJPSw'] = { 'define' : '#define LJ_POT_SWITCH\n/* Use full LJ combination matrix */' }
115 VdwTreatmentDict['VdwLJEwCombGeom'] = { 'define' : '#define LJ_CUT\n#define LJ_EWALD_GEOM\n/* Use full LJ combination matrix + geometric rule for the grid correction */' }
117 # This is OK as an unordered dict
118 EnergiesComputationDict = {
120 'function type' : 'nbk_func_noener',
121 'define' : '/* Will not calculate energies */',
124 'function type' : 'nbk_func_ener',
125 'define' : '#define CALC_ENERGIES',
128 'function type' : 'nbk_func_ener',
129 'define' : '#define CALC_ENERGIES\n#define ENERGY_GROUPS',
133 # This is OK as an unordered dict
134 VerletKernelTypeDict = {
136 'Define' : 'GMX_NBNXN_SIMD_2XNN',
137 'WidthSetup' : '/* Include the full-width SIMD macros */\n',
138 'WidthCheck' : ('#if !(GMX_SIMD_REAL_WIDTH == 8 || GMX_SIMD_REAL_WIDTH == 16)\n' \
139 '#error "unsupported SIMD width"\n' \
144 'Define' : 'GMX_NBNXN_SIMD_4XN',
146 'WidthCheck' : ('#if !(GMX_SIMD_REAL_WIDTH == 2 || GMX_SIMD_REAL_WIDTH == 4 || GMX_SIMD_REAL_WIDTH == 8)\n' \
147 '#error "unsupported SIMD width"\n' \
153 KernelDispatcherTemplate = read_kernel_template("nbnxn_kernel_simd_template.c.pre")
154 KernelsHeaderTemplate = read_kernel_template("nbnxn_kernel_simd_template.h.pre")
156 # For each Verlet kernel type, write three kinds of files:
157 # a header file defining the functions for all the kernels,
158 # a code file containing the kernel function lookup table and
159 # the kernel dispatcher function
160 # for each kernel, a file defining the single C function for that kernel
161 for type in VerletKernelTypeDict:
162 DirName = "../simd_{0}".format(type)
163 KernelNamePrefix = 'nbnxn_kernel'
164 KernelsName = "{0}_simd_{1}".format(KernelNamePrefix,type)
165 KernelsHeaderFileName = "{0}.h".format(KernelsName,type)
166 KernelFunctionLookupTable = {}
167 KernelDeclarations = ''
168 KernelTemplate = read_kernel_template("{0}_kernel.c.pre".format(KernelsName))
170 # Loop over all kernels
171 for ener in EnergiesComputationDict:
172 KernelFunctionLookupTable[ener] = '{\n'
173 for elec in ElectrostaticsDict:
174 KernelFunctionLookupTable[ener] += ' {\n'
175 for ljtreat in VdwTreatmentDict:
176 KernelName = ('{0}_{1}_{2}_{3}_{4}'
177 .format(KernelNamePrefix,elec,ljtreat,ener,type))
179 # Declare the kernel function
180 KernelDeclarations += ('{1:21} {0};\n'
182 EnergiesComputationDict[ener]['function type']))
184 # Write the file with the kernel definition
185 with open('{0}/{1}.c'.format(DirName,KernelName), 'w') as kernelfp:
186 kernelfp.write(FileHeader.format(type))
187 kernelfp.write(KernelTemplate
188 .format(VerletKernelTypeDict[type]['Define'],
189 ElectrostaticsDict[elec]['define'],
190 VdwTreatmentDict[ljtreat]['define'],
191 EnergiesComputationDict[ener]['define'],
192 KernelsHeaderFileName,
194 " " * (len(KernelName) + 1),
195 VerletKernelTypeDict[type]['UnrollSize'],
199 # Enter the kernel function in the lookup table
200 KernelFunctionLookupTable[ener] += ' {0},\n'.format(KernelName)
202 KernelFunctionLookupTable[ener] += ' },\n'
203 KernelFunctionLookupTable[ener] += '};\n'
204 KernelDeclarations += '\n'
206 # Write the header file that declares all the kernel
207 # functions for this type
208 with open('{0}/{1}'.format(DirName,KernelsHeaderFileName),'w') as fp:
209 fp.write(FileHeader.format(type))
210 fp.write(KernelsHeaderTemplate
212 " " * (len(KernelsName) + 1),
215 # Write the file defining the kernel dispatcher
216 # function for this type
217 with open('{0}/{1}'.format(DirName,"{0}.c".format(KernelsName,type)),'w') as fp:
218 fp.write(FileHeader.format(type))
219 fp.write(KernelDispatcherTemplate
220 .format(VerletKernelTypeDict[type]['Define'],
221 VerletKernelTypeDict[type]['WidthSetup'],
222 VerletKernelTypeDict[type]['WidthCheck'],
223 VerletKernelTypeDict[type]['UnrollSize'],
224 KernelsHeaderFileName,
226 ' ' * (len(KernelsName)+1),
227 KernelFunctionLookupTable['F'],
228 KernelFunctionLookupTable['VF'],
229 KernelFunctionLookupTable['VgrpF'],