Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / make_nb_kernel_avx_128_fma_single.py
1 #!/usr/bin/python
2 #
3 # This file is part of the GROMACS molecular simulation package.
4 #
5 # Copyright (c) 2012,2013, by the GROMACS development team, led by
6 # David van der Spoel, Berk Hess, Erik Lindahl, and including many
7 # others, as listed in the AUTHORS file in the top-level source
8 # 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 import sys
37 import os
38 sys.path.append ( "../preprocessor" )
39 from gmxpreprocess import gmxpreprocess
40
41 # "The happiest programs are programs that write other programs."
42 #
43 #
44 # This script controls the generation of Gromacs nonbonded kernels.
45 #
46 # We no longer generate kernels on-the-fly, so this file is not run
47 # during a Gromacs compile - only when we need to update the kernels (=rarely).
48 #
49 # To maximize performance, each combination of interactions in Gromacs
50 # has a separate nonbonded kernel without conditionals in the code.
51 # To avoid writing hundreds of different routines for each architecture,
52 # we instead use a custom preprocessor so we can encode the conditionals
53 # and expand for-loops (e.g, for water-water interactions)
54 # from a general kernel template. While that file will contain quite a
55 # few preprocessor directives, it is still an order of magnitude easier
56 # to maintain than ~200 different kernels (not to mention it avoids bugs).
57 #
58 # To actually generate the kernels, this program iteratively calls the
59 # preprocessor with different define settings corresponding to all
60 # combinations of coulomb/van-der-Waals/geometry options.
61 #
62 # A main goal in the design was to make this new generator _general_. For
63 # this reason we have used a lot of different fields to identify a particular
64 # kernel and interaction. Basically, each kernel will have a name like
65 #
66 # nbkernel_ElecXX_VdwYY_GeomZZ_VF_QQ()
67 #
68 # Where XX/YY/ZZ/VF are strings to identify what the kernel computes.
69 #
70 # Elec/Vdw describe the type of interaction for electrostatics and van der Waals.
71 # The geometry settings correspond e.g. to water-water or water-particle kernels,
72 # and finally the VF setting is V,F,or VF depending on whether we calculate
73 # only the potential, only the force, or both of them. The final string (QQ)
74 # is the architecture/language/optimization of the kernel.
75 #
76 Arch       = 'avx_128_fma_single'
77
78 # Explanation of the 'properties':
79 #
80 # It is cheap to compute r^2, and the kernels require various other functions of r for
81 # different kinds of interaction. Depending on the needs of the kernel and the available
82 # processor instructions, this will be done in different ways.
83 #
84 # 'rinv' means we need 1/r, which is calculated as 1/sqrt(r^2).
85 # 'rinvsq' means we need 1/(r*r). This is calculated as rinv*rinv if we already did rinv, otherwise 1/r^2.
86 # 'r' is similarly calculated as r^2*rinv when needed
87 # 'table' means the interaction is tabulated, in which case we will calculate a table index before the interaction
88 # 'shift' means the interaction will be modified by a constant to make it zero at the cutoff.
89 # 'cutoff' means the interaction is set to 0.0 outside the cutoff
90 #
91
92 FileHeader = \
93 '/*\n' \
94 ' * This file is part of the GROMACS molecular simulation package.\n' \
95 ' *\n' \
96 ' * Copyright (c) 2012,2013, by the GROMACS development team, led by\n' \
97 ' * David van der Spoel, Berk Hess, Erik Lindahl, and including many\n' \
98 ' * others, as listed in the AUTHORS file in the top-level source\n' \
99 ' * directory and at http://www.gromacs.org.\n' \
100 ' *\n' \
101 ' * GROMACS is free software; you can redistribute it and/or\n' \
102 ' * modify it under the terms of the GNU Lesser General Public License\n' \
103 ' * as published by the Free Software Foundation; either version 2.1\n' \
104 ' * of the License, or (at your option) any later version.\n' \
105 ' *\n' \
106 ' * GROMACS is distributed in the hope that it will be useful,\n' \
107 ' * but WITHOUT ANY WARRANTY; without even the implied warranty of\n' \
108 ' * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU\n' \
109 ' * Lesser General Public License for more details.\n' \
110 ' *\n' \
111 ' * You should have received a copy of the GNU Lesser General Public\n' \
112 ' * License along with GROMACS; if not, see\n' \
113 ' * http://www.gnu.org/licenses, or write to the Free Software Foundation,\n' \
114 ' * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.\n' \
115 ' *\n' \
116 ' * If you want to redistribute modifications to GROMACS, please\n' \
117 ' * consider that scientific software is very special. Version\n' \
118 ' * control is crucial - bugs must be traceable. We will be happy to\n' \
119 ' * consider code for inclusion in the official distribution, but\n' \
120 ' * derived work must not be called official GROMACS. Details are found\n' \
121 ' * in the README & COPYING files - if they are missing, get the\n' \
122 ' * official version at http://www.gromacs.org.\n' \
123 ' *\n' \
124 ' * To help us fund GROMACS development, we humbly ask that you cite\n' \
125 ' * the research papers on the package. Check out http://www.gromacs.org.\n' \
126 ' */\n' \
127 '/*\n' \
128 ' * Note: this file was generated by the GROMACS '+Arch+' kernel generator.\n' \
129 ' */\n'
130
131 ###############################################
132 # ELECTROSTATICS
133 # Interactions and flags for them
134 ###############################################
135 ElectrostaticsList = {
136     'None'                    : [],
137     'Coulomb'                 : ['rinv','rinvsq'],
138     'ReactionField'           : ['rinv','rinvsq'],
139     'GeneralizedBorn'         : ['rinv','r'],
140     'CubicSplineTable'        : ['rinv','r','table'],
141     'Ewald'                   : ['rinv','rinvsq','r'],
142 }
143
144
145 ###############################################
146 # VAN DER WAALS
147 # Interactions and flags for them
148 ###############################################
149 VdwList = {
150     'None'                    : [],
151     'LennardJones'            : ['rinvsq'],
152 #    'Buckingham'              : ['rinv','rinvsq','r'], # Disabled for sse4.1 to reduce number of kernels and simply the template
153     'CubicSplineTable'        : ['rinv','r','table'],
154 }
155
156
157 ###############################################
158 # MODIFIERS
159 # Different ways to adjust/modify interactions to conserve energy
160 ###############################################
161 ModifierList = {
162     'None'                    : [],
163     'ExactCutoff'             : ['exactcutoff'],        # Zero the interaction outside the cutoff, used for reaction-field-zero
164     'PotentialShift'          : ['shift','exactcutoff'],
165     'PotentialSwitch'         : ['rinv','r','switch','exactcutoff']
166 }
167
168
169 ###############################################
170 # GEOMETRY COMBINATIONS
171 ###############################################
172 GeometryNameList = [
173     [ 'Particle' , 'Particle' ],
174     [ 'Water3'   , 'Particle' ],
175     [ 'Water3'   , 'Water3'   ],
176     [ 'Water4'   , 'Particle' ],
177     [ 'Water4'   , 'Water4'   ]
178 ]
179
180
181 ###############################################
182 # POTENTIAL / FORCE
183 ###############################################
184 VFList = [
185     'PotentialAndForce',
186 # 'Potential',   # Not used yet
187     'Force'
188 ]
189
190
191 ###############################################
192 # GEOMETRY PROPERTIES
193 ###############################################
194 # Dictionaries with lists telling which interactions are present
195 # 1,2,3 means particles 1,2,3 (but not 0) have electrostatics!
196 GeometryElectrostatics = {
197     'Particle'  : [ 0 ],
198     'Particle2' : [ 0 , 1 ],
199     'Particle3' : [ 0 , 1 , 2 ],
200     'Particle4' : [ 0 , 1 , 2 , 3 ],
201     'Water3'    : [ 0 , 1 , 2 ],
202     'Water4'    : [ 1 , 2 , 3 ]
203 }
204
205 GeometryVdw = {
206     'Particle'  : [ 0 ],
207     'Particle2' : [ 0 , 1 ],
208     'Particle3' : [ 0 , 1 , 2 ],
209     'Particle4' : [ 0 , 1 , 2 , 3 ],
210     'Water3'    : [ 0 ],
211     'Water4'    : [ 0 ]
212 }
213
214
215
216
217 # Dictionary to abbreviate all strings (mixed from all the lists)
218 Abbreviation = {
219     'None'                    : 'None',
220     'Coulomb'                 : 'Coul',
221     'Ewald'                   : 'Ew',
222     'ReactionField'           : 'RF',
223     'GeneralizedBorn'         : 'GB',
224     'CubicSplineTable'        : 'CSTab',
225     'LennardJones'            : 'LJ',
226     'Buckingham'              : 'Bham',
227     'PotentialShift'          : 'Sh',
228     'PotentialSwitch'         : 'Sw',
229     'ExactCutoff'             : 'Cut',
230     'PotentialAndForce'       : 'VF',
231     'Potential'               : 'V',
232     'Force'                   : 'F',
233     'Water3'                  : 'W3',
234     'Water4'                  : 'W4',
235     'Particle'                : 'P1',
236     'Particle2'               : 'P2',
237     'Particle3'               : 'P3',
238     'Particle4'               : 'P4'
239 }
240
241
242 ###############################################
243 # Functions
244 ###############################################
245
246 # Return a string with the kernel name from current settings
247 def MakeKernelFileName(KernelElec,KernelElecMod,KernelVdw,KernelVdwMod,KernelGeom):
248     ElecStr = 'Elec' + Abbreviation[KernelElec]
249     if(KernelElecMod!='None'):
250         ElecStr = ElecStr + Abbreviation[KernelElecMod]
251     VdwStr  = 'Vdw'  + Abbreviation[KernelVdw]
252     if(KernelVdwMod!='None'):
253         VdwStr = VdwStr + Abbreviation[KernelVdwMod]
254     GeomStr = 'Geom' + Abbreviation[KernelGeom[0]] + Abbreviation[KernelGeom[1]]
255     return 'nb_kernel_' + ElecStr + '_' + VdwStr + '_' + GeomStr + '_' + Arch
256
257 def MakeKernelName(KernelElec,KernelElecMod,KernelVdw,KernelVdwMod,KernelGeom,KernelVF):
258     ElecStr = 'Elec' + Abbreviation[KernelElec]
259     if(KernelElecMod!='None'):
260         ElecStr = ElecStr + Abbreviation[KernelElecMod]
261     VdwStr  = 'Vdw'  + Abbreviation[KernelVdw]
262     if(KernelVdwMod!='None'):
263         VdwStr = VdwStr + Abbreviation[KernelVdwMod]
264     GeomStr = 'Geom' + Abbreviation[KernelGeom[0]] + Abbreviation[KernelGeom[1]]
265     VFStr   = Abbreviation[KernelVF]
266     return 'nb_kernel_' + ElecStr + '_' + VdwStr + '_' + GeomStr + '_' + VFStr + '_' + Arch
267
268 # Return a string with a declaration to use for the kernel;
269 # this will be a sequence of string combinations as well as the actual function name
270 # Dont worry about field widths - that is just pretty-printing for the header!
271 def MakeKernelDecl(KernelName,KernelElec,KernelElecMod,KernelVdw,KernelVdwMod,KernelGeom,KernelOther,KernelVF):
272     KernelStr   = '\"'+KernelName+'\"'
273     ArchStr     = '\"'+Arch+'\"'
274     ElecStr     = '\"'+KernelElec+'\"'
275     ElecModStr  = '\"'+KernelElecMod+'\"'
276     VdwStr      = '\"'+KernelVdw+'\"'
277     VdwModStr   = '\"'+KernelVdwMod+'\"'
278     GeomStr     = '\"'+KernelGeom[0]+KernelGeom[1]+'\"'
279     OtherStr    = '\"'+KernelOther+'\"'
280     VFStr       = '\"'+KernelVF+'\"'
281
282     ThisSpec = ArchStr+', '+ElecStr+', '+ElecModStr+', '+VdwStr+', '+VdwModStr+', '+GeomStr+', '+OtherStr+', '+VFStr
283     ThisDecl = '    { '+KernelName+', '+KernelStr+', '+ThisSpec+' }'
284     return ThisDecl
285
286
287 # Returns 1 if this kernel should be created, 0 if we should skip it
288 # This routine is not critical - it is not the end of the world if we create more kernels,
289 # but since the number is pretty large we save both space and compile-time by reducing it a bit.
290 def KeepKernel(KernelElec,KernelElecMod,KernelVdw,KernelVdwMod,KernelGeom,KernelVF):
291
292     # No need for kernels without interactions
293     if(KernelElec=='None' and KernelVdw=='None'):
294         return 0
295
296     # No need for modifiers without interactions
297     if((KernelElec=='None' and KernelElecMod!='None') or (KernelVdw=='None' and KernelVdwMod!='None')):
298         return 0
299
300     # No need for LJ-only water optimization, or water optimization with implicit solvent.
301     if('Water' in KernelGeom[0] and (KernelElec=='None' or 'GeneralizedBorn' in KernelElec)):
302         return 0
303
304     # Non-matching table settings are pointless
305     if( ('Table' in KernelElec) and ('Table' in KernelVdw) and KernelElec!=KernelVdw ):
306         return 0
307
308     # Try to reduce the number of different switch/shift options to get a reasonable number of kernels
309     # For electrostatics, reaction-field can use 'exactcutoff', and ewald can use switch or shift.
310     if(KernelElecMod=='ExactCutoff' and KernelElec!='ReactionField'):
311         return 0
312     if(KernelElecMod in ['PotentialShift','PotentialSwitch'] and KernelElec!='Ewald'):
313         return 0
314     # For Vdw, we support switch and shift for Lennard-Jones/Buckingham
315     if((KernelVdwMod=='ExactCutoff') or
316        (KernelVdwMod in ['PotentialShift','PotentialSwitch'] and KernelVdw not in ['LennardJones','Buckingham'])):
317         return 0
318
319     # Choose either switch or shift and don't mix them...
320     if((KernelElecMod=='PotentialShift' and KernelVdwMod=='PotentialSwitch') or
321        (KernelElecMod=='PotentialSwitch' and KernelVdwMod=='PotentialShift')):
322         return 0
323
324     # Don't use a Vdw kernel with a modifier if the electrostatics one does not have one
325     if(KernelElec!='None' and KernelElecMod=='None' and KernelVdwMod!='None'):
326         return 0
327
328     # Don't use an electrostatics kernel with a modifier if the vdw one does not have one,
329     # unless the electrostatics one is reaction-field with exact cutoff.
330     if(KernelVdw!='None' and KernelVdwMod=='None' and KernelElecMod!='None'):
331         if(KernelElec=='ReactionField' and KernelVdw!='CubicSplineTable'):
332             return 0
333         elif(KernelElec!='ReactionField'):
334             return 0
335
336     return 1
337
338
339
340 #
341 # The preprocessor will automatically expand the interactions for water and other
342 # geometries inside the kernel, but to get this right we need to setup a couple
343 # of defines - we do them in a separate routine to keep the main loop clean.
344 #
345 # While this routine might look a bit complex it is actually quite straightforward,
346 # and the best news is that you wont have to modify _anything_ for a new geometry
347 # as long as you correctly define its Electrostatics/Vdw geometry in the lists above!
348 #
349 def SetDefines(KernelElec,KernelElecMod,KernelVdw,KernelVdwMod,KernelGeom,KernelVF,defines):
350     # What is the _name_ for the i/j group geometry?
351     igeometry            = KernelGeom[0]
352     jgeometry            = KernelGeom[1]
353     # define so we can access it in the source when the preprocessor runs
354     defines['GEOMETRY_I'] = igeometry
355     defines['GEOMETRY_J'] = jgeometry
356
357     # For the i/j groups, extract a python list of which sites have electrostatics
358     # For SPC/TIP3p this will be [1,1,1], while TIP4p (no elec on first site) will be [0,1,1,1]
359     ielec                = GeometryElectrostatics[igeometry]
360     jelec                = GeometryElectrostatics[jgeometry]
361     # Zero out the corresponding lists in case we dont do Elec
362     if(KernelElec=='None'):
363         ielec = []
364         jelec = []
365
366     # Extract similar interaction lists for Vdw interactions (example for SPC: [1,0,0])
367     iVdw                 = GeometryVdw[igeometry]
368     jVdw                 = GeometryVdw[jgeometry]
369
370     # Zero out the corresponding lists in case we dont do Vdw
371     if(KernelVdw=='None'):
372         iVdw = []
373         jVdw = []
374
375     # iany[] and jany[] contains lists of the particles actually used (for interactions) in this kernel
376     iany = list(set(ielec+iVdw))  # convert to+from set to make elements unique
377     jany = list(set(jelec+jVdw))
378
379     defines['PARTICLES_ELEC_I'] = ielec
380     defines['PARTICLES_ELEC_J'] = jelec
381     defines['PARTICLES_VDW_I']  = iVdw
382     defines['PARTICLES_VDW_J']  = jVdw
383     defines['PARTICLES_I']      = iany
384     defines['PARTICLES_J']      = jany
385
386     # elecij,Vdwij are sets with pairs of particles for which the corresponding interaction is done
387     # (and anyij again corresponds to either electrostatics or Vdw)
388     elecij = []
389     Vdwij  = []
390     anyij  = []
391
392     for i in ielec:
393         for j in jelec:
394             elecij.append([i,j])
395
396     for i in iVdw:
397         for j in jVdw:
398             Vdwij.append([i,j])
399
400     for i in iany:
401         for j in jany:
402             if [i,j] in elecij or [i,j] in Vdwij:
403                 anyij.append([i,j])
404
405     defines['PAIRS_IJ']     = anyij
406
407     # Make an 2d list-of-distance-properties-to-calculate for i,j
408     ni = max(iany)+1
409     nj = max(jany)+1
410     # Each element properties[i][j] is an empty list
411     properties = [ [ [] for j in range(0,nj) ] for i in range (0,ni) ]
412     # Add properties to each set
413     for i in range(0,ni):
414         for j in range(0,nj):
415             if [i,j] in elecij:
416                 properties[i][j] = properties[i][j] + ['electrostatics'] + ElectrostaticsList[KernelElec] + ModifierList[KernelElecMod]
417             if [i,j] in Vdwij:
418                 properties[i][j] = properties[i][j] + ['vdw'] + VdwList[KernelVdw] + ModifierList[KernelVdwMod]
419             # Add rinv if we need r
420             if 'r' in properties[i][j]:
421                 properties[i][j] = properties[i][j] + ['rinv']
422             # Add rsq if we need rinv or rinsq
423             if 'rinv' in properties[i][j] or 'rinvsq' in properties[i][j]:
424                 properties[i][j] = properties[i][j] + ['rsq']
425
426     defines['INTERACTION_FLAGS']    = properties
427
428
429
430 def PrintStatistics(ratio):
431     ratio = 100.0*ratio
432     print '\rGenerating %s nonbonded kernels... %5.1f%%' % (Arch,ratio),
433     sys.stdout.flush()
434
435
436
437 defines = {}
438 kerneldecl = []
439
440 cnt     = 0.0
441 nelec   = len(ElectrostaticsList)
442 nVdw    = len(VdwList)
443 nmod    = len(ModifierList)
444 ngeom   = len(GeometryNameList)
445
446 ntot    = nelec*nmod*nVdw*nmod*ngeom
447
448 numKernels = 0
449
450 fpdecl = open('nb_kernel_' + Arch + '.c','w')
451 fpdecl.write( FileHeader )
452 fpdecl.write( '#ifndef nb_kernel_' + Arch + '_h\n' )
453 fpdecl.write( '#define nb_kernel_' + Arch + '_h\n\n' )
454 fpdecl.write( '#include "../nb_kernel.h"\n\n' )
455
456 for KernelElec in ElectrostaticsList:
457     defines['KERNEL_ELEC'] = KernelElec
458
459     for KernelElecMod in ModifierList:
460         defines['KERNEL_MOD_ELEC'] = KernelElecMod
461
462         for KernelVdw in VdwList:
463             defines['KERNEL_VDW'] = KernelVdw
464
465             for KernelVdwMod in ModifierList:
466                 defines['KERNEL_MOD_VDW'] = KernelVdwMod
467
468                 for KernelGeom in GeometryNameList:
469
470                     cnt += 1
471                     KernelFilename = MakeKernelFileName(KernelElec,KernelElecMod,KernelVdw,KernelVdwMod,KernelGeom) + '.c'
472                     fpkernel = open(KernelFilename,'w')
473                     defines['INCLUDE_HEADER'] = 1  # Include header first time in new file
474                     DoHeader = 1
475
476                     for KernelVF in VFList:
477
478                         KernelName = MakeKernelName(KernelElec,KernelElecMod,KernelVdw,KernelVdwMod,KernelGeom,KernelVF)
479
480                         defines['KERNEL_NAME'] = KernelName
481                         defines['KERNEL_VF']   = KernelVF
482
483                         # Check if this is a valid/sane/usable combination
484                         if not KeepKernel(KernelElec,KernelElecMod,KernelVdw,KernelVdwMod,KernelGeom,KernelVF):
485                             continue;
486
487                         # The overall kernel settings determine what the _kernel_ calculates, but for the water
488                         # kernels this does not mean that every pairwise interaction has e.g. Vdw interactions.
489                         # This routine sets defines of what to calculate for each pair of particles in those cases.
490                         SetDefines(KernelElec,KernelElecMod,KernelVdw,KernelVdwMod,KernelGeom,KernelVF,defines)
491
492                         if(DoHeader==1):
493                             fpkernel.write( FileHeader )
494
495                         gmxpreprocess('nb_kernel_template_' + Arch + '.pre', KernelName+'.tmp' , defines, force=1,contentType='C')
496                         numKernels = numKernels + 1
497
498                         defines['INCLUDE_HEADER'] = 0   # Header has been included once now
499                         DoHeader=0
500
501                         # Append temp file contents to the common kernelfile
502                         fptmp = open(KernelName+'.tmp','r')
503                         fpkernel.writelines(fptmp.readlines())
504                         fptmp.close()
505                         os.remove(KernelName+'.tmp')
506
507                         # Add a declaration for this kernel
508                         fpdecl.write('nb_kernel_t ' + KernelName + ';\n');
509
510                         # Add declaration to the buffer
511                         KernelOther=''
512                         kerneldecl.append(MakeKernelDecl(KernelName,KernelElec,KernelElecMod,KernelVdw,KernelVdwMod,KernelGeom,KernelOther,KernelVF))
513
514                     filesize = fpkernel.tell()
515                     fpkernel.close()
516                     if(filesize==0):
517                         os.remove(KernelFilename)
518
519                     PrintStatistics(cnt/ntot)
520                 pass
521             pass
522         pass
523     pass
524 pass
525
526 # Write out the list of settings and corresponding kernels to the declaration file
527 fpdecl.write( '\n\n' )
528 fpdecl.write( 'nb_kernel_info_t\n' )
529 fpdecl.write( 'kernellist_'+Arch+'[] =\n' )
530 fpdecl.write( '{\n' )
531 for decl in kerneldecl[0:-1]:
532     fpdecl.write( decl + ',\n' )
533 fpdecl.write( kerneldecl[-1] + '\n' )
534 fpdecl.write( '};\n\n' )
535 fpdecl.write( 'int\n' )
536 fpdecl.write( 'kernellist_'+Arch+'_size = sizeof(kernellist_'+Arch+')/sizeof(kernellist_'+Arch+'[0]);\n\n')
537 fpdecl.write( '#endif\n')
538 fpdecl.close()