pyslha.py

Tue, 14 Sep 2010 18:54:12 +0100

author
Andy Buckley <andy@insectnation.org>
date
Tue, 14 Sep 2010 18:54:12 +0100
changeset 48
8f42d73776cd
parent 47
9fa86bba2fe9
child 49
1f1923985bb3
permissions
-rw-r--r--

Moving ISAWIG writer functions into module

andy@1 1 #! /usr/bin/env python
andy@1 2
andy@8 3 """\
andy@8 4 A simple but flexible parser of SUSY Les Houches Accord (SLHA) model and decay files.
andy@25 5
andy@39 6 TODO: What to do about unlisted PIDs?
andy@33 7 TODO: Split writeSLHA into writeSLHA{Blocks,Decays}
andy@8 8 """
andy@8 9
andy@8 10 from __future__ import with_statement
andy@26 11
andy@8 12 __author__ = "Andy Buckley <andy.buckley@cern.ch"
andy@47 13 __version__ = "0.2.2"
andy@8 14
andy@1 15
andy@4 16 def _autotype(var):
andy@30 17 """Automatically convert strings to numerical types if possible."""
andy@30 18 if type(var) is not str:
andy@4 19 return var
andy@36 20 if var.isdigit() or (var.startswith("-") and var[1:].isdigit()):
andy@4 21 return int(var)
andy@4 22 try:
andy@4 23 f = float(var)
andy@4 24 return f
andy@4 25 except ValueError:
andy@4 26 return var
andy@4 27
andy@30 28 def _autostr(var):
andy@30 29 """Automatically numerical types to the right sort of string."""
andy@30 30 if type(var) is float:
andy@30 31 return "%e" % var
andy@30 32 return str(var)
andy@30 33
andy@30 34
andy@4 35
andy@12 36 class Block(object):
andy@8 37 """
andy@8 38 Object representation of any BLOCK elements read from the SLHA file. Blocks
andy@8 39 have a name, may have an associated Q value, and then a collection of data
andy@8 40 entries, stored as a recursive dictionary. Types in the dictionary are
andy@8 41 numeric (int or float) when a cast from the string in the file has been
andy@8 42 possible.
andy@8 43 """
andy@5 44 def __init__(self, name, q=None):
andy@1 45 self.name = name
andy@1 46 self.entries = {}
andy@5 47 self.q = _autotype(q)
andy@1 48
andy@1 49 def add_entry(self, entry):
andy@1 50 #print entry
andy@1 51 nextparent = self.entries
andy@1 52 if len(entry) < 2:
andy@1 53 raise Exception("Block entries must be at least a 2-tuple")
andy@4 54 #print "in", entry
andy@4 55 entry = map(_autotype, entry)
andy@4 56 #print "out", entry
andy@1 57 for e in entry[:-2]:
andy@1 58 if e is not entry[-1]:
andy@1 59 nextparent = nextparent.setdefault(e, {})
andy@1 60 nextparent[entry[-2]] = entry[-1]
andy@1 61 #print self.entries
andy@1 62
andy@1 63 def __cmp__(self, other):
andy@31 64 return cmp(self.name, other.name)
andy@1 65
andy@1 66 def __str__(self):
andy@1 67 s = self.name
andy@5 68 if self.q is not None:
andy@5 69 s += " (Q=%s)" % self.q
andy@1 70 s += "\n"
andy@1 71 s += str(self.entries)
andy@1 72 return s
andy@1 73
andy@1 74
andy@12 75 class Decay(object):
andy@8 76 """
andy@8 77 Object representing a decay entry on a particle decribed by the SLHA file.
andy@8 78 'Decay' objects are not a direct representation of a DECAY block in an SLHA
andy@8 79 file... that role, somewhat confusingly, is taken by the Particle class.
andy@8 80
andy@8 81 Decay objects have three properties: a branching ratio, br, an nda number
andy@12 82 (number of daughters == len(ids)), and a tuple of PDG PIDs to which the
andy@12 83 decay occurs. The PDG ID of the particle whose decay this represents may
andy@12 84 also be stored, but this is normally known via the Particle in which the
andy@12 85 decay is stored.
andy@8 86 """
andy@8 87 def __init__(self, br, nda, ids, parentid=None):
andy@8 88 self.parentid = parentid
andy@6 89 self.br = br
andy@6 90 self.nda = nda
andy@6 91 self.ids = ids
andy@29 92 assert(self.nda == len(self.ids))
andy@6 93
andy@6 94 def __cmp__(self, other):
andy@31 95 return cmp(other.br, self.br)
andy@6 96
andy@6 97 def __str__(self):
andy@12 98 return "%e %s" % (self.br, self.ids)
andy@6 99
andy@6 100
andy@12 101 class Particle(object):
andy@8 102 """
andy@8 103 Representation of a single, specific particle, decay block from an SLHA
andy@8 104 file. These objects are not themselves called 'Decay', since that concept
andy@8 105 applies more naturally to the various decays found inside this
andy@8 106 object. Particle classes store the PDG ID (pid) of the particle being
andy@8 107 represented, and optionally the mass (mass) and total decay width
andy@8 108 (totalwidth) of that particle in the SLHA scenario. Masses may also be found
andy@8 109 via the MASS block, from which the Particle.mass property is filled, if at
andy@8 110 all. They also store a list of Decay objects (decays) which are probably the
andy@8 111 item of most interest.
andy@8 112 """
andy@6 113 def __init__(self, pid, totalwidth=None, mass=None):
andy@6 114 self.pid = pid
andy@6 115 self.totalwidth = totalwidth
andy@6 116 self.mass = mass
andy@6 117 self.decays = []
andy@6 118
andy@6 119 def add_decay(self, br, nda, ids):
andy@6 120 self.decays.append(Decay(br, nda, ids))
andy@6 121 self.decays.sort()
andy@6 122
andy@6 123 def __cmp__(self, other):
andy@6 124 if abs(self.pid) == abs(other.pid):
andy@31 125 return cmp(self.pid, other.pid)
andy@31 126 return cmp(abs(self.pid), abs(other.pid))
andy@6 127
andy@6 128 def __str__(self):
andy@6 129 s = str(self.pid)
andy@7 130 if self.mass is not None:
andy@7 131 s += " : mass = %e GeV" % self.mass
andy@6 132 if self.totalwidth is not None:
andy@7 133 s += " : total width = %e GeV" % self.totalwidth
andy@6 134 for d in self.decays:
andy@12 135 if d.br > 0.0:
andy@12 136 s += "\n %s" % d
andy@6 137 return s
andy@1 138
andy@1 139
andy@31 140 def readSLHAFile(spcfilename, **kwargs):
andy@21 141 """
andy@21 142 Read an SLHA file, returning dictionaries of blocks and decays.
andy@31 143
andy@31 144 Other keyword parameters are passed to readSLHA.
andy@21 145 """
andy@21 146 with open(spcfilename, "r") as f:
andy@31 147 return readSLHA(f.read(), kwargs)
andy@21 148
andy@21 149
andy@31 150 def readSLHA(spcstr, ignorenobr=False):
andy@21 151 """
andy@31 152 Read an SLHA definition from a string, returning dictionaries of blocks and
andy@31 153 decays.
andy@31 154
andy@31 155 If the ignorenobr parameter is True, do not store decay entries with a
andy@31 156 branching ratio of zero.
andy@21 157 """
andy@1 158 blocks = {}
andy@1 159 decays = {}
andy@21 160 #
andy@34 161 import re
andy@21 162 currentblock = None
andy@21 163 currentdecay = None
andy@21 164 for line in spcstr.splitlines():
andy@21 165 ## Handle (ignore) comment lines
andy@21 166 if line.startswith("#"):
andy@21 167 continue
andy@21 168 if "#" in line:
andy@21 169 line = line[:line.index("#")]
andy@21 170
andy@21 171 ## Handle BLOCK/DECAY start lines
andy@21 172 if line.upper().startswith("BLOCK"):
andy@47 173 #print line
andy@47 174 match = re.match(r"BLOCK\s+(\w+)(\s+Q=\s*.+)?", line.upper())
andy@21 175 if not match:
andy@8 176 continue
andy@21 177 blockname = match.group(1)
andy@21 178 qstr = match.group(2)
andy@21 179 if qstr is not None:
andy@21 180 qstr = qstr[2:].strip()
andy@21 181 currentblock = blockname
andy@21 182 currentdecay = None
andy@21 183 blocks[blockname] = Block(blockname, q=qstr)
andy@21 184 elif line.upper().startswith("DECAY"):
andy@21 185 match = re.match(r"DECAY\s+(\d+)\s+([\d\.E+-]+).*", line.upper())
andy@21 186 if not match:
andy@21 187 continue
andy@21 188 pdgid = int(match.group(1))
andy@21 189 width = float(match.group(2))
andy@21 190 currentblock = "DECAY"
andy@21 191 currentdecay = pdgid
andy@21 192 decays[pdgid] = Particle(pdgid, width)
andy@21 193 else:
andy@21 194 ## In-block line
andy@21 195 if currentblock is not None:
andy@21 196 items = line.split()
andy@21 197 if len(items) < 1:
andy@6 198 continue
andy@21 199 if currentblock != "DECAY":
andy@21 200 if len(items) < 2:
andy@21 201 ## Treat the ALPHA block differently
andy@21 202 blocks[currentblock].value = _autotype(items[0])
andy@33 203 blocks[currentblock].entries = _autotype(items[0])
andy@8 204 else:
andy@21 205 blocks[currentblock].add_entry(items)
andy@21 206 else:
andy@21 207 br = float(items[0])
andy@21 208 nda = int(items[1])
andy@21 209 ids = map(int, items[2:])
andy@31 210 if br > 0.0 or not ignorenobr:
andy@31 211 decays[currentdecay].add_decay(br, nda, ids)
andy@1 212
andy@8 213 ## Try to populate Particle masses from the MASS block
andy@47 214 # print blocks.keys()
andy@47 215 try:
andy@47 216 for pid in blocks["MASS"].entries.keys():
andy@47 217 if decays.has_key(pid):
andy@47 218 decays[pid].mass = blocks["MASS"].entries[pid]
andy@47 219 except:
andy@47 220 raise Exception("No MASS block found, from which to populate particle masses")
andy@8 221
andy@1 222 return blocks, decays
andy@1 223
andy@1 224
andy@33 225 def readISAWIGFile(isafilename, **kwargs):
andy@33 226 """
andy@33 227 Read a spectrum definition from a file in the ISAWIG format, returning
andy@33 228 dictionaries of blocks and decays. While this is not an SLHA format, it is
andy@33 229 informally supported as a useful mechanism for converting ISAWIG spectra to
andy@33 230 SLHA.
andy@33 231
andy@33 232 Other keyword parameters are passed to readSLHA.
andy@33 233 """
andy@33 234 with open(isafilename, "r") as f:
andy@33 235 return readISAWIG(f.read(), kwargs)
andy@33 236
andy@33 237
andy@33 238 def readISAWIG(isastr, ignorenobr=False):
andy@33 239 """
andy@33 240 Read a spectrum definition from a string in the ISAWIG format, returning
andy@33 241 dictionaries of blocks and decays. While this is not an SLHA format, it is
andy@33 242 informally supported as a useful mechanism for converting ISAWIG spectra to
andy@33 243 SLHA.
andy@33 244
andy@33 245 ISAWIG parsing based on the HERWIG SUSY specification format, from
andy@33 246 http://www.hep.phy.cam.ac.uk/~richardn/HERWIG/ISAWIG/file.html
andy@33 247
andy@33 248 If the ignorenobr parameter is True, do not store decay entries with a
andy@33 249 branching ratio of zero.
andy@33 250 """
andy@33 251
andy@33 252 ## PDG MC ID codes mapped to HERWIG SUSY ID codes, based on
andy@33 253 ## http://www.hep.phy.cam.ac.uk/~richardn/HERWIG/ISAWIG/susycodes.html
andy@33 254 HERWIGID2PDGID = {}
andy@33 255 HERWIGID2PDGID[203] = 25 ## HIGGSL0
andy@33 256 HERWIGID2PDGID[204] = 35 ## HIGGSH0
andy@33 257 HERWIGID2PDGID[205] = 36 ## HIGGSA0
andy@33 258 HERWIGID2PDGID[206] = 37 ## HIGGS+
andy@33 259 HERWIGID2PDGID[207] = -37 ## HIGGS-
andy@33 260 HERWIGID2PDGID[401] = 1000001 ## SSDLBR
andy@33 261 HERWIGID2PDGID[407] = -1000001 ## SSDLBR
andy@33 262 HERWIGID2PDGID[402] = 1000002 ## SSULBR
andy@33 263 HERWIGID2PDGID[408] = -1000002 ## SSUL
andy@33 264 HERWIGID2PDGID[403] = 1000003 ## SSSLBR
andy@33 265 HERWIGID2PDGID[409] = -1000003 ## SSSL
andy@33 266 HERWIGID2PDGID[404] = 1000004 ## SSCLBR
andy@33 267 HERWIGID2PDGID[410] = -1000004 ## SSCL
andy@33 268 HERWIGID2PDGID[405] = 1000005 ## SSB1BR
andy@33 269 HERWIGID2PDGID[411] = -1000005 ## SSB1
andy@33 270 HERWIGID2PDGID[406] = 1000006 ## SST1BR
andy@33 271 HERWIGID2PDGID[412] = -1000006 ## SST1
andy@33 272 HERWIGID2PDGID[413] = 2000001 ## SSDR
andy@33 273 HERWIGID2PDGID[419] = -2000001 ## SSDRBR
andy@33 274 HERWIGID2PDGID[414] = 2000002 ## SSUR
andy@33 275 HERWIGID2PDGID[420] = -2000002 ## SSURBR
andy@33 276 HERWIGID2PDGID[415] = 2000003 ## SSSR
andy@33 277 HERWIGID2PDGID[421] = -2000003 ## SSSRBR
andy@33 278 HERWIGID2PDGID[416] = 2000004 ## SSCR
andy@33 279 HERWIGID2PDGID[422] = -2000004 ## SSCRBR
andy@33 280 HERWIGID2PDGID[417] = 2000005 ## SSB2
andy@33 281 HERWIGID2PDGID[423] = -2000005 ## SSB2BR
andy@33 282 HERWIGID2PDGID[418] = 2000006 ## SST2
andy@33 283 HERWIGID2PDGID[424] = -2000006 ## SST2BR
andy@33 284 HERWIGID2PDGID[425] = 1000011 ## SSEL-
andy@33 285 HERWIGID2PDGID[431] = -1000011 ## SSEL+
andy@33 286 HERWIGID2PDGID[426] = 1000012 ## SSNUEL
andy@33 287 HERWIGID2PDGID[432] = -1000012 ## SSNUELBR
andy@33 288 HERWIGID2PDGID[427] = 1000013 ## SSMUL-
andy@33 289 HERWIGID2PDGID[433] = -1000013 ## SSMUL+
andy@33 290 HERWIGID2PDGID[428] = 1000014 ## SSNUMUL
andy@33 291 HERWIGID2PDGID[434] = -1000014 ## SSNUMLBR
andy@33 292 HERWIGID2PDGID[429] = 1000015 ## SSTAU1-
andy@33 293 HERWIGID2PDGID[435] = -1000015 ## SSTAU1+
andy@33 294 HERWIGID2PDGID[430] = 1000016 ## SSNUTL
andy@33 295 HERWIGID2PDGID[436] = -1000016 ## SSNUTLBR
andy@33 296 HERWIGID2PDGID[437] = 2000011 ## SSEL-
andy@33 297 HERWIGID2PDGID[443] = -2000011 ## SSEL+
andy@33 298 HERWIGID2PDGID[438] = 2000012 ## SSNUEL
andy@33 299 HERWIGID2PDGID[444] = -2000012 ## SSNUELBR
andy@33 300 HERWIGID2PDGID[439] = 2000013 ## SSMUL-
andy@33 301 HERWIGID2PDGID[445] = -2000013 ## SSMUL+
andy@33 302 HERWIGID2PDGID[440] = 2000014 ## SSNUMUL
andy@33 303 HERWIGID2PDGID[446] = -2000014 ## SSNUMLBR
andy@33 304 HERWIGID2PDGID[441] = 2000015 ## SSTAU1-
andy@33 305 HERWIGID2PDGID[447] = -2000015 ## SSTAU1+
andy@33 306 HERWIGID2PDGID[442] = 2000016 ## SSNUTL
andy@33 307 HERWIGID2PDGID[448] = -2000016 ## SSNUTLBR
andy@33 308 HERWIGID2PDGID[449] = 1000021 ## GLUINO
andy@33 309 HERWIGID2PDGID[450] = 1000022 ## NTLINO1
andy@33 310 HERWIGID2PDGID[451] = 1000023 ## NTLINO2
andy@33 311 HERWIGID2PDGID[452] = 1000025 ## NTLINO3
andy@33 312 HERWIGID2PDGID[453] = 1000035 ## NTLINO4
andy@33 313 HERWIGID2PDGID[454] = 1000024 ## CHGINO1+
andy@33 314 HERWIGID2PDGID[456] = -1000024 ## CHGINO1-
andy@33 315 HERWIGID2PDGID[455] = 1000037 ## CHGINO2+
andy@33 316 HERWIGID2PDGID[457] = -1000037 ## CHGINO2-
andy@33 317 HERWIGID2PDGID[458] = 1000039 ## GRAVTINO
andy@33 318
andy@33 319 blocks = {}
andy@33 320 decays = {}
andy@35 321 LINES = isastr.splitlines()
andy@33 322
andy@33 323 def getnextvalidline():
andy@35 324 while LINES:
andy@35 325 s = LINES.pop(0).strip()
andy@33 326 ## Return None if EOF reached
andy@33 327 if len(s) == 0:
andy@33 328 continue
andy@33 329 ## Strip comments
andy@33 330 if "#" in s:
andy@33 331 s = s[:s.index("#")].strip()
andy@33 332 ## Return if non-empty
andy@33 333 if len(s) > 0:
andy@33 334 return s
andy@33 335
andy@34 336 def getnextvalidlineitems():
andy@34 337 return map(_autotype, getnextvalidline().split())
andy@34 338
andy@34 339 ## Populate MASS block and create decaying particle objects
andy@35 340 masses = Block("MASS")
andy@33 341 numentries = int(getnextvalidline())
andy@33 342 for i in xrange(numentries):
andy@34 343 hwid, mass, lifetime = getnextvalidlineitems()
andy@34 344 width = 1.0/(lifetime * 1.51926778e24) ## width in GeV == hbar/lifetime in seconds
andy@34 345 pdgid = HERWIGID2PDGID.get(hwid, hwid)
andy@34 346 masses.add_entry((pdgid, mass))
andy@34 347 decays[pdgid] = Particle(pdgid, width, mass)
andy@34 348 #print pdgid, mass, width
andy@34 349 blocks["MASS"] = masses
andy@33 350
andy@34 351 ## Populate decays
andy@34 352 for n in xrange(numentries):
andy@34 353 numdecays = int(getnextvalidline())
andy@34 354 for d in xrange(numdecays):
andy@40 355 #print n, numentries-1, d, numdecays-1
andy@34 356 decayitems = getnextvalidlineitems()
andy@34 357 hwid = decayitems[0]
andy@34 358 pdgid = HERWIGID2PDGID.get(hwid, hwid)
andy@34 359 br = decayitems[1]
andy@34 360 nme = decayitems[2]
andy@34 361 daughter_hwids = decayitems[3:]
andy@34 362 daughter_pdgids = []
andy@34 363 for hw in daughter_hwids:
andy@34 364 if hw != 0:
andy@34 365 daughter_pdgids.append(HERWIGID2PDGID.get(hw, hw))
andy@35 366 if not decays.has_key(pdgid):
andy@40 367 #print "Decay for unlisted particle %d, %d" % (hwid, pdgid)
andy@38 368 decays[pdgid] = Particle(pdgid)
andy@38 369 decays[pdgid].add_decay(br, len(daughter_pdgids), daughter_pdgids)
andy@33 370
andy@33 371
andy@34 372 ## Now the SUSY parameters
andy@34 373 TANB, ALPHAH = getnextvalidlineitems()
andy@34 374 blocks["MINPAR"] = Block("MINPAR")
andy@34 375 blocks["MINPAR"].add_entry((3, TANB))
andy@34 376 blocks["ALPHA"] = Block("ALPHA")
andy@34 377 blocks["ALPHA"].entries = ALPHAH
andy@34 378 #
andy@34 379 ## Neutralino mixing matrix
andy@34 380 blocks["NMIX"] = Block("NMIX")
andy@34 381 for i in xrange(1, 5):
andy@34 382 nmix_i = getnextvalidlineitems()
andy@34 383 for j, v in enumerate(nmix_i):
andy@34 384 blocks["NMIX"].add_entry((i, j+1, v))
andy@34 385 #
andy@34 386 ## Chargino mixing matrices V and U
andy@34 387 blocks["VMIX"] = Block("VMIX")
andy@34 388 vmix = getnextvalidlineitems()
andy@34 389 blocks["VMIX"].add_entry((1, 1, vmix[0]))
andy@34 390 blocks["VMIX"].add_entry((1, 2, vmix[1]))
andy@34 391 blocks["VMIX"].add_entry((2, 1, vmix[2]))
andy@34 392 blocks["VMIX"].add_entry((2, 2, vmix[3]))
andy@34 393 blocks["UMIX"] = Block("UMIX")
andy@34 394 umix = getnextvalidlineitems()
andy@34 395 blocks["UMIX"].add_entry((1, 1, umix[0]))
andy@34 396 blocks["UMIX"].add_entry((1, 2, umix[1]))
andy@34 397 blocks["UMIX"].add_entry((2, 1, umix[2]))
andy@34 398 blocks["UMIX"].add_entry((2, 2, umix[3]))
andy@34 399 #
andy@34 400 THETAT, THETAB, THETAL = getnextvalidlineitems()
andy@34 401 import math
andy@34 402 blocks["STOPMIX"] = Block("STOPMIX")
andy@34 403 blocks["STOPMIX"].add_entry((1, 1, math.cos(THETAT)))
andy@34 404 blocks["STOPMIX"].add_entry((1, 2, -math.sin(THETAT)))
andy@34 405 blocks["STOPMIX"].add_entry((2, 1, math.sin(THETAT)))
andy@34 406 blocks["STOPMIX"].add_entry((2, 2, math.cos(THETAT)))
andy@34 407 blocks["SBOTMIX"] = Block("SBOTMIX")
andy@34 408 blocks["SBOTMIX"].add_entry((1, 1, math.cos(THETAB)))
andy@34 409 blocks["SBOTMIX"].add_entry((1, 2, -math.sin(THETAB)))
andy@34 410 blocks["SBOTMIX"].add_entry((2, 1, math.sin(THETAB)))
andy@34 411 blocks["SBOTMIX"].add_entry((2, 2, math.cos(THETAB)))
andy@34 412 blocks["STAUMIX"] = Block("STAUMIX")
andy@34 413 blocks["STAUMIX"].add_entry((1, 1, math.cos(THETAL)))
andy@34 414 blocks["STAUMIX"].add_entry((1, 2, -math.sin(THETAL)))
andy@34 415 blocks["STAUMIX"].add_entry((2, 1, math.sin(THETAL)))
andy@34 416 blocks["STAUMIX"].add_entry((2, 2, math.cos(THETAL)))
andy@34 417 #
andy@34 418 ATSS, ABSS, ALSS = getnextvalidlineitems()
andy@34 419 blocks["AU"] = Block("AU")
andy@34 420 blocks["AU"].add_entry((3, 3, ATSS))
andy@34 421 blocks["AD"] = Block("AD")
andy@34 422 blocks["AD"].add_entry((3, 3, ABSS))
andy@34 423 blocks["AE"] = Block("AE")
andy@34 424 blocks["AE"].add_entry((3, 3, ALSS))
andy@34 425 #
andy@34 426 MUSS = getnextvalidlineitems()[0]
andy@34 427 blocks["MINPAR"].add_entry((4, MUSS))
andy@34 428 #
andy@34 429 return blocks, decays
andy@33 430
andy@21 431
andy@31 432 def writeSLHAFile(spcfilename, blocks, decays, **kwargs):
andy@29 433 """
andy@29 434 Write an SLHA file from the supplied blocks and decays dicts.
andy@31 435
andy@31 436 Other keyword parameters are passed to writeSLHA.
andy@29 437 """
andy@29 438 with open(spcfilename, "w") as f:
andy@31 439 f.write(writeSLHA(blocks, decays, kwargs))
andy@29 440
andy@29 441
andy@33 442 ## TODO: Split writeSLHA into writeSLHA{Blocks,Decays}
andy@33 443
andy@31 444 def writeSLHA(blocks, decays, ignorenobr=False):
andy@29 445 """
andy@29 446 Return an SLHA definition as a string, from the supplied blocks and decays dicts.
andy@29 447 """
andy@31 448 sep = " "
andy@29 449 out = ""
andy@30 450 def dict_hier_strs(d, s=""):
andy@30 451 if type(d) is dict:
andy@30 452 for k, v in sorted(d.iteritems()):
andy@31 453 for s2 in dict_hier_strs(v, s + sep + _autostr(k)):
andy@30 454 yield s2
andy@30 455 else:
andy@31 456 yield s + sep + _autostr(d)
andy@30 457 ## Blocks
andy@30 458 for bname, b in sorted(blocks.iteritems()):
andy@29 459 namestr = b.name
andy@29 460 if b.q is not None:
andy@29 461 namestr += " Q= %e" % b.q
andy@29 462 out += "BLOCK %s\n" % namestr
andy@30 463 for s in dict_hier_strs(b.entries):
andy@31 464 out += sep + s + "\n"
andy@29 465 out += "\n"
andy@30 466 ## Decays
andy@30 467 for pid, particle in sorted(decays.iteritems()):
andy@40 468 out += "DECAY %d %e\n" % (particle.pid, particle.totalwidth or -1)
andy@30 469 for d in sorted(particle.decays):
andy@31 470 if d.br > 0.0 or not ignorenobr:
andy@31 471 products_str = " ".join(map(str, d.ids))
andy@31 472 out += sep + "%e" % d.br + sep + "%d" % len(d.ids) + sep + products_str + "\n"
andy@29 473 out += "\n"
andy@29 474 return out
andy@29 475
andy@29 476
andy@29 477
andy@48 478 def writeISAWIGFile(isafilename, blocks, decays, **kwargs):
andy@48 479 """
andy@48 480 Write an ISAWIG file from the supplied blocks and decays dicts.
andy@48 481
andy@48 482 Other keyword parameters are passed to writeISAWIG.
andy@48 483
andy@48 484 TODO: Handle RPV SUSY
andy@48 485 """
andy@48 486 with open(isafilename, "w") as f:
andy@48 487 f.write(writeISAWIG(blocks, decays, kwargs))
andy@48 488
andy@48 489
andy@48 490 def writeISAWIG(blocks, decays, ignorenobr=False):
andy@48 491 """
andy@48 492 Return an ISAWIG definition as a string, from the supplied blocks and decays dicts.
andy@48 493
andy@48 494 ISAWIG parsing based on the HERWIG SUSY specification format, from
andy@48 495 http://www.hep.phy.cam.ac.uk/~richardn/HERWIG/ISAWIG/file.html
andy@48 496
andy@48 497 If the ignorenobr parameter is True, do not write decay entries with a
andy@48 498 branching ratio of zero.
andy@48 499 """
andy@48 500 ## PDG MC ID codes mapped to HERWIG SUSY ID codes, based on
andy@48 501 ## http://www.hep.phy.cam.ac.uk/~richardn/HERWIG/ISAWIG/susycodes.html
andy@48 502 PDGID2HERWIGID = {}
andy@48 503 PDGID2HERWIGID[ 25] = 203 ## HIGGSL0 (ADDED)
andy@48 504 PDGID2HERWIGID[ 26] = 203 ## HIGGSL0
andy@48 505 PDGID2HERWIGID[ 35] = 204 ## HIGGSH0
andy@48 506 PDGID2HERWIGID[ 36] = 205 ## HIGGSA0
andy@48 507 PDGID2HERWIGID[ 37] = 206 ## HIGGS+
andy@48 508 PDGID2HERWIGID[ -37] = 207 ## HIGGS-
andy@48 509 PDGID2HERWIGID[ 1000001] = 401 ## SSDLBR
andy@48 510 PDGID2HERWIGID[-1000001] = 407 ## SSDLBR
andy@48 511 PDGID2HERWIGID[ 1000002] = 402 ## SSULBR
andy@48 512 PDGID2HERWIGID[-1000002] = 408 ## SSUL
andy@48 513 PDGID2HERWIGID[ 1000003] = 403 ## SSSLBR
andy@48 514 PDGID2HERWIGID[-1000003] = 409 ## SSSL
andy@48 515 PDGID2HERWIGID[ 1000004] = 404 ## SSCLBR
andy@48 516 PDGID2HERWIGID[-1000004] = 410 ## SSCL
andy@48 517 PDGID2HERWIGID[ 1000005] = 405 ## SSB1BR
andy@48 518 PDGID2HERWIGID[-1000005] = 411 ## SSB1
andy@48 519 PDGID2HERWIGID[ 1000006] = 406 ## SST1BR
andy@48 520 PDGID2HERWIGID[-1000006] = 412 ## SST1
andy@48 521 PDGID2HERWIGID[ 2000001] = 413 ## SSDR
andy@48 522 PDGID2HERWIGID[-2000001] = 419 ## SSDRBR
andy@48 523 PDGID2HERWIGID[ 2000002] = 414 ## SSUR
andy@48 524 PDGID2HERWIGID[-2000002] = 420 ## SSURBR
andy@48 525 PDGID2HERWIGID[ 2000003] = 415 ## SSSR
andy@48 526 PDGID2HERWIGID[-2000003] = 421 ## SSSRBR
andy@48 527 PDGID2HERWIGID[ 2000004] = 416 ## SSCR
andy@48 528 PDGID2HERWIGID[-2000004] = 422 ## SSCRBR
andy@48 529 PDGID2HERWIGID[ 2000005] = 417 ## SSB2
andy@48 530 PDGID2HERWIGID[-2000005] = 423 ## SSB2BR
andy@48 531 PDGID2HERWIGID[ 2000006] = 418 ## SST2
andy@48 532 PDGID2HERWIGID[-2000006] = 424 ## SST2BR
andy@48 533 PDGID2HERWIGID[ 1000011] = 425 ## SSEL-
andy@48 534 PDGID2HERWIGID[-1000011] = 431 ## SSEL+
andy@48 535 PDGID2HERWIGID[ 1000012] = 426 ## SSNUEL
andy@48 536 PDGID2HERWIGID[-1000012] = 432 ## SSNUELBR
andy@48 537 PDGID2HERWIGID[ 1000013] = 427 ## SSMUL-
andy@48 538 PDGID2HERWIGID[-1000013] = 433 ## SSMUL+
andy@48 539 PDGID2HERWIGID[ 1000014] = 428 ## SSNUMUL
andy@48 540 PDGID2HERWIGID[-1000014] = 434 ## SSNUMLBR
andy@48 541 PDGID2HERWIGID[ 1000015] = 429 ## SSTAU1-
andy@48 542 PDGID2HERWIGID[-1000015] = 435 ## SSTAU1+
andy@48 543 PDGID2HERWIGID[ 1000016] = 430 ## SSNUTL
andy@48 544 PDGID2HERWIGID[-1000016] = 436 ## SSNUTLBR
andy@48 545 PDGID2HERWIGID[ 2000011] = 437 ## SSEL-
andy@48 546 PDGID2HERWIGID[-2000011] = 443 ## SSEL+
andy@48 547 PDGID2HERWIGID[ 2000012] = 438 ## SSNUEL
andy@48 548 PDGID2HERWIGID[-2000012] = 444 ## SSNUELBR
andy@48 549 PDGID2HERWIGID[ 2000013] = 439 ## SSMUL-
andy@48 550 PDGID2HERWIGID[-2000013] = 445 ## SSMUL+
andy@48 551 PDGID2HERWIGID[ 2000014] = 440 ## SSNUMUL
andy@48 552 PDGID2HERWIGID[-2000014] = 446 ## SSNUMLBR
andy@48 553 PDGID2HERWIGID[ 2000015] = 441 ## SSTAU1-
andy@48 554 PDGID2HERWIGID[-2000015] = 447 ## SSTAU1+
andy@48 555 PDGID2HERWIGID[ 2000016] = 442 ## SSNUTL
andy@48 556 PDGID2HERWIGID[-2000016] = 448 ## SSNUTLBR
andy@48 557 PDGID2HERWIGID[ 1000021] = 449 ## GLUINO
andy@48 558 PDGID2HERWIGID[ 1000022] = 450 ## NTLINO1
andy@48 559 PDGID2HERWIGID[ 1000023] = 451 ## NTLINO2
andy@48 560 PDGID2HERWIGID[ 1000025] = 452 ## NTLINO3
andy@48 561 PDGID2HERWIGID[ 1000035] = 453 ## NTLINO4
andy@48 562 PDGID2HERWIGID[ 1000024] = 454 ## CHGINO1+
andy@48 563 PDGID2HERWIGID[-1000024] = 456 ## CHGINO1-
andy@48 564 PDGID2HERWIGID[ 1000037] = 455 ## CHGINO2+
andy@48 565 PDGID2HERWIGID[-1000037] = 457 ## CHGINO2-
andy@48 566 PDGID2HERWIGID[ 1000039] = 458 ## GRAVTINO
andy@48 567
andy@48 568 masses = blocks["MASS"].entries
andy@48 569
andy@48 570 ## Init output string
andy@48 571 out = ""
andy@48 572
andy@48 573 ## First write out masses section:
andy@48 574 ## Number of SUSY + top particles
andy@48 575 ## IDHW, RMASS(IDHW), RLTIM(IDHW)
andy@48 576 ## repeated for each particle
andy@48 577 ## IDHW is the HERWIG identity code.
andy@48 578 ## RMASS and RTLIM are the mass in GeV, and lifetime in seconds respectively.
andy@48 579 massout = ""
andy@48 580 for pid in masses.keys():
andy@48 581 lifetime = -1
andy@48 582 try:
andy@48 583 width = decays[pid].totalwidth
andy@48 584 if width and width > 0:
andy@48 585 lifetime = 1.0/(width * 1.51926778e24) ## lifetime in seconds == hbar/width in GeV
andy@48 586 except:
andy@48 587 pass
andy@48 588 massout += "%d %e %e\n" % (PDGID2HERWIGID.get(pid, pid), masses[pid], lifetime)
andy@48 589 out += "%d\n" % massout.count("\n")
andy@48 590 out += massout
andy@48 591
andy@48 592 assert(len(masses) == len(decays))
andy@48 593
andy@48 594 ## Next each particles decay modes together with their branching ratios and matrix element codes
andy@48 595 ## Number of decay modes for a given particle (IDK)
andy@48 596 ## IDK(*), BRFRAC(*), NME(*) & IDKPRD(1-5,*)
andy@48 597 ## repeated for each mode.
andy@48 598 ## Repeated for each particle.
andy@48 599 ## IDK is the HERWIG code for the decaying particle, BRFRAC is the branching ratio of
andy@48 600 ## the decay mode. NME is a code for the matrix element to be used, either from the
andy@48 601 ## SUSY elements or the main HERWIG MEs. IDKPRD are the HERWIG identity codes of the decay products.
andy@48 602 for i, pid in enumerate(decays.keys()):
andy@48 603 # if not decays.has_key(pid):
andy@48 604 # continue
andy@48 605 hwid = PDGID2HERWIGID.get(pid, pid)
andy@48 606 decayout = ""
andy@48 607 #decayout += "@@@@ %d %d %d\n" % (i, pid, hwid)
andy@48 608 for i_d, d in enumerate(decays[pid].decays):
andy@48 609 ## Skip decay if it has no branching ratio
andy@48 610 if ignorenobr and d.br == 0:
andy@48 611 continue
andy@48 612 ## TODO: Identify decay matrix element to use
andy@48 613 #decayout += str(i_d + 1) + " "
andy@48 614 decayout += "%d %e 0 " % (hwid, d.br)
andy@48 615 ## TODO: Order decay products as required
andy@48 616 ids = [0,0,0,0,0]
andy@48 617 for i, pid in enumerate(d.ids):
andy@48 618 ids[i] = pid
andy@48 619 ids = map(str, ids)
andy@48 620 decayout += " ".join(ids) + "\n"
andy@48 621 decayout = "%d\n" % decayout.count("\n") + decayout
andy@48 622 out += decayout
andy@48 623
andy@48 624 ## Now the SUSY parameters
andy@48 625 ## TANB, ALPHAH:
andy@48 626 out += "%e %e\n" % (blocks["MINPAR"].entries[3], blocks["ALPHA"].entries)
andy@48 627 ## Neutralino mixing matrix
andy@48 628 nmix = blocks["NMIX"].entries
andy@48 629 for i in xrange(1, 5):
andy@48 630 out += "%e %e %e %e\n" % (nmix[i][1], nmix[i][2], nmix[i][3], nmix[i][4])
andy@48 631 ## Chargino mixing matrices V and U
andy@48 632 vmix = blocks["VMIX"].entries
andy@48 633 out += "%e %e %e %e\n" % (vmix[1][1], vmix[1][2], vmix[2][1], vmix[2][2])
andy@48 634 umix = blocks["UMIX"].entries
andy@48 635 out += "%e %e %e %e\n" % (umix[1][1], umix[1][2], umix[2][1], umix[2][2])
andy@48 636 # THETAT,THETAB,THETAL
andy@48 637 import math
andy@48 638 out += "%e %e %e\n" % (math.acos(blocks["STOPMIX"].entries[1][1]),
andy@48 639 math.acos(blocks["SBOTMIX"].entries[1][1]),
andy@48 640 math.acos(blocks["STAUMIX"].entries[1][1]))
andy@48 641 # ATSS,ABSS,ALSS
andy@48 642 out += "%e %e %e\n" % (blocks["AU"].entries[3][3],
andy@48 643 blocks["AD"].entries[3][3],
andy@48 644 blocks["AE"].entries[3][3])
andy@48 645 # MUSS == sign(mu)
andy@48 646 out += "%f\n" % blocks["MINPAR"].entries[4]
andy@48 647
andy@48 648 ## TODO: Handle RPV SUSY
andy@48 649
andy@48 650 return out
andy@48 651
andy@48 652
andy@48 653
andy@1 654 if __name__ == "__main__":
andy@1 655 import sys
andy@1 656 for a in sys.argv[1:]:
andy@35 657 if a.endswith(".isa"):
andy@35 658 blocks, decays = readISAWIGFile(a)
andy@35 659 else:
andy@35 660 blocks, decays = readSLHAFile(a)
andy@3 661
andy@5 662 for bname, b in sorted(blocks.iteritems()):
andy@5 663 print b
andy@5 664 print
andy@3 665
andy@47 666 print blocks.keys()
andy@47 667
andy@4 668 print blocks["MASS"].entries[25]
andy@6 669 print
andy@6 670
andy@6 671 for p in sorted(decays.values()):
andy@6 672 print p
andy@6 673 print
andy@29 674
andy@31 675 print writeSLHA(blocks, decays, ignorenobr=True)

mercurial