Converted 4xn kernel to C++
[gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_file_generator / make_verlet_simd_kernel_files.py
blob0b541e752d20174f60e5fcac29b69464d81d1038
1 #!/usr/bin/python
3 # This file is part of the GROMACS molecular simulation package.
5 # Copyright (c) 2013,2014,2015, 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
42 # are:
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
54 # parallelism.
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
59 # contain:
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
66 # and kernel type.
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).
75 import re
76 import sys
77 import os
78 os.chdir(os.path.dirname(os.path.abspath(__file__)))
79 import collections # Requires Python 2.7
80 sys.path.append('../../../../../admin')
81 from copyright import create_copyright_header
83 FileHeader = create_copyright_header('2012,2013,2014,2015')
84 FileHeader += """/*
85 * Note: this file was generated by the Verlet kernel generator for
86 * kernel type {0}.
89 """
91 def read_kernel_template(filename):
92 with open(filename, "r") as TemplateFile:
93 TemplateText = TemplateFile.read()
94 copyright_re = r'/\*\n \* This file is part of the GROMACS molecular simulation package\.\n( \*.*\n)* \*/\n'
95 match = re.match(copyright_re, TemplateText)
96 if match:
97 TemplateText = TemplateText[match.end():]
98 return TemplateText
100 # The dict order must match the order of an enumeration in
101 # nbnxn_kernel_simd_template.c.pre
102 ElectrostaticsDict = collections.OrderedDict()
103 ElectrostaticsDict['ElecRF'] = { 'define' : '#define CALC_COUL_RF' }
104 ElectrostaticsDict['ElecQSTab'] = { 'define' : '#define CALC_COUL_TAB' }
105 ElectrostaticsDict['ElecQSTabTwinCut'] = { 'define' : '#define CALC_COUL_TAB\n#define VDW_CUTOFF_CHECK /* Use twin-range cut-off */' }
106 ElectrostaticsDict['ElecEw'] = { 'define' : '#define CALC_COUL_EWALD' }
107 ElectrostaticsDict['ElecEwTwinCut'] = { 'define' : '#define CALC_COUL_EWALD\n#define VDW_CUTOFF_CHECK /* Use twin-range cut-off */' }
109 # The dict order must match the order of a C enumeration.
110 VdwTreatmentDict = collections.OrderedDict()
111 VdwTreatmentDict['VdwLJCombGeom'] = { 'define' : '#define LJ_CUT\n#define LJ_COMB_GEOM' }
112 VdwTreatmentDict['VdwLJCombLB'] = { 'define' : '#define LJ_CUT\n#define LJ_COMB_LB' }
113 VdwTreatmentDict['VdwLJ'] = { 'define' : '#define LJ_CUT\n/* Use full LJ combination matrix */' }
114 VdwTreatmentDict['VdwLJFSw'] = { 'define' : '#define LJ_FORCE_SWITCH\n/* Use full LJ combination matrix */' }
115 VdwTreatmentDict['VdwLJPSw'] = { 'define' : '#define LJ_POT_SWITCH\n/* Use full LJ combination matrix */' }
116 VdwTreatmentDict['VdwLJEwCombGeom'] = { 'define' : '#define LJ_CUT\n#define LJ_EWALD_GEOM\n/* Use full LJ combination matrix + geometric rule for the grid correction */' }
118 # This is OK as an unordered dict
119 EnergiesComputationDict = {
120 'F' : {
121 'function type' : 'nbk_func_noener',
122 'define' : '/* Will not calculate energies */',
124 'VF' : {
125 'function type' : 'nbk_func_ener',
126 'define' : '#define CALC_ENERGIES',
128 'VgrpF' : {
129 'function type' : 'nbk_func_ener',
130 'define' : '#define CALC_ENERGIES\n#define ENERGY_GROUPS',
134 # This is OK as an unordered dict
135 VerletKernelTypeDict = {
136 '2xnn' : {
137 'Define' : 'GMX_NBNXN_SIMD_2XNN',
138 'WidthSetup' : '/* Include the full-width SIMD macros */\n',
139 'WidthCheck' : ('#if !(GMX_SIMD_REAL_WIDTH == 8 || GMX_SIMD_REAL_WIDTH == 16)\n' \
140 '#error "unsupported SIMD width"\n' \
141 '#endif\n'),
142 'UnrollSize' : 2,
144 '4xn' : {
145 'Define' : 'GMX_NBNXN_SIMD_4XN',
146 'WidthSetup' : (''),
147 'WidthCheck' : ('#if !(GMX_SIMD_REAL_WIDTH == 2 || GMX_SIMD_REAL_WIDTH == 4 || GMX_SIMD_REAL_WIDTH == 8)\n' \
148 '#error "unsupported SIMD width"\n' \
149 '#endif\n'),
150 'UnrollSize' : 1,
154 KernelDispatcherTemplate = read_kernel_template("nbnxn_kernel_simd_template.cpp.pre")
155 KernelsHeaderTemplate = read_kernel_template("nbnxn_kernel_simd_template.h.pre")
157 # For each Verlet kernel type, write three kinds of files:
158 # a header file defining the functions for all the kernels,
159 # a code file containing the kernel function lookup table and
160 # the kernel dispatcher function
161 # for each kernel, a file defining the single C function for that kernel
162 for type in VerletKernelTypeDict:
163 DirName = "../simd_{0}".format(type)
164 KernelNamePrefix = 'nbnxn_kernel'
165 KernelsName = "{0}_simd_{1}".format(KernelNamePrefix,type)
166 KernelsHeaderFileName = "{0}.h".format(KernelsName,type)
167 KernelsHeaderPathName = "gromacs/mdlib/nbnxn_kernels/simd_{0}/{1}".format(type,KernelsHeaderFileName)
168 KernelFunctionLookupTable = {}
169 KernelDeclarations = ''
170 KernelTemplate = read_kernel_template("{0}_kernel.cpp.pre".format(KernelsName))
172 # Loop over all kernels
173 for ener in EnergiesComputationDict:
174 KernelFunctionLookupTable[ener] = '{\n'
175 for elec in ElectrostaticsDict:
176 KernelFunctionLookupTable[ener] += ' {\n'
177 for ljtreat in VdwTreatmentDict:
178 KernelName = ('{0}_{1}_{2}_{3}_{4}'
179 .format(KernelNamePrefix,elec,ljtreat,ener,type))
181 # Declare the kernel function
182 KernelDeclarations += ('{1:21} {0};\n'
183 .format(KernelName,
184 EnergiesComputationDict[ener]['function type']))
186 # Write the file with the kernel definition
187 with open('{0}/{1}.cpp'.format(DirName,KernelName), 'w') as kernelfp:
188 kernelfp.write(FileHeader.format(type))
189 kernelfp.write(KernelTemplate
190 .format(VerletKernelTypeDict[type]['Define'],
191 ElectrostaticsDict[elec]['define'],
192 VdwTreatmentDict[ljtreat]['define'],
193 EnergiesComputationDict[ener]['define'],
194 KernelsHeaderPathName,
195 KernelName,
196 " " * (len(KernelName) + 1),
197 VerletKernelTypeDict[type]['UnrollSize'],
201 # Enter the kernel function in the lookup table
202 KernelFunctionLookupTable[ener] += ' {0},\n'.format(KernelName)
204 KernelFunctionLookupTable[ener] += ' },\n'
205 KernelFunctionLookupTable[ener] += '};\n'
206 KernelDeclarations += '\n'
208 # Write the header file that declares all the kernel
209 # functions for this type
210 with open('{0}/{1}'.format(DirName,KernelsHeaderFileName),'w') as fp:
211 fp.write(FileHeader.format(type))
212 fp.write(KernelsHeaderTemplate
213 .format(KernelsName,
214 " " * (len(KernelsName) + 1),
215 KernelDeclarations))
217 # Write the file defining the kernel dispatcher
218 # function for this type
219 with open('{0}/{1}'.format(DirName,"{0}.cpp".format(KernelsName,type)),'w') as fp:
220 fp.write(FileHeader.format(type))
221 fp.write(KernelDispatcherTemplate
222 .format(VerletKernelTypeDict[type]['Define'],
223 VerletKernelTypeDict[type]['WidthSetup'],
224 VerletKernelTypeDict[type]['WidthCheck'],
225 VerletKernelTypeDict[type]['UnrollSize'],
226 KernelsHeaderFileName,
227 KernelsName,
228 ' ' * (len(KernelsName)+1),
229 KernelFunctionLookupTable['F'],
230 KernelFunctionLookupTable['VF'],
231 KernelFunctionLookupTable['VgrpF'],
235 sys.exit()