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