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