0f435c2b899e82a010ea7b2cb20d117d365a1f7f
[alexxy/gromacs.git] / src / gromacs / nbnxm / kernel_file_generator / make_verlet_simd_kernel_files.py
1 #!/usr/bin/env python3
2 #
3 # This file is part of the GROMACS molecular simulation package.
4 #
5 # Copyright (c) 2013,2014,2015,2017,2018,2019, 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.
9 #
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.
14 #
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.
19 #
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.
24 #
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.
32 #
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.
35
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 2xMM and 4xM. 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:
43 #
44 #   A single header file that declares all the kernel functions for
45 #   this nbnxn kernel structure type and a function pointer table.
46 #
47 #   Many C kernel files, each defining a single kernel function. These
48 #   functions can take a noticeable time to compile, and should tend
49 #   to be in seperate files to take advantage of make-time
50 #   parallelism.
51 #
52 # This script should be run from the directory in which it is
53 # located. The generated files are located in ../kernels_simd_<type>.
54 # There are three other files in those locations that are not generated.
55 # These contain:
56 #
57 #   setup logic peculiar to the kernel structure type but common to
58 #   all the kernels within that type, and
59 #
60 #   the logic for the outer and inner loops of the kernels, as
61 #   customized by numerous preprocessor defines to suit the hardware
62 #   and kernel type.
63 #
64 # Note that while functions for both nbnxn kernel structures are
65 # compiled and built into an mdrun executable, because that executable
66 # is not portable, only the functions for the useful nbnxn kernel
67 # structure for the hardware selected at CMake time contain real
68 # kernel logic. A run-time error occurs if an inappropriate kernel
69 # dispatcher function is called (but that is normally impossible).
70
71 import re
72 import sys
73 import os
74 os.chdir(os.path.dirname(os.path.abspath(__file__)))
75 import collections # Requires Python 2.7
76 sys.path.append('../../../../admin')
77 from copyright import create_copyright_header
78
79 FileHeader = create_copyright_header('2012,2013,2014,2015,2019')
80 FileHeader += """/*
81  * Note: this file was generated by the Verlet kernel generator for
82  * kernel type {0}.
83  */
84
85 """
86
87 def read_kernel_template(filename):
88     with open(filename, "r") as TemplateFile:
89         TemplateText = TemplateFile.read()
90     copyright_re = r'/\*\n \* This file is part of the GROMACS molecular simulation package\.\n( \*.*\n)* \*/\n'
91     match = re.match(copyright_re, TemplateText)
92     if match:
93         TemplateText = TemplateText[match.end():]
94     return TemplateText
95
96 # The dict order must match the order of an enumeration in
97 # kernel_simd_template.c.pre
98 ElectrostaticsDict = collections.OrderedDict()
99 ElectrostaticsDict['ElecRF'] = { 'define' : '#define CALC_COUL_RF' }
100 ElectrostaticsDict['ElecQSTab'] = { 'define' : '#define CALC_COUL_TAB' }
101 ElectrostaticsDict['ElecQSTabTwinCut'] = { 'define' : '#define CALC_COUL_TAB\n#define VDW_CUTOFF_CHECK /* Use twin-range cut-off */' }
102 ElectrostaticsDict['ElecEw'] = { 'define' : '#define CALC_COUL_EWALD' }
103 ElectrostaticsDict['ElecEwTwinCut'] = { 'define' : '#define CALC_COUL_EWALD\n#define VDW_CUTOFF_CHECK /* Use twin-range cut-off */' }
104  
105 # The dict order must match the order of a C enumeration.
106 VdwTreatmentDict = collections.OrderedDict()
107 VdwTreatmentDict['VdwLJCombGeom'] = { 'define' : '#define LJ_CUT\n#define LJ_COMB_GEOM' }
108 VdwTreatmentDict['VdwLJCombLB'] = { 'define' : '#define LJ_CUT\n#define LJ_COMB_LB' }
109 VdwTreatmentDict['VdwLJ'] = { 'define' : '#define LJ_CUT\n/* Use full LJ combination matrix */' }
110 VdwTreatmentDict['VdwLJFSw'] = { 'define' : '#define LJ_FORCE_SWITCH\n/* Use full LJ combination matrix */' }
111 VdwTreatmentDict['VdwLJPSw'] = { 'define' : '#define LJ_POT_SWITCH\n/* Use full LJ combination matrix */' }
112 VdwTreatmentDict['VdwLJEwCombGeom'] = { 'define' : '#define LJ_CUT\n#define LJ_EWALD_GEOM\n/* Use full LJ combination matrix + geometric rule for the grid correction */' }
113
114 # This is OK as an unordered dict
115 EnergiesComputationDict = {
116     'F'  : {
117         'function type' : 'nbk_func_noener',
118         'define' : '/* Will not calculate energies */',
119     },
120     'VF'    : {
121         'function type' : 'nbk_func_ener',
122         'define' : '#define CALC_ENERGIES',
123     },
124     'VgrpF' : {
125         'function type' : 'nbk_func_ener',
126         'define' : '#define CALC_ENERGIES\n#define ENERGY_GROUPS',
127     },
128 }
129
130 # This is OK as an unordered dict
131 VerletKernelTypeDict = {
132     '2xmm' : {
133         'Define' : 'GMX_NBNXN_SIMD_2XNN',
134         'WidthSetup' : '/* Include the full-width SIMD macros */\n',
135         'WidthCheck' : ('#if !(GMX_SIMD_REAL_WIDTH == 8 || GMX_SIMD_REAL_WIDTH == 16)\n' \
136                         '#error "unsupported SIMD width"\n' \
137                         '#endif\n'),
138         'UnrollSize' : 2,
139     },
140     '4xm' : {
141         'Define' : 'GMX_NBNXN_SIMD_4XN',
142         'WidthSetup' : (''),
143         'WidthCheck' : ('#if !(GMX_SIMD_REAL_WIDTH == 2 || GMX_SIMD_REAL_WIDTH == 4 || GMX_SIMD_REAL_WIDTH == 8)\n' \
144                         '#error "unsupported SIMD width"\n' \
145                         '#endif\n'),
146         'UnrollSize' : 1,
147     },
148 }
149
150 KernelsHeaderTemplate = read_kernel_template("kernel_simd_template.h.pre")
151
152 # For each Verlet kernel type, write two kinds of files:
153 #   a header file defining the functions for all the kernels and
154 #     the kernel function lookup table
155 #   for each kernel, a file defining the single C function for that kernel
156 for type in VerletKernelTypeDict:
157     DirName = "../kernels_simd_{0}".format(type)
158     KernelNamePrefix = 'nbnxm_kernel'
159     KernelFileNamePrefix = 'kernel'
160     KernelsName = "{0}_simd_{1}".format(KernelNamePrefix,type)
161     KernelsFileName = "{0}_simd_{1}".format(KernelFileNamePrefix,type)
162     KernelsHeaderFileName = "kernels.h"
163     KernelsHeaderPathName = "{1}".format(type,KernelsHeaderFileName)
164     KernelFunctionLookupTable = {}
165     KernelDeclarations = ''
166     KernelTemplate = read_kernel_template("{0}_kernel.cpp.pre".format(KernelsFileName))
167
168     # Loop over all kernels
169     for ener in EnergiesComputationDict:
170         KernelFunctionLookupTable[ener] = '{\n'
171         for elec in ElectrostaticsDict:
172             KernelFunctionLookupTable[ener] += '    {\n'
173             for ljtreat in VdwTreatmentDict:
174                 KernelName = ('{0}_{1}_{2}_{3}_{4}'
175                               .format(KernelNamePrefix,elec,ljtreat,ener,type))
176                 KernelFileName = ('{0}_{1}_{2}_{3}'
177                                   .format(KernelFileNamePrefix,elec,ljtreat,ener,type))
178
179                 # Declare the kernel function
180                 KernelDeclarations += ('{1:21} {0};\n'
181                                        .format(KernelName,
182                                                EnergiesComputationDict[ener]['function type']))
183
184                 # Write the file with the kernel definition
185                 with open('{0}/{1}.cpp'.format(DirName,KernelFileName), '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                                            KernelsHeaderPathName,
193                                            KernelName,
194                                            " " * (len(KernelName) + 1),
195                                            VerletKernelTypeDict[type]['UnrollSize'],
196                                        )
197                                )
198
199                 # Enter the kernel function in the lookup table
200                 KernelFunctionLookupTable[ener] += '        {0},\n'.format(KernelName)
201
202             KernelFunctionLookupTable[ener] += '    },\n'
203         KernelFunctionLookupTable[ener] += '};\n'
204         KernelDeclarations += '\n'
205
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
211                  .format(KernelDeclarations,
212                          type,
213                          KernelFunctionLookupTable['F'],
214                          KernelFunctionLookupTable['VF'],
215                          KernelFunctionLookupTable['VgrpF'])
216              )
217
218 sys.exit()