Package madgraph :: Package various :: Module lhe_parser
[hide private]
[frames] | no frames]

Source Code for Module madgraph.various.lhe_parser

   1  from __future__ import division 
   2  import collections 
   3  import random 
   4  import re 
   5  import operator 
   6  import numbers 
   7  import math 
   8  import time 
   9  import os 
  10  import shutil 
  11  import sys 
  12   
  13  pjoin = os.path.join 
  14   
  15  if '__main__' == __name__: 
  16      import sys 
  17      sys.path.append('../../') 
  18  import misc 
  19  import logging 
  20  import gzip 
  21  import banner as banner_mod 
  22  logger = logging.getLogger("madgraph.lhe_parser") 
23 24 -class Particle(object):
25 """ """ 26 # regular expression not use anymore to speed up the computation 27 #pattern=re.compile(r'''^\s* 28 # (?P<pid>-?\d+)\s+ #PID 29 # (?P<status>-?\d+)\s+ #status (1 for output particle) 30 # (?P<mother1>-?\d+)\s+ #mother 31 # (?P<mother2>-?\d+)\s+ #mother 32 # (?P<color1>[+-e.\d]*)\s+ #color1 33 # (?P<color2>[+-e.\d]*)\s+ #color2 34 # (?P<px>[+-e.\d]*)\s+ #px 35 # (?P<py>[+-e.\d]*)\s+ #py 36 # (?P<pz>[+-e.\d]*)\s+ #pz 37 # (?P<E>[+-e.\d]*)\s+ #E 38 # (?P<mass>[+-e.\d]*)\s+ #mass 39 # (?P<vtim>[+-e.\d]*)\s+ #displace vertex 40 # (?P<helicity>[+-e.\d]*)\s* #helicity 41 # ($|(?P<comment>\#[\d|D]*)) #comment/end of string 42 # ''',66) #verbose+ignore case 43 44 45
46 - def __init__(self, line=None, event=None):
47 """ """ 48 49 if isinstance(line, Particle): 50 for key in line.__dict__: 51 setattr(self, key, getattr(line, key)) 52 if event: 53 self.event = event 54 return 55 56 self.event = event 57 self.event_id = len(event) #not yet in the event 58 # LHE information 59 self.pid = 0 60 self.status = 0 # -1:initial. 1:final. 2: propagator 61 self.mother1 = None 62 self.mother2 = None 63 self.color1 = 0 64 self.color2 = None 65 self.px = 0 66 self.py = 0 67 self.pz = 0 68 self.E = 0 69 self.mass = 0 70 self.vtim = 0 71 self.helicity = 9 72 self.rwgt = 0 73 self.comment = '' 74 75 if line: 76 self.parse(line)
77 78 @property
79 - def pdg(self):
80 "convenient alias" 81 return self.pid
82
83 - def parse(self, line):
84 """parse the line""" 85 86 args = line.split() 87 keys = ['pid', 'status','mother1','mother2','color1', 'color2', 'px','py','pz','E', 88 'mass','vtim','helicity'] 89 90 for key,value in zip(keys,args): 91 setattr(self, key, float(value)) 92 self.pid = int(self.pid) 93 94 self.comment = ' '.join(args[len(keys):]) 95 if self.comment.startswith(('|','#')): 96 self.comment = self.comment[1:]
97 98 # Note that mother1/mother2 will be modified by the Event parse function to replace the 99 # integer by a pointer to the actual particle object. 100
101 - def __str__(self):
102 """string representing the particles""" 103 return " %8d %2d %4d %4d %4d %4d %+13.10e %+13.10e %+13.10e %14.10e %14.10e %10.4e %10.4e" \ 104 % (self.pid, 105 self.status, 106 (self.mother1 if isinstance(self.mother1, numbers.Number) else self.mother1.event_id+1) if self.mother1 else 0, 107 (self.mother2 if isinstance(self.mother2, numbers.Number) else self.mother2.event_id+1) if self.mother2 else 0, 108 self.color1, 109 self.color2, 110 self.px, 111 self.py, 112 self.pz, 113 self.E, 114 self.mass, 115 self.vtim, 116 self.helicity)
117
118 - def __eq__(self, other):
119 120 if not isinstance(other, Particle): 121 return False 122 if self.pid == other.pid and \ 123 self.status == other.status and \ 124 self.mother1 == other.mother1 and \ 125 self.mother2 == other.mother2 and \ 126 self.color1 == other.color1 and \ 127 self.color2 == other.color2 and \ 128 self.px == other.px and \ 129 self.py == other.py and \ 130 self.pz == other.pz and \ 131 self.E == other.E and \ 132 self.mass == other.mass and \ 133 self.vtim == other.vtim and \ 134 self.helicity == other.helicity: 135 return True 136 return False
137
138 - def set_momentum(self, momentum):
139 140 self.E = momentum.E 141 self.px = momentum.px 142 self.py = momentum.py 143 self.pz = momentum.pz
144
145 - def add_decay(self, decay_event):
146 """associate to this particle the decay in the associate event""" 147 148 return self.event.add_decay_to_particle(self.event_id, decay_event)
149
150 - def __repr__(self):
151 return 'Particle("%s", event=%s)' % (str(self), self.event)
152
153 154 -class EventFile(object):
155 """A class to allow to read both gzip and not gzip file""" 156 157
158 - def __new__(self, path, mode='r', *args, **opt):
159 160 if not path.endswith(".gz"): 161 return file.__new__(EventFileNoGzip, path, mode, *args, **opt) 162 elif mode == 'r' and not os.path.exists(path) and os.path.exists(path[:-3]): 163 return EventFile.__new__(EventFileNoGzip, path[:-3], mode, *args, **opt) 164 else: 165 try: 166 return gzip.GzipFile.__new__(EventFileGzip, path, mode, *args, **opt) 167 except IOError, error: 168 raise 169 except Exception, error: 170 if mode == 'r': 171 misc.gunzip(path) 172 return file.__new__(EventFileNoGzip, path[:-3], mode, *args, **opt)
173 174
175 - def __init__(self, path, mode='r', *args, **opt):
176 """open file and read the banner [if in read mode]""" 177 178 self.to_zip = False 179 if path.endswith('.gz') and mode == 'w' and\ 180 isinstance(self, EventFileNoGzip): 181 path = path[:-3] 182 self.to_zip = True # force to zip the file at close() with misc.gzip 183 184 self.parsing = True # check if/when we need to parse the event. 185 self.eventgroup = False 186 try: 187 super(EventFile, self).__init__(path, mode, *args, **opt) 188 except IOError: 189 if '.gz' in path and isinstance(self, EventFileNoGzip) and\ 190 mode == 'r' and os.path.exists(path[:-3]): 191 super(EventFile, self).__init__(path[:-3], mode, *args, **opt) 192 else: 193 raise 194 195 self.banner = '' 196 if mode == 'r': 197 line = '' 198 while '</init>' not in line.lower(): 199 try: 200 line = super(EventFile, self).next() 201 except StopIteration: 202 self.seek(0) 203 self.banner = '' 204 break 205 if "<event" in line.lower(): 206 self.seek(0) 207 self.banner = '' 208 break 209 210 self.banner += line
211
212 - def get_banner(self):
213 """return a banner object""" 214 import madgraph.various.banner as banner 215 if isinstance(self.banner, banner.Banner): 216 return self.banner 217 218 output = banner.Banner() 219 output.read_banner(self.banner) 220 return output
221 222 @property
223 - def cross(self):
224 """return the cross-section of the file #from the banner""" 225 try: 226 return self._cross 227 except Exception: 228 pass 229 230 onebanner = self.get_banner() 231 self._cross = onebanner.get_cross() 232 return self._cross
233
234 - def __len__(self):
235 if self.closed: 236 return 0 237 if hasattr(self,"len"): 238 return self.len 239 self.seek(0) 240 nb_event=0 241 with misc.TMP_variable(self, 'parsing', False): 242 for _ in self: 243 nb_event +=1 244 self.len = nb_event 245 self.seek(0) 246 return self.len
247
248 - def next(self):
249 """get next event""" 250 251 if not self.eventgroup: 252 text = '' 253 line = '' 254 mode = 0 255 while '</event>' not in line: 256 line = super(EventFile, self).next() 257 if '<event' in line: 258 mode = 1 259 text = '' 260 if mode: 261 text += line 262 if self.parsing: 263 out = Event(text) 264 if len(out) == 0: 265 raise Exception 266 return out 267 else: 268 return text 269 else: 270 events = [] 271 text = '' 272 line = '' 273 mode = 0 274 while '</eventgroup>' not in line: 275 line = super(EventFile, self).next() 276 if '<eventgroup' in line: 277 events=[] 278 text = '' 279 elif '<event' in line: 280 text='' 281 mode=1 282 elif '</event>' in line: 283 if self.parsing: 284 events.append(Event(text)) 285 else: 286 events.append(text) 287 text = '' 288 mode = 0 289 if mode: 290 text += line 291 if len(events) == 0: 292 return self.next() 293 return events
294 295
296 - def initialize_unweighting(self, get_wgt, trunc_error):
297 """ scan once the file to return 298 - the list of the hightest weight (of size trunc_error*NB_EVENT 299 - the cross-section by type of process 300 - the total number of events in the file 301 """ 302 303 # We need to loop over the event file to get some information about the 304 # new cross-section/ wgt of event. 305 self.seek(0) 306 all_wgt = [] 307 cross = collections.defaultdict(int) 308 nb_event = 0 309 for event in self: 310 nb_event +=1 311 wgt = get_wgt(event) 312 cross['all'] += wgt 313 cross['abs'] += abs(wgt) 314 cross[event.ievent] += wgt 315 all_wgt.append(abs(wgt)) 316 # avoid all_wgt to be too large 317 if nb_event % 20000 == 0: 318 all_wgt.sort() 319 # drop the lowest weight 320 nb_keep = max(20, int(nb_event*trunc_error*15)) 321 all_wgt = all_wgt[-nb_keep:] 322 323 #final selection of the interesting weight to keep 324 all_wgt.sort() 325 # drop the lowest weight 326 nb_keep = max(20, int(nb_event*trunc_error*10)) 327 all_wgt = all_wgt[-nb_keep:] 328 self.seek(0) 329 return all_wgt, cross, nb_event
330
331 - def write_events(self, event):
332 """ write a single events or a list of event 333 if self.eventgroup is ON, then add <eventgroup> around the lists of events 334 """ 335 if isinstance(event, Event): 336 if self.eventgroup: 337 self.write('<eventgroup>\n%s\n</eventgroup>\n' % event) 338 else: 339 self.write(str(event)) 340 elif isinstance(event, list): 341 if self.eventgroup: 342 self.write('<eventgroup>\n') 343 for evt in event: 344 self.write(str(evt)) 345 if self.eventgroup: 346 self.write('</eventgroup>\n')
347
348 - def unweight(self, outputpath, get_wgt=None, max_wgt=0, trunc_error=0, 349 event_target=0, log_level=logging.INFO, normalization='average'):
350 """unweight the current file according to wgt information wgt. 351 which can either be a fct of the event or a tag in the rwgt list. 352 max_wgt allow to do partial unweighting. 353 trunc_error allow for dynamical partial unweighting 354 event_target reweight for that many event with maximal trunc_error. 355 (stop to write event when target is reached) 356 """ 357 if not get_wgt: 358 def weight(event): 359 return event.wgt
360 get_wgt = weight 361 unwgt_name = "central weight" 362 elif isinstance(get_wgt, str): 363 unwgt_name =get_wgt 364 def get_wgt(event): 365 event.parse_reweight() 366 return event.reweight_data[unwgt_name]
367 else: 368 unwgt_name = get_wgt.func_name 369 370 # check which weight to write 371 if hasattr(self, "written_weight"): 372 written_weight = lambda x: math.copysign(self.written_weight,float(x)) 373 else: 374 written_weight = lambda x: x 375 376 all_wgt, cross, nb_event = self.initialize_unweighting(get_wgt, trunc_error) 377 378 # function that need to be define on the flight 379 def max_wgt_for_trunc(trunc): 380 """find the weight with the maximal truncation.""" 381 382 xsum = 0 383 i=1 384 while (xsum - all_wgt[-i] * (i-1) <= cross['abs'] * trunc): 385 max_wgt = all_wgt[-i] 386 xsum += all_wgt[-i] 387 i +=1 388 if i == len(all_wgt): 389 break 390 391 return max_wgt 392 # end of the function 393 394 # choose the max_weight 395 if not max_wgt: 396 if trunc_error == 0 or len(all_wgt)<2 or event_target: 397 max_wgt = all_wgt[-1] 398 else: 399 max_wgt = max_wgt_for_trunc(trunc_error) 400 401 # need to modify the banner so load it to an object 402 if self.banner: 403 try: 404 import internal 405 except: 406 import madgraph.various.banner as banner_module 407 else: 408 import internal.banner as banner_module 409 if not isinstance(self.banner, banner_module.Banner): 410 banner = self.get_banner() 411 # 1. modify the cross-section 412 banner.modify_init_cross(cross) 413 # 3. add information about change in weight 414 banner["unweight"] = "unweighted by %s" % unwgt_name 415 else: 416 banner = self.banner 417 # modify the lha strategy 418 curr_strategy = banner.get_lha_strategy() 419 if normalization in ['unit', 'sum']: 420 strategy = 3 421 else: 422 strategy = 4 423 if curr_strategy >0: 424 banner.set_lha_strategy(abs(strategy)) 425 else: 426 banner.set_lha_strategy(-1*abs(strategy)) 427 428 # Do the reweighting (up to 20 times if we have target_event) 429 nb_try = 20 430 nb_keep = 0 431 for i in range(nb_try): 432 self.seek(0) 433 if event_target: 434 if i==0: 435 max_wgt = max_wgt_for_trunc(0) 436 else: 437 #guess the correct max_wgt based on last iteration 438 efficiency = nb_keep/nb_event 439 needed_efficiency = event_target/nb_event 440 last_max_wgt = max_wgt 441 needed_max_wgt = last_max_wgt * efficiency / needed_efficiency 442 443 min_max_wgt = max_wgt_for_trunc(trunc_error) 444 max_wgt = max(min_max_wgt, needed_max_wgt) 445 max_wgt = min(max_wgt, all_wgt[-1]) 446 if max_wgt == last_max_wgt: 447 if nb_keep < event_target and log_level>=10: 448 logger.log(log_level+10,"fail to reach target %s", event_target) 449 break 450 else: 451 break 452 453 #create output file (here since we are sure that we have to rewrite it) 454 if outputpath: 455 outfile = EventFile(outputpath, "w") 456 # need to write banner information 457 # need to see what to do with rwgt information! 458 if self.banner and outputpath: 459 banner.write(outfile, close_tag=False) 460 461 # scan the file 462 nb_keep = 0 463 trunc_cross = 0 464 for event in self: 465 r = random.random() 466 wgt = get_wgt(event) 467 if abs(wgt) < r * max_wgt: 468 continue 469 elif wgt > 0: 470 nb_keep += 1 471 event.wgt = written_weight(max(wgt, max_wgt)) 472 if abs(wgt) > max_wgt: 473 trunc_cross += abs(wgt) - max_wgt 474 if event_target ==0 or nb_keep <= event_target: 475 if outputpath: 476 outfile.write(str(event)) 477 478 elif wgt < 0: 479 nb_keep += 1 480 event.wgt = -1* written_weight(max(abs(wgt), max_wgt)) 481 if abs(wgt) > max_wgt: 482 trunc_cross += abs(wgt) - max_wgt 483 if outputpath and (event_target ==0 or nb_keep <= event_target): 484 outfile.write(str(event)) 485 486 if event_target and nb_keep > event_target: 487 if not outputpath: 488 #no outputpath define -> wants only the nb of unweighted events 489 continue 490 elif event_target and i != nb_try-1 and nb_keep >= event_target *1.05: 491 outfile.write("</LesHouchesEvents>\n") 492 outfile.close() 493 #logger.log(log_level, "Found Too much event %s. Try to reduce truncation" % nb_keep) 494 continue 495 else: 496 outfile.write("</LesHouchesEvents>\n") 497 outfile.close() 498 break 499 elif event_target == 0: 500 if outputpath: 501 outfile.write("</LesHouchesEvents>\n") 502 outfile.close() 503 break 504 elif outputpath: 505 outfile.write("</LesHouchesEvents>\n") 506 outfile.close() 507 # logger.log(log_level, "Found only %s event. Reduce max_wgt" % nb_keep) 508 509 else: 510 # pass here if event_target > 0 and all the attempt fail. 511 logger.log(log_level+10,"fail to reach target event %s (iteration=%s)", event_target,i) 512 513 # logger.log(log_level, "Final maximum weight used for final "+\ 514 # "unweighting is %s yielding %s events." % (max_wgt,nb_keep)) 515 516 if event_target: 517 nb_events_unweighted = nb_keep 518 nb_keep = min( event_target, nb_keep) 519 else: 520 nb_events_unweighted = nb_keep 521 522 logger.log(log_level, "write %i event (efficiency %.2g %%, truncation %.2g %%) after %i iteration(s)", 523 nb_keep, nb_events_unweighted/nb_event*100, trunc_cross/cross['abs']*100, i) 524 525 #correct the weight in the file if not the correct number of event 526 if nb_keep != event_target and hasattr(self, "written_weight") and strategy !=4: 527 written_weight = lambda x: math.copysign(self.written_weight*event_target/nb_keep, float(x)) 528 startfile = EventFile(outputpath) 529 tmpname = pjoin(os.path.dirname(outputpath), "wgtcorrected_"+ os.path.basename(outputpath)) 530 outfile = EventFile(tmpname, "w") 531 outfile.write(startfile.banner) 532 for event in startfile: 533 event.wgt = written_weight(event.wgt) 534 outfile.write(str(event)) 535 outfile.write("</LesHouchesEvents>\n") 536 startfile.close() 537 outfile.close() 538 shutil.move(tmpname, outputpath) 539 540 541 542 543 self.max_wgt = max_wgt 544 return nb_keep 545
546 - def apply_fct_on_event(self, *fcts, **opts):
547 """ apply one or more fct on all event. """ 548 549 opt= {"print_step": 5000, "maxevent":float("inf"),'no_output':False} 550 opt.update(opts) 551 start = time.time() 552 nb_fct = len(fcts) 553 out = [] 554 for i in range(nb_fct): 555 out.append([]) 556 self.seek(0) 557 nb_event = 0 558 for event in self: 559 nb_event += 1 560 if opt["print_step"] and (nb_event % opt["print_step"]) == 0: 561 if hasattr(self,"len"): 562 print("currently at %s/%s event [%is]" % (nb_event, self.len, time.time()-start)) 563 else: 564 print("currently at %s event [%is]" % (nb_event, time.time()-start)) 565 for i in range(nb_fct): 566 value = fcts[i](event) 567 if not opt['no_output']: 568 out[i].append(value) 569 if nb_event > opt['maxevent']: 570 break 571 if nb_fct == 1: 572 return out[0] 573 else: 574 return out
575
576 - def split(self, nb_event=0, partition=None, cwd=os.path.curdir, zip=False):
577 """split the file in multiple file. Do not change the weight!""" 578 579 nb_file = -1 580 for i, event in enumerate(self): 581 if (not (partition is None) and i==sum(partition[:nb_file+1])) or \ 582 (partition is None and i % nb_event == 0): 583 if i: 584 #close previous file 585 current.write('</LesHouchesEvent>\n') 586 current.close() 587 # create the new file 588 nb_file +=1 589 # If end of partition then finish writing events here. 590 if not partition is None and (nb_file+1>len(partition)): 591 return nb_file 592 if zip: 593 current = EventFile(pjoin(cwd,'%s_%s.lhe.gz' % (self.name, nb_file)),'w') 594 else: 595 current = open(pjoin(cwd,'%s_%s.lhe' % (self.name, nb_file)),'w') 596 current.write(self.banner) 597 current.write(str(event)) 598 if i!=0: 599 current.write('</LesHouchesEvent>\n') 600 current.close() 601 602 return nb_file +1
603
604 - def update_HwU(self, hwu, fct, name='lhe', keep_wgt=False, maxevents=sys.maxint):
605 """take a HwU and add this event file for the function fct""" 606 607 if not isinstance(hwu, list): 608 hwu = [hwu] 609 610 class HwUUpdater(object): 611 612 def __init__(self, fct, keep_wgt): 613 614 self.fct = fct 615 self.first = True 616 self.keep_wgt = keep_wgt
617 618 def add(self, event): 619 620 value = self.fct(event) 621 # initialise the curve for the first call 622 if self.first: 623 for h in hwu: 624 # register the variables 625 if isinstance(value, dict): 626 h.add_line(value.keys()) 627 else: 628 629 h.add_line(name) 630 if self.keep_wgt is True: 631 event.parse_reweight() 632 h.add_line(['%s_%s' % (name, key) 633 for key in event.reweight_data]) 634 elif self.keep_wgt: 635 h.add_line(self.keep_wgt.values()) 636 self.first = False 637 # Fill the histograms 638 for h in hwu: 639 if isinstance(value, tuple): 640 h.addEvent(value[0], value[1]) 641 else: 642 h.addEvent(value,{name:event.wgt}) 643 if self.keep_wgt: 644 event.parse_reweight() 645 if self.keep_wgt is True: 646 data = dict(('%s_%s' % (name, key),event.reweight_data[key]) 647 for key in event.reweight_data) 648 h.addEvent(value, data) 649 else: 650 data = dict(( value,event.reweight_data[key]) 651 for key,value in self.keep_wgt.items()) 652 h.addEvent(value, data) 653 654 655 656 self.apply_fct_on_event(HwUUpdater(fct,keep_wgt).add, no_output=True,maxevent=maxevents) 657 return hwu 658
659 - def create_syscalc_data(self, out_path, pythia_input=None):
660 """take the lhe file and add the matchscale from the pythia_input file""" 661 662 if pythia_input: 663 def next_data(): 664 for line in open(pythia_input): 665 if line.startswith('#'): 666 continue 667 data = line.split() 668 print (int(data[0]), data[-3], data[-2], data[-1]) 669 yield (int(data[0]), data[-3], data[-2], data[-1])
670 else: 671 def next_data(): 672 i=0 673 while 1: 674 yield [i,0,0,0] 675 i+=1 676 sys_iterator = next_data() 677 #ensure that we are at the beginning of the file 678 self.seek(0) 679 out = open(out_path,'w') 680 681 pdf_pattern = re.compile(r'''<init>(.*)</init>''', re.M+re.S) 682 init = pdf_pattern.findall(self.banner)[0].split('\n',2)[1] 683 id1, id2, _, _, _, _, pdf1,pdf2,_,_ = init.split() 684 id = [int(id1), int(id2)] 685 type = [] 686 for i in range(2): 687 if abs(id[i]) == 2212: 688 if i > 0: 689 type.append(1) 690 else: 691 type.append(-1) 692 else: 693 type.append(0) 694 pdf = max(int(pdf1),int(pdf2)) 695 696 out.write("<header>\n" + \ 697 "<orgpdf>%i</orgpdf>\n" % pdf + \ 698 "<beams> %s %s</beams>\n" % tuple(type) + \ 699 "</header>\n") 700 701 702 nevt, smin, smax, scomp = sys_iterator.next() 703 for i, orig_event in enumerate(self): 704 if i < nevt: 705 continue 706 new_event = Event() 707 sys = orig_event.parse_syscalc_info() 708 new_event.syscalc_data = sys 709 if smin: 710 new_event.syscalc_data['matchscale'] = "%s %s %s" % (smin, scomp, smax) 711 out.write(str(new_event), nevt) 712 try: 713 nevt, smin, smax, scomp = sys_iterator.next() 714 except StopIteration: 715 break 716
717 - def get_alphas(self, scale, lhapdf_config='lhapdf-config'):
718 """return the alphas value associated to a given scale""" 719 720 if hasattr(self, 'alpsrunner'): 721 return self.alpsrunner(scale) 722 723 # 724 banner = banner_mod.Banner(self.banner) 725 run_card = banner.charge_card('run_card') 726 use_runner = False 727 if abs(run_card['lpp1']) != 1 and abs(run_card['lpp2']) != 1: 728 # no pdf use. -> use Runner 729 use_runner = True 730 else: 731 # try to use lhapdf 732 lhapdf = misc.import_python_lhapdf(lhapdf_config) 733 if not lhapdf: 734 logger.warning('fail to link to lhapdf for the alphas-running. Use Two loop computation') 735 use_runner = True 736 try: 737 self.pdf = lhapdf.mkPDF(int(self.banner.run_card.get_lhapdf_id())) 738 except Exception: 739 logger.warning('fail to link to lhapdf for the alphas-running. Use Two loop computation') 740 use_runner = True 741 742 if not use_runner: 743 self.alpsrunner = lambda scale: self.pdf.alphasQ(scale) 744 else: 745 try: 746 from models.model_reader import Alphas_Runner 747 except ImportError: 748 root = os.path.dirname(__file__) 749 root_path = pjoin(root, os.pardir, os.pardir) 750 try: 751 import internal.madevent_interface as me_int 752 cmd = me_int.MadEventCmd(root_path,force_run=True) 753 except ImportError: 754 import internal.amcnlo_run_interface as me_int 755 cmd = me_int.Cmd(root_path,force_run=True) 756 if 'mg5_path' in cmd.options and cmd.options['mg5_path']: 757 sys.path.append(cmd.options['mg5_path']) 758 from models.model_reader import Alphas_Runner 759 760 if not hasattr(banner, 'param_card'): 761 param_card = banner.charge_card('param_card') 762 else: 763 param_card = banner.param_card 764 765 asmz = param_card.get_value('sminputs', 3, 0.13) 766 nloop =2 767 zmass = param_card.get_value('mass', 23, 91.188) 768 cmass = param_card.get_value('mass', 4, 1.4) 769 if cmass == 0: 770 cmass = 1.4 771 bmass = param_card.get_value('mass', 5, 4.7) 772 if bmass == 0: 773 bmass = 4.7 774 self.alpsrunner = Alphas_Runner(asmz, nloop, zmass, cmass, bmass) 775 776 777 778 return self.alpsrunner(scale)
779
780 781 782 783 784 785 786 787 -class EventFileGzip(EventFile, gzip.GzipFile):
788 """A way to read/write a gzipped lhef event"""
789
790 -class EventFileNoGzip(EventFile, file):
791 """A way to read a standard event file""" 792
793 - def close(self,*args, **opts):
794 795 out = super(EventFileNoGzip, self).close(*args, **opts) 796 if self.to_zip: 797 misc.gzip(self.name)
798
799 -class MultiEventFile(EventFile):
800 """a class to read simultaneously multiple file and read them in mixing them. 801 Unweighting can be done at the same time. 802 The number of events in each file need to be provide in advance 803 (if not provide the file is first read to find that number""" 804
805 - def __new__(cls, start_list=[],parse=True):
806 return object.__new__(MultiEventFile)
807
808 - def __init__(self, start_list=[], parse=True):
809 """if trunc_error is define here then this allow 810 to only read all the files twice and not three times.""" 811 self.eventgroup = False 812 self.files = [] 813 self.parsefile = parse #if self.files is formatted or just the path 814 self.banner = '' 815 self.initial_nb_events = [] 816 self.total_event_in_files = 0 817 self.curr_nb_events = [] 818 self.allcross = [] 819 self.error = [] 820 self.across = [] 821 self.scales = [] 822 if start_list: 823 if parse: 824 for p in start_list: 825 self.add(p) 826 else: 827 self.files = start_list 828 self._configure = False
829
830 - def close(self,*args,**opts):
831 for f in self.files: 832 f.close(*args, **opts)
833
834 - def add(self, path, cross, error, across, nb_event=0, scale=1):
835 """ add a file to the pool, across allow to reweight the sum of weight 836 in the file to the given cross-section 837 """ 838 839 if across == 0: 840 # No event linked to this channel -> so no need to include it 841 return 842 843 obj = EventFile(path) 844 obj.eventgroup = self.eventgroup 845 if len(self.files) == 0 and not self.banner: 846 self.banner = obj.banner 847 self.curr_nb_events.append(0) 848 self.initial_nb_events.append(0) 849 self.allcross.append(cross) 850 self.across.append(across) 851 self.error.append(error) 852 self.scales.append(scale) 853 self.files.append(obj) 854 if nb_event: 855 obj.len = nb_event 856 self._configure = False 857 return obj
858
859 - def __iter__(self):
860 return self
861
862 - def next(self):
863 864 if not self._configure: 865 self.configure() 866 867 remaining_event = self.total_event_in_files - sum(self.curr_nb_events) 868 if remaining_event == 0: 869 raise StopIteration 870 # determine which file need to be read 871 nb_event = random.randint(1, remaining_event) 872 sum_nb=0 873 for i, obj in enumerate(self.files): 874 sum_nb += self.initial_nb_events[i] - self.curr_nb_events[i] 875 if nb_event <= sum_nb: 876 self.curr_nb_events[i] += 1 877 event = obj.next() 878 if not self.eventgroup: 879 event.sample_scale = self.scales[i] # for file reweighting 880 else: 881 for evt in event: 882 evt.sample_scale = self.scales[i] 883 return event 884 else: 885 raise Exception
886 887
888 - def define_init_banner(self, wgt, lha_strategy, proc_charac=None):
889 """define the part of the init_banner""" 890 891 if not self.banner: 892 return 893 894 # compute the cross-section of each splitted channel 895 grouped_cross = {} 896 grouped_error = {} 897 for i,ff in enumerate(self.files): 898 filename = ff.name 899 from_init = False 900 Pdir = [P for P in filename.split(os.path.sep) if P.startswith('P')] 901 if Pdir: 902 Pdir = Pdir[-1] 903 group = Pdir.split("_")[0][1:] 904 if not group.isdigit(): 905 from_init = True 906 else: 907 from_init = True 908 909 if not from_init: 910 if group in grouped_cross: 911 grouped_cross[group] += self.allcross[i] 912 grouped_error[group] += self.error[i]**2 913 else: 914 grouped_cross[group] = self.allcross[i] 915 grouped_error[group] = self.error[i]**2 916 else: 917 ban = banner_mod.Banner(ff.banner) 918 for line in ban['init'].split('\n'): 919 splitline = line.split() 920 if len(splitline)==4: 921 cross, error, _, group = splitline 922 if int(group) in grouped_cross: 923 grouped_cross[group] += float(cross) 924 grouped_error[group] += float(error)**2 925 else: 926 grouped_cross[group] = float(cross) 927 grouped_error[group] = float(error)**2 928 nb_group = len(grouped_cross) 929 930 # compute the information for the first line 931 try: 932 run_card = self.banner.run_card 933 except: 934 run_card = self.banner.charge_card("run_card") 935 936 init_information = run_card.get_banner_init_information() 937 #correct for special case 938 if proc_charac and proc_charac['ninitial'] == 1: 939 #special case for 1>N 940 init_information = run_card.get_banner_init_information() 941 event = self.next() 942 init_information["idbmup1"] = event[0].pdg 943 init_information["ebmup1"] = event[0].mass 944 init_information["idbmup2"] = 0 945 init_information["ebmup2"] = 0 946 self.seek(0) 947 else: 948 # check special case without PDF for one (or both) beam 949 if init_information["idbmup1"] == 0: 950 event = self.next() 951 init_information["idbmup1"]= event[0].pdg 952 if init_information["idbmup2"] == 0: 953 init_information["idbmup2"]= event[1].pdg 954 self.seek(0) 955 if init_information["idbmup2"] == 0: 956 event = self.next() 957 init_information["idbmup2"] = event[1].pdg 958 self.seek(0) 959 960 init_information["nprup"] = nb_group 961 962 if run_card["lhe_version"] < 3: 963 init_information["generator_info"] = "" 964 else: 965 init_information["generator_info"] = "<generator name='MadGraph5_aMC@NLO' version='%s'>please cite 1405.0301 </generator>\n" \ 966 % misc.get_pkg_info()['version'] 967 968 # cross_information: 969 cross_info = "%(cross)e %(error)e %(wgt)e %(id)i" 970 init_information["cross_info"] = [] 971 for id in grouped_cross: 972 conv = {"id": int(id), "cross": grouped_cross[id], "error": math.sqrt(grouped_error[id]), 973 "wgt": wgt} 974 init_information["cross_info"].append( cross_info % conv) 975 init_information["cross_info"] = '\n'.join(init_information["cross_info"]) 976 init_information['lha_stra'] = -1 * abs(lha_strategy) 977 978 template_init =\ 979 """ %(idbmup1)i %(idbmup2)i %(ebmup1)e %(ebmup2)e %(pdfgup1)i %(pdfgup2)i %(pdfsup1)i %(pdfsup2)i %(lha_stra)i %(nprup)i 980 %(cross_info)s 981 %(generator_info)s 982 """ 983 984 self.banner["init"] = template_init % init_information
985 986 987
988 - def initialize_unweighting(self, getwgt, trunc_error):
989 """ scan once the file to return 990 - the list of the hightest weight (of size trunc_error*NB_EVENT 991 - the cross-section by type of process 992 - the total number of events in the files 993 In top of that it initialise the information for the next routine 994 to determine how to choose which file to read 995 """ 996 self.seek(0) 997 all_wgt = [] 998 total_event = 0 999 sum_cross = collections.defaultdict(int) 1000 for i,f in enumerate(self.files): 1001 nb_event = 0 1002 # We need to loop over the event file to get some information about the 1003 # new cross-section/ wgt of event. 1004 cross = collections.defaultdict(int) 1005 new_wgt =[] 1006 for event in f: 1007 nb_event += 1 1008 total_event += 1 1009 event.sample_scale = 1 1010 wgt = getwgt(event) 1011 cross['all'] += wgt 1012 cross['abs'] += abs(wgt) 1013 cross[event.ievent] += wgt 1014 new_wgt.append(abs(wgt)) 1015 # avoid all_wgt to be too large 1016 if nb_event % 20000 == 0: 1017 new_wgt.sort() 1018 # drop the lowest weight 1019 nb_keep = max(20, int(nb_event*trunc_error*15)) 1020 new_wgt = new_wgt[-nb_keep:] 1021 if nb_event == 0: 1022 raise Exception 1023 # store the information 1024 self.initial_nb_events[i] = nb_event 1025 self.scales[i] = self.across[i]/cross['abs'] if self.across[i] else 1 1026 #misc.sprint("sum of wgt in event %s is %s. Should be %s => scale %s (nb_event: %s)" 1027 # % (i, cross['all'], self.allcross[i], self.scales[i], nb_event)) 1028 for key in cross: 1029 sum_cross[key] += cross[key]* self.scales[i] 1030 all_wgt +=[self.scales[i] * w for w in new_wgt] 1031 all_wgt.sort() 1032 nb_keep = max(20, int(total_event*trunc_error*10)) 1033 all_wgt = all_wgt[-nb_keep:] 1034 1035 self.total_event_in_files = total_event 1036 #final selection of the interesting weight to keep 1037 all_wgt.sort() 1038 # drop the lowest weight 1039 nb_keep = max(20, int(total_event*trunc_error*10)) 1040 all_wgt = all_wgt[-nb_keep:] 1041 self.seek(0) 1042 self._configure = True 1043 return all_wgt, sum_cross, total_event
1044
1045 - def configure(self):
1046 1047 self._configure = True 1048 for i,f in enumerate(self.files): 1049 self.initial_nb_events[i] = len(f) 1050 self.total_event_in_files = sum(self.initial_nb_events)
1051
1052 - def __len__(self):
1053 1054 return len(self.files)
1055
1056 - def seek(self, pos):
1057 """ """ 1058 1059 if pos !=0: 1060 raise Exception 1061 for i in range(len(self)): 1062 self.curr_nb_events[i] = 0 1063 for f in self.files: 1064 f.seek(pos)
1065
1066 - def unweight(self, outputpath, get_wgt, **opts):
1067 """unweight the current file according to wgt information wgt. 1068 which can either be a fct of the event or a tag in the rwgt list. 1069 max_wgt allow to do partial unweighting. 1070 trunc_error allow for dynamical partial unweighting 1071 event_target reweight for that many event with maximal trunc_error. 1072 (stop to write event when target is reached) 1073 """ 1074 1075 if isinstance(get_wgt, str): 1076 unwgt_name =get_wgt 1077 def get_wgt_multi(event): 1078 event.parse_reweight() 1079 return event.reweight_data[unwgt_name] * event.sample_scale
1080 else: 1081 unwgt_name = get_wgt.func_name 1082 get_wgt_multi = lambda event: get_wgt(event) * event.sample_scale 1083 #define the weighting such that we have built-in the scaling 1084 1085 if 'proc_charac' in opts: 1086 if opts['proc_charac']: 1087 proc_charac = opts['proc_charac'] 1088 else: 1089 proc_charac=None 1090 del opts['proc_charac'] 1091 else: 1092 proc_charac = None 1093 1094 if 'event_target' in opts and opts['event_target']: 1095 if 'normalization' in opts: 1096 if opts['normalization'] == 'sum': 1097 new_wgt = sum(self.across)/opts['event_target'] 1098 strategy = 3 1099 elif opts['normalization'] == 'average': 1100 strategy = 4 1101 new_wgt = sum(self.across) 1102 elif opts['normalization'] == 'unit': 1103 strategy =3 1104 new_wgt = 1. 1105 else: 1106 strategy = 4 1107 new_wgt = sum(self.across) 1108 self.define_init_banner(new_wgt, strategy, proc_charac=proc_charac) 1109 self.written_weight = new_wgt 1110 elif 'write_init' in opts and opts['write_init']: 1111 self.define_init_banner(0,0, proc_charac=proc_charac) 1112 del opts['write_init'] 1113 return super(MultiEventFile, self).unweight(outputpath, get_wgt_multi, **opts)
1114
1115 - def write(self, path, random=False, banner=None, get_info=False):
1116 """ """ 1117 1118 if isinstance(path, str): 1119 out = EventFile(path, 'w') 1120 if self.parsefile and not banner: 1121 banner = self.files[0].banner 1122 elif not banner: 1123 firstlhe = EventFile(self.files[0]) 1124 banner = firstlhe.banner 1125 else: 1126 out = path 1127 if banner: 1128 out.write(banner) 1129 nb_event = 0 1130 info = collections.defaultdict(float) 1131 if random and self.open: 1132 for event in self: 1133 nb_event +=1 1134 out.write(event) 1135 if get_info: 1136 event.parse_reweight() 1137 for key, value in event.reweight_data.items(): 1138 info[key] += value 1139 info['central'] += event.wgt 1140 elif not random: 1141 for i,f in enumerate(self.files): 1142 #check if we need to parse the file or not 1143 if not self.parsefile: 1144 if i==0: 1145 try: 1146 lhe = firstlhe 1147 except: 1148 lhe = EventFile(f) 1149 else: 1150 lhe = EventFile(f) 1151 else: 1152 lhe = f 1153 for event in lhe: 1154 nb_event +=1 1155 if get_info: 1156 event.parse_reweight() 1157 for key, value in event.reweight_data.items(): 1158 info[key] += value 1159 info['central'] += event.wgt 1160 out.write(str(event)) 1161 lhe.close() 1162 out.write("</LesHouchesEvents>\n") 1163 return nb_event, info
1164
1165 - def remove(self):
1166 """ """ 1167 if self.parsefile: 1168 for f in self.files: 1169 os.remove(f.name) 1170 else: 1171 for f in self.files: 1172 os.remove(f)
1173
1174 1175 1176 -class Event(list):
1177 """Class storing a single event information (list of particles + global information)""" 1178 1179 warning_order = True # raise a warning if the order of the particle are not in accordance of child/mother 1180
1181 - def __init__(self, text=None):
1182 """The initialization of an empty Event (or one associate to a text file)""" 1183 list.__init__(self) 1184 1185 # First line information 1186 self.nexternal = 0 1187 self.ievent = 0 1188 self.wgt = 0 1189 self.aqcd = 0 1190 self.scale = 0 1191 self.aqed = 0 1192 self.aqcd = 0 1193 # Weight information 1194 self.tag = '' 1195 self.eventflag = {} # for information in <event > 1196 self.comment = '' 1197 self.reweight_data = {} 1198 self.matched_scale_data = None 1199 self.syscalc_data = {} 1200 if text: 1201 self.parse(text)
1202 1203 1204
1205 - def parse(self, text):
1206 """Take the input file and create the structured information""" 1207 #text = re.sub(r'</?event>', '', text) # remove pointless tag 1208 status = 'first' 1209 for line in text.split('\n'): 1210 line = line.strip() 1211 if not line: 1212 continue 1213 elif line[0] == '#': 1214 self.comment += '%s\n' % line 1215 continue 1216 elif line.startswith('<event'): 1217 if '=' in line: 1218 found = re.findall(r"""(\w*)=(?:(?:['"])([^'"]*)(?=['"])|(\S*))""",line) 1219 #for '<event line=4 value=\'3\' error="5" test=" 1 and 2">\n' 1220 #return [('line', '', '4'), ('value', '3', ''), ('error', '5', ''), ('test', ' 1 and 2', '')] 1221 self.eventflag = dict((n, a1) if a1 else (n,a2) for n,a1,a2 in found) 1222 # return {'test': ' 1 and 2', 'line': '4', 'value': '3', 'error': '5'} 1223 continue 1224 1225 elif 'first' == status: 1226 if '<rwgt>' in line: 1227 status = 'tag' 1228 else: 1229 self.assign_scale_line(line) 1230 status = 'part' 1231 continue 1232 if '<' in line: 1233 status = 'tag' 1234 1235 if 'part' == status: 1236 part = Particle(line, event=self) 1237 if part.E != 0: 1238 self.append(part) 1239 elif self.nexternal: 1240 self.nexternal-=1 1241 else: 1242 if '</event>' in line: 1243 line = line.replace('</event>','',1) 1244 self.tag += '%s\n' % line 1245 1246 self.assign_mother()
1247
1248 - def assign_mother(self):
1249 # assign the mother: 1250 for i,particle in enumerate(self): 1251 if i < particle.mother1 or i < particle.mother2: 1252 if self.warning_order: 1253 logger.warning("Order of particle in the event did not agree with parent/child order. This might be problematic for some code.") 1254 Event.warning_order = False 1255 self.reorder_mother_child() 1256 return self.assign_mother() 1257 1258 if particle.mother1: 1259 try: 1260 particle.mother1 = self[int(particle.mother1) -1] 1261 except Exception: 1262 logger.warning("WRONG MOTHER INFO %s", self) 1263 particle.mother1 = 0 1264 if particle.mother2: 1265 try: 1266 particle.mother2 = self[int(particle.mother2) -1] 1267 except Exception: 1268 logger.warning("WRONG MOTHER INFO %s", self) 1269 particle.mother2 = 0
1270
1271 - def rescale_weights(self, ratio):
1272 """change all the weights by a given ratio""" 1273 1274 self.wgt *= ratio 1275 self.parse_reweight() 1276 for key in self.reweight_data: 1277 self.reweight_data[key] *= ratio 1278 return self.wgt
1279
1280 - def reorder_mother_child(self):
1281 """check and correct the mother/child position. 1282 only correct one order by call (but this is a recursive call)""" 1283 1284 tomove, position = None, None 1285 for i,particle in enumerate(self): 1286 if i < particle.mother1: 1287 # move i after particle.mother1 1288 tomove, position = i, particle.mother1-1 1289 break 1290 if i < particle.mother2: 1291 tomove, position = i, particle.mother2-1 1292 1293 # nothing to change -> we are done 1294 if not tomove: 1295 return 1296 1297 # move the particles: 1298 particle = self.pop(tomove) 1299 self.insert(int(position), particle) 1300 1301 #change the mother id/ event_id in the event. 1302 for i, particle in enumerate(self): 1303 particle.event_id = i 1304 #misc.sprint( i, particle.event_id) 1305 m1, m2 = particle.mother1, particle.mother2 1306 if m1 == tomove +1: 1307 particle.mother1 = position+1 1308 elif tomove < m1 <= position +1: 1309 particle.mother1 -= 1 1310 if m2 == tomove +1: 1311 particle.mother2 = position+1 1312 elif tomove < m2 <= position +1: 1313 particle.mother2 -= 1 1314 # re-call the function for the next potential change 1315 return self.reorder_mother_child()
1316 1317 1318 1319 1320 1321
1322 - def parse_reweight(self):
1323 """Parse the re-weight information in order to return a dictionary 1324 {key: value}. If no group is define group should be '' """ 1325 if self.reweight_data: 1326 return self.reweight_data 1327 self.reweight_data = {} 1328 self.reweight_order = [] 1329 start, stop = self.tag.find('<rwgt>'), self.tag.find('</rwgt>') 1330 if start != -1 != stop : 1331 pattern = re.compile(r'''<\s*wgt id=(?:\'|\")(?P<id>[^\'\"]+)(?:\'|\")\s*>\s*(?P<val>[\ded+-.]*)\s*</wgt>''',re.I) 1332 data = pattern.findall(self.tag[start:stop]) 1333 try: 1334 self.reweight_data = dict([(pid, float(value)) for (pid, value) in data 1335 if not self.reweight_order.append(pid)]) 1336 # the if is to create the order file on the flight 1337 except ValueError, error: 1338 raise Exception, 'Event File has unvalid weight. %s' % error 1339 self.tag = self.tag[:start] + self.tag[stop+7:] 1340 return self.reweight_data
1341
1342 - def parse_nlo_weight(self, real_type=(1,11), threshold=None):
1343 """ """ 1344 if hasattr(self, 'nloweight'): 1345 return self.nloweight 1346 1347 start, stop = self.tag.find('<mgrwgt>'), self.tag.find('</mgrwgt>') 1348 if start != -1 != stop : 1349 1350 text = self.tag[start+8:stop] 1351 self.nloweight = NLO_PARTIALWEIGHT(text, self, real_type=real_type, 1352 threshold=threshold) 1353 return self.nloweight
1354
1355 - def rewrite_nlo_weight(self, wgt=None):
1356 """get the string associate to the weight""" 1357 1358 text="""<mgrwgt> 1359 %(total_wgt).10e %(nb_wgt)i %(nb_event)i 0 1360 %(event)s 1361 %(wgt)s 1362 </mgrwgt>""" 1363 1364 1365 if not wgt: 1366 if not hasattr(self, 'nloweight'): 1367 return 1368 wgt = self.nloweight 1369 1370 data = {'total_wgt': wgt.total_wgt, #need to check name and meaning, 1371 'nb_wgt': wgt.nb_wgt, 1372 'nb_event': wgt.nb_event, 1373 'event': '\n'.join(p.__str__(mode='fortran') for p in wgt.momenta), 1374 'wgt':'\n'.join(w.__str__(mode='formatted') 1375 for e in wgt.cevents for w in e.wgts)} 1376 1377 data['total_wgt'] = sum([w.ref_wgt for e in wgt.cevents for w in e.wgts]) 1378 start, stop = self.tag.find('<mgrwgt>'), self.tag.find('</mgrwgt>') 1379 1380 self.tag = self.tag[:start] + text % data + self.tag[stop+9:]
1381 1382
1383 - def parse_lo_weight(self):
1384 """ """ 1385 1386 1387 if hasattr(self, 'loweight'): 1388 return self.loweight 1389 1390 if not hasattr(Event, 'loweight_pattern'): 1391 Event.loweight_pattern = re.compile('''<rscale>\s*(?P<nqcd>\d+)\s+(?P<ren_scale>[\d.e+-]+)\s*</rscale>\s*\n\s* 1392 <asrwt>\s*(?P<asrwt>[\s\d.+-e]+)\s*</asrwt>\s*\n\s* 1393 <pdfrwt\s+beam=["']?1["']?\>\s*(?P<beam1>[\s\d.e+-]*)\s*</pdfrwt>\s*\n\s* 1394 <pdfrwt\s+beam=["']?2["']?\>\s*(?P<beam2>[\s\d.e+-]*)\s*</pdfrwt>\s*\n\s* 1395 <totfact>\s*(?P<totfact>[\d.e+-]*)\s*</totfact> 1396 ''',re.X+re.I+re.M) 1397 1398 start, stop = self.tag.find('<mgrwt>'), self.tag.find('</mgrwt>') 1399 1400 if start != -1 != stop : 1401 text = self.tag[start+8:stop] 1402 1403 info = Event.loweight_pattern.search(text) 1404 if not info: 1405 raise Exception, '%s not parsed'% text 1406 self.loweight={} 1407 self.loweight['n_qcd'] = int(info.group('nqcd')) 1408 self.loweight['ren_scale'] = float(info.group('ren_scale')) 1409 self.loweight['asrwt'] =[float(x) for x in info.group('asrwt').split()[1:]] 1410 self.loweight['tot_fact'] = float(info.group('totfact')) 1411 1412 args = info.group('beam1').split() 1413 npdf = int(args[0]) 1414 self.loweight['n_pdfrw1'] = npdf 1415 self.loweight['pdf_pdg_code1'] = [int(i) for i in args[1:1+npdf]] 1416 self.loweight['pdf_x1'] = [float(i) for i in args[1+npdf:1+2*npdf]] 1417 self.loweight['pdf_q1'] = [float(i) for i in args[1+2*npdf:1+3*npdf]] 1418 args = info.group('beam2').split() 1419 npdf = int(args[0]) 1420 self.loweight['n_pdfrw2'] = npdf 1421 self.loweight['pdf_pdg_code2'] = [int(i) for i in args[1:1+npdf]] 1422 self.loweight['pdf_x2'] = [float(i) for i in args[1+npdf:1+2*npdf]] 1423 self.loweight['pdf_q2'] = [float(i) for i in args[1+2*npdf:1+3*npdf]] 1424 1425 else: 1426 return None 1427 return self.loweight
1428 1429
1430 - def parse_matching_scale(self):
1431 """Parse the line containing the starting scale for the shower""" 1432 1433 if self.matched_scale_data is not None: 1434 return self.matched_scale_data 1435 1436 self.matched_scale_data = [] 1437 1438 1439 pattern = re.compile("<scales\s|</scales>") 1440 data = re.split(pattern,self.tag) 1441 if len(data) == 1: 1442 return [] 1443 else: 1444 tmp = {} 1445 start,content, end = data 1446 self.tag = "%s%s" % (start, end) 1447 pattern = re.compile("pt_clust_(\d*)=\"([\de+-.]*)\"") 1448 for id,value in pattern.findall(content): 1449 tmp[int(id)] = float(value) 1450 for i in range(1, len(self)+1): 1451 if i in tmp: 1452 self.matched_scale_data.append(tmp[i]) 1453 else: 1454 self.matched_scale_data.append(-1) 1455 return self.matched_scale_data
1456
1457 - def parse_syscalc_info(self):
1458 """ parse the flag for syscalc between <mgrwt></mgrwt> 1459 <mgrwt> 1460 <rscale> 3 0.26552898E+03</rscale> 1461 <asrwt>0</asrwt> 1462 <pdfrwt beam="1"> 1 21 0.14527945E+00 0.26552898E+03</pdfrwt> 1463 <pdfrwt beam="2"> 1 21 0.15249110E-01 0.26552898E+03</pdfrwt> 1464 <totfact> 0.10344054E+04</totfact> 1465 </mgrwt> 1466 """ 1467 if self.syscalc_data: 1468 return self.syscalc_data 1469 1470 pattern = re.compile("<mgrwt>|</mgrwt>") 1471 pattern2 = re.compile("<(?P<tag>[\w]*)(?:\s*(\w*)=[\"'](.*)[\"']\s*|\s*)>(.*)</(?P=tag)>") 1472 data = re.split(pattern,self.tag) 1473 if len(data) == 1: 1474 return [] 1475 else: 1476 tmp = {} 1477 start,content, end = data 1478 self.tag = "%s%s" % (start, end) 1479 for tag, key, keyval, tagval in pattern2.findall(content): 1480 if key: 1481 self.syscalc_data[(tag, key, keyval)] = tagval 1482 else: 1483 self.syscalc_data[tag] = tagval 1484 return self.syscalc_data
1485 1486
1487 - def add_decay_to_particle(self, position, decay_event):
1488 """define the decay of the particle id by the event pass in argument""" 1489 1490 this_particle = self[position] 1491 #change the status to internal particle 1492 this_particle.status = 2 1493 this_particle.helicity = 0 1494 1495 # some usefull information 1496 decay_particle = decay_event[0] 1497 this_4mom = FourMomentum(this_particle) 1498 nb_part = len(self) #original number of particle 1499 1500 thres = decay_particle.E*1e-10 1501 assert max(decay_particle.px, decay_particle.py, decay_particle.pz) < thres,\ 1502 "not on rest particle %s %s %s %s" % (decay_particle.E, decay_particle.px,decay_particle.py,decay_particle.pz) 1503 1504 self.nexternal += decay_event.nexternal -1 1505 old_scales = list(self.parse_matching_scale()) 1506 if old_scales: 1507 jet_position = sum(1 for i in range(position) if self[i].status==1) 1508 initial_pos = sum(1 for i in range(position) if self[i].status==-1) 1509 self.matched_scale_data.pop(initial_pos+jet_position) 1510 # add the particle with only handling the 4-momenta/mother 1511 # color information will be corrected later. 1512 for particle in decay_event[1:]: 1513 # duplicate particle to avoid border effect 1514 new_particle = Particle(particle, self) 1515 new_particle.event_id = len(self) 1516 self.append(new_particle) 1517 if old_scales: 1518 self.matched_scale_data.append(old_scales[initial_pos+jet_position]) 1519 # compute and assign the new four_momenta 1520 new_momentum = this_4mom.boost(FourMomentum(new_particle)) 1521 new_particle.set_momentum(new_momentum) 1522 # compute the new mother 1523 for tag in ['mother1', 'mother2']: 1524 mother = getattr(particle, tag) 1525 if isinstance(mother, Particle): 1526 mother_id = getattr(particle, tag).event_id 1527 if mother_id == 0: 1528 setattr(new_particle, tag, this_particle) 1529 else: 1530 try: 1531 setattr(new_particle, tag, self[nb_part + mother_id -1]) 1532 except Exception, error: 1533 print error 1534 misc.sprint( self) 1535 misc.sprint(nb_part + mother_id -1) 1536 misc.sprint(tag) 1537 misc.sprint(position, decay_event) 1538 misc.sprint(particle) 1539 misc.sprint(len(self), nb_part + mother_id -1) 1540 raise 1541 elif tag == "mother2" and isinstance(particle.mother1, Particle): 1542 new_particle.mother2 = this_particle 1543 else: 1544 raise Exception, "Something weird happens. Please report it for investigation" 1545 # Need to correct the color information of the particle 1546 # first find the first available color index 1547 max_color=501 1548 for particle in self[:nb_part]: 1549 max_color=max(max_color, particle.color1, particle.color2) 1550 1551 # define a color mapping and assign it: 1552 color_mapping = {} 1553 color_mapping[decay_particle.color1] = this_particle.color1 1554 color_mapping[decay_particle.color2] = this_particle.color2 1555 for particle in self[nb_part:]: 1556 if particle.color1: 1557 if particle.color1 not in color_mapping: 1558 max_color +=1 1559 color_mapping[particle.color1] = max_color 1560 particle.color1 = max_color 1561 else: 1562 particle.color1 = color_mapping[particle.color1] 1563 if particle.color2: 1564 if particle.color2 not in color_mapping: 1565 max_color +=1 1566 color_mapping[particle.color2] = max_color 1567 particle.color2 = max_color 1568 else: 1569 particle.color2 = color_mapping[particle.color2]
1570
1571 - def add_decays(self, pdg_to_decay):
1572 """use auto-recursion""" 1573 1574 pdg_to_decay = dict(pdg_to_decay) 1575 1576 for i,particle in enumerate(self): 1577 if particle.status != 1: 1578 continue 1579 if particle.pdg in pdg_to_decay and pdg_to_decay[particle.pdg]: 1580 one_decay = pdg_to_decay[particle.pdg].pop() 1581 self.add_decay_to_particle(i, one_decay) 1582 return self.add_decays(pdg_to_decay) 1583 return self
1584 1585 1586
1587 - def remove_decay(self, pdg_code=0, event_id=None):
1588 1589 to_remove = [] 1590 if event_id is not None: 1591 to_remove.append(self[event_id]) 1592 1593 if pdg_code: 1594 for particle in self: 1595 if particle.pid == pdg_code: 1596 to_remove.append(particle) 1597 1598 new_event = Event() 1599 # copy first line information + ... 1600 for tag in ['nexternal', 'ievent', 'wgt', 'aqcd', 'scale', 'aqed','tag','comment']: 1601 setattr(new_event, tag, getattr(self, tag)) 1602 1603 for particle in self: 1604 if isinstance(particle.mother1, Particle) and particle.mother1 in to_remove: 1605 to_remove.append(particle) 1606 if particle.status == 1: 1607 new_event.nexternal -= 1 1608 continue 1609 elif isinstance(particle.mother2, Particle) and particle.mother2 in to_remove: 1610 to_remove.append(particle) 1611 if particle.status == 1: 1612 new_event.nexternal -= 1 1613 continue 1614 else: 1615 new_event.append(Particle(particle)) 1616 1617 #ensure that the event_id is correct for all_particle 1618 # and put the status to 1 for removed particle 1619 for pos, particle in enumerate(new_event): 1620 particle.event_id = pos 1621 if particle in to_remove: 1622 particle.status = 1 1623 return new_event
1624
1625 - def get_decay(self, pdg_code=0, event_id=None):
1626 1627 to_start = [] 1628 if event_id is not None: 1629 to_start.append(self[event_id]) 1630 1631 elif pdg_code: 1632 for particle in self: 1633 if particle.pid == pdg_code: 1634 to_start.append(particle) 1635 break 1636 1637 new_event = Event() 1638 # copy first line information + ... 1639 for tag in ['ievent', 'wgt', 'aqcd', 'scale', 'aqed','tag','comment']: 1640 setattr(new_event, tag, getattr(self, tag)) 1641 1642 # Add the decaying particle 1643 old2new = {} 1644 new_decay_part = Particle(to_start[0]) 1645 new_decay_part.mother1 = None 1646 new_decay_part.mother2 = None 1647 new_decay_part.status = -1 1648 old2new[new_decay_part.event_id] = len(old2new) 1649 new_event.append(new_decay_part) 1650 1651 1652 # add the other particle 1653 for particle in self: 1654 if isinstance(particle.mother1, Particle) and particle.mother1.event_id in old2new\ 1655 or isinstance(particle.mother2, Particle) and particle.mother2.event_id in old2new: 1656 old2new[particle.event_id] = len(old2new) 1657 new_event.append(Particle(particle)) 1658 1659 #ensure that the event_id is correct for all_particle 1660 # and correct the mother1/mother2 by the new reference 1661 nexternal = 0 1662 for pos, particle in enumerate(new_event): 1663 particle.event_id = pos 1664 if particle.mother1: 1665 particle.mother1 = new_event[old2new[particle.mother1.event_id]] 1666 if particle.mother2: 1667 particle.mother2 = new_event[old2new[particle.mother2.event_id]] 1668 if particle.status in [-1,1]: 1669 nexternal +=1 1670 new_event.nexternal = nexternal 1671 1672 return new_event
1673 1674
1675 - def check(self):
1676 """check various property of the events""" 1677 1678 #1. Check that the 4-momenta are conserved 1679 E, px, py, pz = 0,0,0,0 1680 absE, abspx, abspy, abspz = 0,0,0,0 1681 for particle in self: 1682 coeff = 1 1683 if particle.status == -1: 1684 coeff = -1 1685 elif particle.status != 1: 1686 continue 1687 E += coeff * particle.E 1688 absE += abs(particle.E) 1689 px += coeff * particle.px 1690 py += coeff * particle.py 1691 pz += coeff * particle.pz 1692 abspx += abs(particle.px) 1693 abspy += abs(particle.py) 1694 abspz += abs(particle.pz) 1695 # check that relative error is under control 1696 threshold = 5e-7 1697 if E/absE > threshold: 1698 logger.critical(self) 1699 raise Exception, "Do not conserve Energy %s, %s" % (E/absE, E) 1700 if px/abspx > threshold: 1701 logger.critical(self) 1702 raise Exception, "Do not conserve Px %s, %s" % (px/abspx, px) 1703 if py/abspy > threshold: 1704 logger.critical(self) 1705 raise Exception, "Do not conserve Py %s, %s" % (py/abspy, py) 1706 if pz/abspz > threshold: 1707 logger.critical(self) 1708 raise Exception, "Do not conserve Pz %s, %s" % (pz/abspz, pz) 1709 1710 #2. check the color of the event 1711 self.check_color_structure()
1712
1713 - def assign_scale_line(self, line):
1714 """read the line corresponding to global event line 1715 format of the line is: 1716 Nexternal IEVENT WEIGHT SCALE AEW AS 1717 """ 1718 inputs = line.split() 1719 assert len(inputs) == 6 1720 self.nexternal=int(inputs[0]) 1721 self.ievent=int(inputs[1]) 1722 self.wgt=float(inputs[2]) 1723 self.scale=float(inputs[3]) 1724 self.aqed=float(inputs[4]) 1725 self.aqcd=float(inputs[5])
1726
1727 - def get_tag_and_order(self):
1728 """Return the unique tag identifying the SubProcesses for the generation. 1729 Usefull for program like MadSpin and Reweight module.""" 1730 1731 initial, final, order = [], [], [[], []] 1732 for particle in self: 1733 if particle.status == -1: 1734 initial.append(particle.pid) 1735 order[0].append(particle.pid) 1736 elif particle.status == 1: 1737 final.append(particle.pid) 1738 order[1].append(particle.pid) 1739 initial.sort(), final.sort() 1740 tag = (tuple(initial), tuple(final)) 1741 return tag, order
1742 1743 @staticmethod
1744 - def mass_shuffle(momenta, sqrts, new_mass, new_sqrts=None):
1745 """use the RAMBO method to shuffle the PS. initial sqrts is preserved.""" 1746 1747 if not new_sqrts: 1748 new_sqrts = sqrts 1749 1750 oldm = [p.mass_sqr for p in momenta] 1751 newm = [m**2 for m in new_mass] 1752 tot_mom = sum(momenta, FourMomentum()) 1753 if tot_mom.pt2 > 1e-5: 1754 boost_back = FourMomentum(tot_mom.mass,0,0,0).boost_to_restframe(tot_mom) 1755 for i,m in enumerate(momenta): 1756 momenta[i] = m.boost_to_restframe(tot_mom) 1757 1758 # this is the equation 4.3 of RAMBO paper 1759 f = lambda chi: new_sqrts - sum(math.sqrt(max(0, M + chi**2*(p.E**2-m))) 1760 for M,p,m in zip(newm, momenta,oldm)) 1761 # this is the derivation of the function 1762 df = lambda chi: -1* sum(chi*(p.E**2-m)/math.sqrt(max(0,(p.E**2-m)*chi**2+M)) 1763 for M,p,m in zip(newm, momenta,oldm)) 1764 1765 if sum(new_mass) > new_sqrts: 1766 return momenta, 0 1767 try: 1768 chi = misc.newtonmethod(f, df, 1.0, error=1e-7,maxiter=1000) 1769 except: 1770 return momenta, 0 1771 # create the new set of momenta # eq. (4.2) 1772 new_momenta = [] 1773 for i,p in enumerate(momenta): 1774 new_momenta.append( 1775 FourMomentum(math.sqrt(newm[i]+chi**2*(p.E**2-oldm[i])), 1776 chi*p.px, chi*p.py, chi*p.pz)) 1777 1778 #if __debug__: 1779 # for i,p in enumerate(new_momenta): 1780 # misc.sprint(p.mass_sqr, new_mass[i]**2, i,p, momenta[i]) 1781 # assert p.mass_sqr == new_mass[i]**2 1782 1783 # compute the jacobian factor (eq. 4.9) 1784 jac = chi**(3*len(momenta)-3) 1785 jac *= reduce(operator.mul,[p.E/k.E for p,k in zip(momenta, new_momenta)],1) 1786 jac *= sum(p.norm_sq/p.E for p in momenta) 1787 jac /= sum(k.norm_sq/k.E for k in new_momenta) 1788 1789 # boost back the events in the lab-frame 1790 if tot_mom.pt2 > 1e-5: 1791 for i,m in enumerate(new_momenta): 1792 new_momenta[i] = m.boost_to_restframe(boost_back) 1793 return new_momenta, jac
1794 1795 1796 1797
1798 - def change_ext_mass(self, new_param_card):
1799 """routine to rescale the mass via RAMBO method. no internal mass preserve. 1800 sqrts is preserve (RAMBO algo) 1801 """ 1802 1803 old_momenta = [] 1804 new_masses = [] 1805 change_mass = False # check if we need to change the mass 1806 for part in self: 1807 if part.status == 1: 1808 old_momenta.append(FourMomentum(part)) 1809 new_masses.append(new_param_card.get_value('mass', abs(part.pid))) 1810 if not misc.equal(part.mass, new_masses[-1], 4, zero_limit=10): 1811 change_mass = True 1812 1813 if not change_mass: 1814 return 1 1815 1816 sqrts = self.sqrts 1817 1818 # apply the RAMBO algo 1819 new_mom, jac = self.mass_shuffle(old_momenta, sqrts, new_masses) 1820 1821 #modify the momenta of the particles: 1822 ind =0 1823 for part in self: 1824 if part.status==1: 1825 part.E, part.px, part.py, part.pz, part.mass = \ 1826 new_mom[ind].E, new_mom[ind].px, new_mom[ind].py, new_mom[ind].pz,new_mom[ind].mass 1827 ind+=1 1828 return jac
1829
1830 - def change_sqrts(self, new_sqrts):
1831 """routine to rescale the momenta to change the invariant mass""" 1832 1833 old_momenta = [] 1834 incoming = [] 1835 masses = [] 1836 for part in self: 1837 if part.status == -1: 1838 incoming.append(FourMomentum(part)) 1839 if part.status == 1: 1840 old_momenta.append(FourMomentum(part)) 1841 masses.append(part.mass) 1842 1843 p_init = FourMomentum() 1844 p_inits = [] 1845 n_init = 0 1846 for p in incoming: 1847 n_init +=1 1848 p_init += p 1849 p_inits.append(p) 1850 old_sqrts = p_init.mass 1851 1852 new_mom, jac = self.mass_shuffle(old_momenta, old_sqrts, masses, new_sqrts=new_sqrts) 1853 1854 #modify the momenta of the particles: 1855 ind =0 1856 for part in self: 1857 if part.status==1: 1858 part.E, part.px, part.py, part.pz, part.mass = \ 1859 new_mom[ind].E, new_mom[ind].px, new_mom[ind].py, new_mom[ind].pz,new_mom[ind].mass 1860 ind+=1 1861 1862 #change the initial state 1863 p_init = FourMomentum() 1864 for part in self: 1865 if part.status==1: 1866 p_init += part 1867 if n_init == 1: 1868 for part in self: 1869 if part.status == -1: 1870 part.E, part.px, part.py, part.pz = \ 1871 p_init.E, p_init.px, p_init.py, p_init.pz 1872 elif n_init ==2: 1873 if not misc.equal(p_init.px, 0) or not misc.equal(p_init.py, 0): 1874 raise Exception 1875 if not misc.equal(p_inits[0].px, 0) or not misc.equal(p_inits[0].py, 0): 1876 raise Exception 1877 #assume that initial energy is written as 1878 # p1 = (sqrts/2*exp(eta), 0, 0 , E1) 1879 # p2 = (sqrts/2*exp(-eta), 0, 0 , -E2) 1880 # keep eta fix 1881 eta = math.log(2*p_inits[0].E/old_sqrts) 1882 new_p = [[new_sqrts/2*math.exp(eta), 0., 0., new_sqrts/2*math.exp(eta)], 1883 [new_sqrts/2*math.exp(-eta), 0., 0., -new_sqrts/2*math.exp(-eta)]] 1884 1885 ind=0 1886 for part in self: 1887 if part.status == -1: 1888 part.E, part.px, part.py, part.pz = new_p[ind] 1889 ind+=1 1890 if ind ==2: 1891 break 1892 else: 1893 raise Exception 1894 1895 return jac
1896 1897
1898 - def get_helicity(self, get_order, allow_reversed=True):
1899 """return a list with the helicities in the order asked for""" 1900 1901 #avoid to modify the input 1902 order = [list(get_order[0]), list(get_order[1])] 1903 out = [9] *(len(order[0])+len(order[1])) 1904 for i, part in enumerate(self): 1905 if part.status == 1: #final 1906 try: 1907 ind = order[1].index(part.pid) 1908 except ValueError, error: 1909 if not allow_reversed: 1910 raise error 1911 else: 1912 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 1913 try: 1914 return self.get_helicity(order, False) 1915 except ValueError: 1916 raise error 1917 position = len(order[0]) + ind 1918 order[1][ind] = 0 1919 elif part.status == -1: 1920 try: 1921 ind = order[0].index(part.pid) 1922 except ValueError, error: 1923 if not allow_reversed: 1924 raise error 1925 else: 1926 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 1927 try: 1928 return self.get_helicity(order, False) 1929 except ValueError: 1930 raise error 1931 1932 position = ind 1933 order[0][ind] = 0 1934 else: #intermediate 1935 continue 1936 out[position] = int(part.helicity) 1937 return out
1938 1939
1940 - def check_color_structure(self):
1941 """check the validity of the color structure""" 1942 1943 #1. check that each color is raised only once. 1944 color_index = collections.defaultdict(int) 1945 for particle in self: 1946 if particle.status in [-1,1]: 1947 if particle.color1: 1948 color_index[particle.color1] +=1 1949 if -7 < particle.pdg < 0: 1950 raise Exception, "anti-quark with color tag" 1951 if particle.color2: 1952 color_index[particle.color2] +=1 1953 if 7 > particle.pdg > 0: 1954 raise Exception, "quark with anti-color tag" 1955 1956 1957 for key,value in color_index.items(): 1958 if value > 2: 1959 print self 1960 print key, value 1961 raise Exception, 'Wrong color_flow' 1962 1963 1964 #2. check that each parent present have coherent color-structure 1965 check = [] 1966 popup_index = [] #check that the popup index are created in a unique way 1967 for particle in self: 1968 mothers = [] 1969 childs = [] 1970 if particle.mother1: 1971 mothers.append(particle.mother1) 1972 if particle.mother2 and particle.mother2 is not particle.mother1: 1973 mothers.append(particle.mother2) 1974 if not mothers: 1975 continue 1976 if (particle.mother1.event_id, particle.mother2.event_id) in check: 1977 continue 1978 check.append((particle.mother1.event_id, particle.mother2.event_id)) 1979 1980 childs = [p for p in self if p.mother1 is particle.mother1 and \ 1981 p.mother2 is particle.mother2] 1982 1983 mcolors = [] 1984 manticolors = [] 1985 for m in mothers: 1986 if m.color1: 1987 if m.color1 in manticolors: 1988 manticolors.remove(m.color1) 1989 else: 1990 mcolors.append(m.color1) 1991 if m.color2: 1992 if m.color2 in mcolors: 1993 mcolors.remove(m.color2) 1994 else: 1995 manticolors.append(m.color2) 1996 ccolors = [] 1997 canticolors = [] 1998 for m in childs: 1999 if m.color1: 2000 if m.color1 in canticolors: 2001 canticolors.remove(m.color1) 2002 else: 2003 ccolors.append(m.color1) 2004 if m.color2: 2005 if m.color2 in ccolors: 2006 ccolors.remove(m.color2) 2007 else: 2008 canticolors.append(m.color2) 2009 for index in mcolors[:]: 2010 if index in ccolors: 2011 mcolors.remove(index) 2012 ccolors.remove(index) 2013 for index in manticolors[:]: 2014 if index in canticolors: 2015 manticolors.remove(index) 2016 canticolors.remove(index) 2017 2018 if mcolors != []: 2019 #only case is a epsilon_ijk structure. 2020 if len(canticolors) + len(mcolors) != 3: 2021 logger.critical(str(self)) 2022 raise Exception, "Wrong color flow for %s -> %s" ([m.pid for m in mothers], [c.pid for c in childs]) 2023 else: 2024 popup_index += canticolors 2025 elif manticolors != []: 2026 #only case is a epsilon_ijk structure. 2027 if len(ccolors) + len(manticolors) != 3: 2028 logger.critical(str(self)) 2029 raise Exception, "Wrong color flow for %s -> %s" ([m.pid for m in mothers], [c.pid for c in childs]) 2030 else: 2031 popup_index += ccolors 2032 2033 # Check that color popup (from epsilon_ijk) are raised only once 2034 if len(popup_index) != len(set(popup_index)): 2035 logger.critical(self) 2036 raise Exception, "Wrong color flow: identical poping-up index, %s" % (popup_index)
2037
2038 - def __eq__(self, other):
2039 """two event are the same if they have the same momentum. other info are ignored""" 2040 2041 if other is None: 2042 return False 2043 2044 for i,p in enumerate(self): 2045 if p.E != other[i].E: 2046 return False 2047 elif p.pz != other[i].pz: 2048 return False 2049 elif p.px != other[i].px: 2050 return False 2051 elif p.py != other[i].py: 2052 return False 2053 return True
2054 2055
2056 - def __str__(self, event_id=''):
2057 """return a correctly formatted LHE event""" 2058 2059 out="""<event%(event_flag)s> 2060 %(scale)s 2061 %(particles)s 2062 %(comments)s 2063 %(tag)s 2064 %(reweight)s 2065 </event> 2066 """ 2067 if event_id not in ['', None]: 2068 self.eventflag['event'] = str(event_id) 2069 2070 if self.eventflag: 2071 event_flag = ' %s' % ' '.join('%s="%s"' % (k,v) for (k,v) in self.eventflag.items()) 2072 else: 2073 event_flag = '' 2074 2075 if self.nexternal: 2076 scale_str = "%2d %6d %+13.7e %14.8e %14.8e %14.8e" % \ 2077 (self.nexternal,self.ievent,self.wgt,self.scale,self.aqed,self.aqcd) 2078 else: 2079 scale_str = '' 2080 2081 if self.reweight_data: 2082 # check that all key have an order if not add them at the end 2083 if set(self.reweight_data.keys()) != set(self.reweight_order): 2084 self.reweight_order += [k for k in self.reweight_data.keys() \ 2085 if k not in self.reweight_order] 2086 2087 reweight_str = '<rwgt>\n%s\n</rwgt>' % '\n'.join( 2088 '<wgt id=\'%s\'> %+13.7e </wgt>' % (i, float(self.reweight_data[i])) 2089 for i in self.reweight_order if i in self.reweight_data) 2090 else: 2091 reweight_str = '' 2092 2093 tag_str = self.tag 2094 if hasattr(self, 'nloweight') and self.nloweight.modified: 2095 self.rewrite_nlo_weight() 2096 tag_str = self.tag 2097 2098 if self.matched_scale_data: 2099 tmp_scale = ' '.join(['pt_clust_%i=\"%s\"' % (i+1,v) 2100 for i,v in enumerate(self.matched_scale_data) 2101 if v!=-1]) 2102 if tmp_scale: 2103 tag_str = "<scales %s></scales>%s" % (tmp_scale, self.tag) 2104 2105 if self.syscalc_data: 2106 keys= ['rscale', 'asrwt', ('pdfrwt', 'beam', '1'), ('pdfrwt', 'beam', '2'), 2107 'matchscale', 'totfact'] 2108 sys_str = "<mgrwt>\n" 2109 template = """<%(key)s%(opts)s>%(values)s</%(key)s>\n""" 2110 for k in keys: 2111 if k not in self.syscalc_data: 2112 continue 2113 replace = {} 2114 replace['values'] = self.syscalc_data[k] 2115 if isinstance(k, str): 2116 replace['key'] = k 2117 replace['opts'] = '' 2118 else: 2119 replace['key'] = k[0] 2120 replace['opts'] = ' %s=\"%s\"' % (k[1],k[2]) 2121 sys_str += template % replace 2122 sys_str += "</mgrwt>\n" 2123 reweight_str = sys_str + reweight_str 2124 2125 out = out % {'event_flag': event_flag, 2126 'scale': scale_str, 2127 'particles': '\n'.join([str(p) for p in self]), 2128 'tag': tag_str, 2129 'comments': self.comment, 2130 'reweight': reweight_str} 2131 2132 return re.sub('[\n]+', '\n', out)
2133
2134 - def get_momenta(self, get_order, allow_reversed=True):
2135 """return the momenta vector in the order asked for""" 2136 2137 #avoid to modify the input 2138 order = [list(get_order[0]), list(get_order[1])] 2139 out = [''] *(len(order[0])+len(order[1])) 2140 for i, part in enumerate(self): 2141 if part.status == 1: #final 2142 try: 2143 ind = order[1].index(part.pid) 2144 except ValueError, error: 2145 if not allow_reversed: 2146 raise error 2147 else: 2148 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 2149 try: 2150 return self.get_momenta_str(order, False) 2151 except ValueError: 2152 raise error 2153 position = len(order[0]) + ind 2154 order[1][ind] = 0 2155 elif part.status == -1: 2156 try: 2157 ind = order[0].index(part.pid) 2158 except ValueError, error: 2159 if not allow_reversed: 2160 raise error 2161 else: 2162 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 2163 try: 2164 return self.get_momenta_str(order, False) 2165 except ValueError: 2166 raise error 2167 2168 position = ind 2169 order[0][ind] = 0 2170 else: #intermediate 2171 continue 2172 2173 out[position] = (part.E, part.px, part.py, part.pz) 2174 2175 return out
2176 2177
2178 - def get_scale(self,type):
2179 2180 if type == 1: 2181 return self.get_et_scale() 2182 elif type == 2: 2183 return self.get_ht_scale() 2184 elif type == 3: 2185 return self.get_ht_scale(prefactor=0.5) 2186 elif type == 4: 2187 return self.get_sqrts_scale() 2188 elif type == -1: 2189 return self.get_ht_scale(prefactor=0.5)
2190 2191
2192 - def get_ht_scale(self, prefactor=1):
2193 2194 scale = 0 2195 for particle in self: 2196 if particle.status != 1: 2197 continue 2198 p=FourMomentum(particle) 2199 scale += math.sqrt(p.mass_sqr + p.pt**2) 2200 2201 return prefactor * scale
2202 2203
2204 - def get_et_scale(self, prefactor=1):
2205 2206 scale = 0 2207 for particle in self: 2208 if particle.status != 1: 2209 continue 2210 p = FourMomentum(particle) 2211 pt = p.pt 2212 if (pt>0): 2213 scale += p.E*pt/math.sqrt(pt**2+p.pz**2) 2214 2215 return prefactor * scale
2216 2217 @property
2218 - def sqrts(self):
2219 return self.get_sqrts_scale(1)
2220
2221 - def get_sqrts_scale(self, prefactor=1):
2222 2223 scale = 0 2224 init = [] 2225 for particle in self: 2226 if particle.status == -1: 2227 init.append(FourMomentum(particle)) 2228 if len(init) == 1: 2229 return init[0].mass 2230 elif len(init)==2: 2231 return math.sqrt((init[0]+init[1])**2)
2232 2233 2234 2235
2236 - def get_momenta_str(self, get_order, allow_reversed=True):
2237 """return the momenta str in the order asked for""" 2238 2239 out = self.get_momenta(get_order, allow_reversed) 2240 #format 2241 format = '%.12f' 2242 format_line = ' '.join([format]*4) + ' \n' 2243 out = [format_line % one for one in out] 2244 out = ''.join(out).replace('e','d') 2245 return out
2246
2247 -class WeightFile(EventFile):
2248 """A class to allow to read both gzip and not gzip file. 2249 containing only weight from pythia --generated by SysCalc""" 2250
2251 - def __new__(self, path, mode='r', *args, **opt):
2252 if path.endswith(".gz"): 2253 try: 2254 return gzip.GzipFile.__new__(WeightFileGzip, path, mode, *args, **opt) 2255 except IOError, error: 2256 raise 2257 except Exception, error: 2258 if mode == 'r': 2259 misc.gunzip(path) 2260 return file.__new__(WeightFileNoGzip, path[:-3], mode, *args, **opt) 2261 else: 2262 return file.__new__(WeightFileNoGzip, path, mode, *args, **opt)
2263 2264
2265 - def __init__(self, path, mode='r', *args, **opt):
2266 """open file and read the banner [if in read mode]""" 2267 2268 super(EventFile, self).__init__(path, mode, *args, **opt) 2269 self.banner = '' 2270 if mode == 'r': 2271 line = '' 2272 while '</header>' not in line.lower(): 2273 try: 2274 line = super(EventFile, self).next() 2275 except StopIteration: 2276 self.seek(0) 2277 self.banner = '' 2278 break 2279 if "<event" in line.lower(): 2280 self.seek(0) 2281 self.banner = '' 2282 break 2283 2284 self.banner += line
2285
2286 2287 -class WeightFileGzip(WeightFile, EventFileGzip):
2288 pass
2289
2290 -class WeightFileNoGzip(WeightFile, EventFileNoGzip):
2291 pass
2292
2293 2294 -class FourMomentum(object):
2295 """a convenient object for 4-momenta operation""" 2296
2297 - def __init__(self, obj=0, px=0, py=0, pz=0, E=0):
2298 """initialize the four momenta""" 2299 2300 if obj is 0 and E: 2301 obj = E 2302 2303 if isinstance(obj, (FourMomentum, Particle)): 2304 px = obj.px 2305 py = obj.py 2306 pz = obj.pz 2307 E = obj.E 2308 elif isinstance(obj, (list, tuple)): 2309 assert len(obj) ==4 2310 E = obj[0] 2311 px = obj[1] 2312 py = obj[2] 2313 pz = obj[3] 2314 elif isinstance(obj, str): 2315 obj = [float(i) for i in obj.split()] 2316 assert len(obj) ==4 2317 E = obj[0] 2318 px = obj[1] 2319 py = obj[2] 2320 pz = obj[3] 2321 else: 2322 E =obj 2323 2324 2325 self.E = float(E) 2326 self.px = float(px) 2327 self.py = float(py) 2328 self.pz = float(pz)
2329 2330 @property
2331 - def mass(self):
2332 """return the mass""" 2333 return math.sqrt(max(self.E**2 - self.px**2 - self.py**2 - self.pz**2,0))
2334 2335 @property
2336 - def mass_sqr(self):
2337 """return the mass square""" 2338 return self.E**2 - self.px**2 - self.py**2 - self.pz**2
2339 2340 @property
2341 - def pt(self):
2342 return math.sqrt(max(0, self.pt2))
2343 2344 @property
2345 - def pseudorapidity(self):
2346 norm = math.sqrt(self.px**2 + self.py**2+self.pz**2) 2347 return 0.5* math.log((norm - self.pz) / (norm + self.pz))
2348 2349 @property
2350 - def rapidity(self):
2351 return 0.5* math.log((self.E +self.pz) / (self.E - self.pz))
2352 2353 2354 @property
2355 - def pt2(self):
2356 """ return the pt square """ 2357 2358 return self.px**2 + self.py**2
2359 2360 @property
2361 - def norm(self):
2362 """ return |\vec p| """ 2363 return math.sqrt(self.px**2 + self.py**2 + self.pz**2)
2364 2365 @property
2366 - def norm_sq(self):
2367 """ return |\vec p|^2 """ 2368 return self.px**2 + self.py**2 + self.pz**2
2369 2370
2371 - def __add__(self, obj):
2372 2373 assert isinstance(obj, FourMomentum) 2374 new = FourMomentum(self.E+obj.E, 2375 self.px + obj.px, 2376 self.py + obj.py, 2377 self.pz + obj.pz) 2378 return new
2379
2380 - def __iadd__(self, obj):
2381 """update the object with the sum""" 2382 self.E += obj.E 2383 self.px += obj.px 2384 self.py += obj.py 2385 self.pz += obj.pz 2386 return self
2387
2388 - def __sub__(self, obj):
2389 2390 assert isinstance(obj, FourMomentum) 2391 new = FourMomentum(self.E-obj.E, 2392 self.px - obj.px, 2393 self.py - obj.py, 2394 self.pz - obj.pz) 2395 return new
2396
2397 - def __isub__(self, obj):
2398 """update the object with the sum""" 2399 self.E -= obj.E 2400 self.px -= obj.px 2401 self.py -= obj.py 2402 self.pz -= obj.pz 2403 return self
2404
2405 - def __mul__(self, obj):
2406 if isinstance(obj, FourMomentum): 2407 return self.E*obj.E - self.px *obj.px - self.py * obj.py - self.pz * obj.pz 2408 elif isinstance(obj, (float, int)): 2409 return FourMomentum(obj*self.E,obj*self.px,obj*self.py,obj*self.pz ) 2410 else: 2411 raise NotImplemented
2412 __rmul__ = __mul__ 2413
2414 - def __pow__(self, power):
2415 assert power in [1,2] 2416 2417 if power == 1: 2418 return FourMomentum(self) 2419 elif power == 2: 2420 return self.mass_sqr
2421
2422 - def __repr__(self):
2423 return 'FourMomentum(%s,%s,%s,%s)' % (self.E, self.px, self.py,self.pz)
2424
2425 - def __str__(self, mode='python'):
2426 if mode == 'python': 2427 return self.__repr__() 2428 elif mode == 'fortran': 2429 return '%.10e %.10e %.10e %.10e' % self.get_tuple()
2430
2431 - def get_tuple(self):
2432 return (self.E, self.px, self.py,self.pz)
2433
2434 - def boost(self, mom):
2435 """mom 4-momenta is suppose to be given in the rest frame of this 4-momenta. 2436 the output is the 4-momenta in the frame of this 4-momenta 2437 function copied from HELAS routine.""" 2438 2439 2440 pt = self.px**2 + self.py**2 + self.pz**2 2441 if pt: 2442 s3product = self.px * mom.px + self.py * mom.py + self.pz * mom.pz 2443 mass = self.mass 2444 lf = (mom.E + (self.E - mass) * s3product / pt ) / mass 2445 return FourMomentum(E=(self.E*mom.E+s3product)/mass, 2446 px=mom.px + self.px * lf, 2447 py=mom.py + self.py * lf, 2448 pz=mom.pz + self.pz * lf) 2449 else: 2450 return FourMomentum(mom)
2451
2452 - def zboost(self, pboost=None, E=0, pz=0):
2453 """Both momenta should be in the same frame. 2454 The boost perform correspond to the boost required to set pboost at 2455 rest (only z boost applied). 2456 """ 2457 if isinstance(pboost, FourMomentum): 2458 E = pboost.E 2459 pz = pboost.pz 2460 2461 #beta = pz/E 2462 gamma = E / math.sqrt(E**2-pz**2) 2463 gammabeta = pz / math.sqrt(E**2-pz**2) 2464 2465 out = FourMomentum([gamma*self.E - gammabeta*self.pz, 2466 self.px, 2467 self.py, 2468 gamma*self.pz - gammabeta*self.E]) 2469 2470 if abs(out.pz) < 1e-6 * out.E: 2471 out.pz = 0 2472 return out
2473
2474 - def boost_to_restframe(self, pboost):
2475 """apply the boost transformation such that pboost is at rest in the new frame. 2476 First apply a rotation to allign the pboost to the z axis and then use 2477 zboost routine (see above) 2478 """ 2479 2480 if pboost.px == 0 == pboost.py: 2481 out = self.zboost(E=pboost.E,pz=pboost.pz) 2482 return out 2483 2484 2485 # write pboost as (E, p cosT sinF, p sinT sinF, p cosF) 2486 # rotation such that it become (E, 0 , 0 , p ) is 2487 # cosT sinF , -sinT , cosT sinF 2488 # sinT cosF , cosT , sinT sinF 2489 # -sinT , 0 , cosF 2490 p = math.sqrt( pboost.px**2 + pboost.py**2+ pboost.pz**2) 2491 cosF = pboost.pz / p 2492 sinF = math.sqrt(1-cosF**2) 2493 sinT = pboost.py/p/sinF 2494 cosT = pboost.px/p/sinF 2495 2496 out=FourMomentum([self.E, 2497 self.px*cosT*cosF + self.py*sinT*cosF-self.pz*sinF, 2498 -self.px*sinT+ self.py*cosT, 2499 self.px*cosT*sinF + self.py*sinT*sinF + self.pz*cosF 2500 ]) 2501 out = out.zboost(E=pboost.E,pz=p) 2502 return out
2503
2504 2505 2506 2507 -class OneNLOWeight(object):
2508
2509 - def __init__(self, input, real_type=(1,11)):
2510 """ """ 2511 2512 self.real_type = real_type 2513 if isinstance(input, str): 2514 self.parse(input)
2515
2516 - def __str__(self, mode='display'):
2517 2518 if mode == 'display': 2519 out = """ pwgt: %(pwgt)s 2520 born, real : %(born)s %(real)s 2521 pdgs : %(pdgs)s 2522 bjks : %(bjks)s 2523 scales**2, gs: %(scales2)s %(gs)s 2524 born/real related : %(born_related)s %(real_related)s 2525 type / nfks : %(type)s %(nfks)s 2526 to merge : %(to_merge_pdg)s in %(merge_new_pdg)s 2527 ref_wgt : %(ref_wgt)s""" % self.__dict__ 2528 return out 2529 elif mode == 'formatted': 2530 format_var = [] 2531 variable = [] 2532 2533 def to_add_full(f, v, format_var, variable): 2534 """ function to add to the formatted output""" 2535 if isinstance(v, list): 2536 format_var += [f]*len(v) 2537 variable += v 2538 else: 2539 format_var.append(f) 2540 variable.append(v)
2541 to_add = lambda x,y: to_add_full(x,y, format_var, variable) 2542 #set the formatting 2543 to_add('%.10e', [p*self.bias_wgt for p in self.pwgt]) 2544 to_add('%.10e', self.born) 2545 to_add('%.10e', self.real) 2546 to_add('%i', self.nexternal) 2547 to_add('%i', self.pdgs) 2548 to_add('%i', self.qcdpower) 2549 to_add('%.10e', self.bjks) 2550 to_add('%.10e', self.scales2) 2551 to_add('%.10e', self.gs) 2552 to_add('%i', [self.born_related, self.real_related]) 2553 to_add('%i' , [self.type, self.nfks]) 2554 to_add('%i' , self.to_merge_pdg) 2555 to_add('%i', self.merge_new_pdg) 2556 to_add('%.10e', self.ref_wgt*self.bias_wgt) 2557 to_add('%.10e', self.bias_wgt) 2558 return ' '.join(format_var) % tuple(variable)
2559 2560
2561 - def parse(self, text, keep_bias=False):
2562 """parse the line and create the related object. 2563 keep bias allow to not systematically correct for the bias in the written information""" 2564 #0.546601845792D+00 0.000000000000D+00 0.000000000000D+00 0.119210435309D+02 0.000000000000D+00 5 -1 2 -11 12 21 0 0.24546101D-01 0.15706890D-02 0.12586055D+04 0.12586055D+04 0.12586055D+04 1 2 2 2 5 2 2 0.539995789976D+04 2565 #0.274922677249D+01 0.000000000000D+00 0.000000000000D+00 0.770516514633D+01 0.113763730192D+00 5 21 2 -11 12 1 2 0.52500539D-02 0.30205908D+00 0.45444066D+04 0.45444066D+04 0.45444066D+04 0.12520062D+01 1 2 1 3 5 1 -1 0.110944218997D+05 2566 # below comment are from Rik description email 2567 data = text.split() 2568 # 1. The first three doubles are, as before, the 'wgt', i.e., the overall event of this 2569 # contribution, and the ones multiplying the log[mu_R/QES] and the log[mu_F/QES] 2570 # stripped of alpha_s and the PDFs. 2571 # from example: 0.274922677249D+01 0.000000000000D+00 0.000000000000D+00 2572 self.pwgt = [float(f) for f in data[:3]] 2573 # 2. The next two doubles are the values of the (corresponding) Born and 2574 # real-emission matrix elements. You can either use these values to check 2575 # that the newly computed original matrix element weights are correct, 2576 # or directly use these so that you don't have to recompute the original weights. 2577 # For contributions for which the real-emission matrix elements were 2578 # not computed, the 2nd of these numbers is zero. The opposite is not true, 2579 # because each real-emission phase-space configuration has an underlying Born one 2580 # (this is not unique, but on our code we made a specific choice here). 2581 # This latter information is useful if the real-emission matrix elements 2582 # are unstable; you can then reweight with the Born instead. 2583 # (see also point 9 below, where the momentum configurations are assigned). 2584 # I don't think this instability is real problem when reweighting the real-emission 2585 # with tree-level matrix elements (as we generally would do), but is important 2586 # when reweighting with loop-squared contributions as we have been doing for gg->H. 2587 # (I'm not sure that reweighting tree-level with loop^2 is something that 2588 # we can do in general, because we don't really know what to do with the 2589 # virtual matrix elements because we cannot generate 2-loop diagrams.) 2590 # from example: 0.770516514633D+01 0.113763730192D+00 2591 self.born = float(data[3]) 2592 self.real = float(data[4]) 2593 # 3. integer: number of external particles of the real-emission configuration (as before) 2594 # from example: 5 2595 self.nexternal = int(data[5]) 2596 # 4. PDG codes corresponding to the real-emission configuration (as before) 2597 # from example: 21 2 -11 12 1 2 2598 self.pdgs = [int(i) for i in data[6:6+self.nexternal]] 2599 flag = 6+self.nexternal # new starting point for the position 2600 # 5. next integer is the power of g_strong in the matrix elements (as before) 2601 # from example: 2 2602 self.qcdpower = int(data[flag]) 2603 # 6. 2 doubles: The bjorken x's used for this contribution (as before) 2604 # from example: 0.52500539D-02 0.30205908D+00 2605 self.bjks = [float(f) for f in data[flag+1:flag+3]] 2606 # 7. 3 doubles: The Ellis-sexton scale, the renormalisation scale and the factorisation scale, all squared, used for this contribution (as before) 2607 # from example: 0.45444066D+04 0.45444066D+04 0.45444066D+04 2608 self.scales2 = [float(f) for f in data[flag+3:flag+6]] 2609 # 8.the value of g_strong 2610 # from example: 0.12520062D+01 2611 self.gs = float(data[flag+6]) 2612 # 9. 2 integers: the corresponding Born and real-emission type kinematics. (in the list of momenta) 2613 # Note that also the Born-kinematics has n+1 particles, with, in general, 2614 # one particle with zero momentum (this is not ALWAYS the case, 2615 # there could also be 2 particles with perfectly collinear momentum). 2616 # To convert this from n+1 to a n particles, you have to sum the momenta 2617 # of the two particles that 'merge', see point 12 below. 2618 # from example: 1 2 2619 self.born_related = int(data[flag+7]) 2620 self.real_related = int(data[flag+8]) 2621 # 10. 1 integer: the 'type'. This is the information you should use to determine 2622 # if to reweight with Born, virtual or real-emission matrix elements. 2623 # (Apart from the possible problems with complicated real-emission matrix elements 2624 # that need to be computed very close to the soft/collinear limits, see point 2 above. 2625 # I guess that for tree-level this is always okay, but when reweighting 2626 # a tree-level contribution with a one-loop squared one, as we do 2627 # for gg->Higgs, this is important). 2628 # type=1 : real-emission: 2629 # type=2 : Born: 2630 # type=3 : integrated counter terms: 2631 # type=4 : soft counter-term: 2632 # type=5 : collinear counter-term: 2633 # type=6 : soft-collinear counter-term: 2634 # type=7 : O(alphaS) expansion of Sudakov factor for NNLL+NLO: 2635 # type=8 : soft counter-term (with n+1-body kin.): 2636 # type=9 : collinear counter-term (with n+1-body kin.): 2637 # type=10: soft-collinear counter-term (with n+1-body kin.): 2638 # type=11: real-emission (with n-body kin.): 2639 # type=12: MC subtraction with n-body kin.: 2640 # type=13: MC subtraction with n+1-body kin.: 2641 # type=14: virtual corrections minus approximate virtual 2642 # type=15: approximate virtual corrections: 2643 # from example: 1 2644 self.type = int(data[flag+9]) 2645 # 11. 1 integer: The FKS configuration for this contribution (not really 2646 # relevant for anything, but is used in checking the reweighting to 2647 # get scale & PDF uncertainties). 2648 # from example: 3 2649 self.nfks = int(data[flag+10]) 2650 # 12. 2 integers: the two particles that should be merged to form the 2651 # born contribution from the real-emission one. Remove these two particles 2652 # from the (ordered) list of PDG codes, and insert a newly created particle 2653 # at the location of the minimum of the two particles removed. 2654 # I.e., if you merge particles 2 and 4, you have to insert the new particle 2655 # as the 2nd particle. And particle 5 and above will be shifted down by one. 2656 # from example: 5 1 2657 self.to_merge_pdg = [int (f) for f in data[flag+11:flag+13]] 2658 # 13. 1 integer: the PDG code of the particle that is created after merging the two particles at point 12. 2659 # from example -1 2660 self.merge_new_pdg = int(data[flag+13]) 2661 # 14. 1 double: the reference number that one should be able to reconstruct 2662 # form the weights (point 1 above) and the rest of the information of this line. 2663 # This is really the contribution to this event as computed by the code 2664 # (and is passed to the integrator). It contains everything. 2665 # from example: 0.110944218997D+05 2666 self.ref_wgt = float(data[flag+14]) 2667 # 15. The bias weight. This weight is included in the self.ref_wgt, as well as in 2668 # the self.pwgt. However, it is already removed from the XWGTUP (and 2669 # scale/pdf weights). That means that in practice this weight is not used. 2670 try: 2671 self.bias_wgt = float(data[flag+15]) 2672 except IndexError: 2673 self.bias_wgt = 1.0 2674 2675 if not keep_bias: 2676 self.ref_wgt /= self.bias_wgt 2677 self.pwgt = [p/self.bias_wgt for p in self.pwgt] 2678 2679 #check the momenta configuration linked to the event 2680 if self.type in self.real_type: 2681 self.momenta_config = self.real_related 2682 else: 2683 self.momenta_config = self.born_related
2684
2685 2686 -class NLO_PARTIALWEIGHT(object):
2687
2688 - class BasicEvent(list):
2689 2690
2691 - def __init__(self, momenta, wgts, event, real_type=(1,11)):
2692 2693 list.__init__(self, momenta) 2694 assert self 2695 self.soft = False 2696 self.wgts = wgts 2697 self.pdgs = list(wgts[0].pdgs) 2698 self.event = event 2699 self.real_type = real_type 2700 2701 if wgts[0].momenta_config == wgts[0].born_related: 2702 # need to remove one momenta. 2703 ind1, ind2 = [ind-1 for ind in wgts[0].to_merge_pdg] 2704 if ind1> ind2: 2705 ind1, ind2 = ind2, ind1 2706 if ind1 >= sum(1 for p in event if p.status==-1): 2707 new_p = self[ind1] + self[ind2] 2708 else: 2709 new_p = self[ind1] - self[ind2] 2710 self.pop(ind1) 2711 self.insert(ind1, new_p) 2712 self.pop(ind2) 2713 self.pdgs.pop(ind1) 2714 self.pdgs.insert(ind1, wgts[0].merge_new_pdg ) 2715 self.pdgs.pop(ind2) 2716 # DO NOT update the pdgs of the partial weight! 2717 2718 elif any(w.type in self.real_type for w in wgts): 2719 if any(w.type not in self.real_type for w in wgts): 2720 raise Exception 2721 # Do nothing !!! 2722 # previously (commented we were checking here if the particle 2723 # were too soft this is done later now 2724 # The comment line below allow to convert this event 2725 # to a born one (old method) 2726 # self.pop(ind1) 2727 # self.insert(ind1, new_p) 2728 # self.pop(ind2) 2729 # self.pdgs.pop(ind1) 2730 # self.pdgs.insert(ind1, wgts[0].merge_new_pdg ) 2731 # self.pdgs.pop(ind2) 2732 # # DO NOT update the pdgs of the partial weight! 2733 else: 2734 raise Exception
2735
2736 - def check_fks_singularity(self, ind1, ind2, nb_init=2, threshold=None):
2737 """check that the propagator associated to ij is not too light 2738 [related to soft-collinear singularity]""" 2739 2740 if threshold is None: 2741 threshold = 1e-8 2742 2743 if ind1> ind2: 2744 ind1, ind2 = ind2, ind1 2745 if ind1 >= nb_init: 2746 new_p = self[ind1] + self[ind2] 2747 else: 2748 new_p = self[ind1] - self[ind2] 2749 2750 inv_mass = new_p.mass_sqr 2751 if nb_init == 2: 2752 shat = (self[0]+self[1]).mass_sqr 2753 else: 2754 shat = self[0].mass_sqr 2755 2756 2757 if (abs(inv_mass)/shat < threshold): 2758 return True 2759 else: 2760 return False
2761 2762
2763 - def get_pdg_code(self):
2764 return self.pdgs
2765
2766 - def get_tag_and_order(self):
2767 """ return the tag and order for this basic event""" 2768 (initial, _), _ = self.event.get_tag_and_order() 2769 order = self.get_pdg_code() 2770 2771 2772 initial, out = order[:len(initial)], order[len(initial):] 2773 initial.sort() 2774 out.sort() 2775 return (tuple(initial), tuple(out)), order
2776
2777 - def get_momenta(self, get_order, allow_reversed=True):
2778 """return the momenta vector in the order asked for""" 2779 2780 #avoid to modify the input 2781 order = [list(get_order[0]), list(get_order[1])] 2782 out = [''] *(len(order[0])+len(order[1])) 2783 pdgs = self.get_pdg_code() 2784 for pos, part in enumerate(self): 2785 if pos < len(get_order[0]): #initial 2786 try: 2787 ind = order[0].index(pdgs[pos]) 2788 except ValueError, error: 2789 if not allow_reversed: 2790 raise error 2791 else: 2792 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 2793 try: 2794 return self.get_momenta(order, False) 2795 except ValueError: 2796 raise error 2797 2798 2799 position = ind 2800 order[0][ind] = 0 2801 else: #final 2802 try: 2803 ind = order[1].index(pdgs[pos]) 2804 except ValueError, error: 2805 if not allow_reversed: 2806 raise error 2807 else: 2808 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 2809 try: 2810 return self.get_momenta(order, False) 2811 except ValueError: 2812 raise error 2813 position = len(order[0]) + ind 2814 order[1][ind] = 0 2815 2816 out[position] = (part.E, part.px, part.py, part.pz) 2817 2818 return out
2819 2820
2821 - def get_helicity(self, *args):
2822 return [9] * len(self)
2823 2824 @property
2825 - def aqcd(self):
2826 return self.event.aqcd
2827
2828 - def get_ht_scale(self, prefactor=1):
2829 2830 scale = 0 2831 for particle in self: 2832 p = particle 2833 scale += math.sqrt(max(0, p.mass_sqr + p.pt**2)) 2834 2835 return prefactor * scale
2836
2837 - def get_et_scale(self, prefactor=1):
2838 2839 scale = 0 2840 for particle in self: 2841 p = particle 2842 pt = p.pt 2843 if (pt>0): 2844 scale += p.E*pt/math.sqrt(pt**2+p.pz**2) 2845 2846 return prefactor * scale
2847 2848
2849 - def get_sqrts_scale(self, event,prefactor=1):
2850 2851 scale = 0 2852 nb_init = 0 2853 for particle in event: 2854 if particle.status == -1: 2855 nb_init+=1 2856 if nb_init == 1: 2857 return self[0].mass 2858 elif nb_init==2: 2859 return math.sqrt((self[0]+self[1])**2)
2860 2861 2862 2863
2864 - def __init__(self, input, event, real_type=(1,11), threshold=None):
2865 2866 self.real_type = real_type 2867 self.event = event 2868 self.total_wgt = 0. 2869 self.nb_event = 0 2870 self.nb_wgts = 0 2871 self.threshold = threshold 2872 self.modified = False #set on True if we decide to change internal infor 2873 # that need to be written in the event file. 2874 #need to be set manually when this is the case 2875 if isinstance(input, str): 2876 self.parse(input)
2877 2878 2879
2880 - def parse(self, text):
2881 """create the object from the string information (see example below)""" 2882 #0.2344688900d+00 8 2 0 2883 #0.4676614699d+02 0.0000000000d+00 0.0000000000d+00 0.4676614699d+02 2884 #0.4676614699d+02 0.0000000000d+00 0.0000000000d+00 -.4676614699d+02 2885 #0.4676614699d+02 0.2256794794d+02 0.4332148227d+01 0.4073073437d+02 2886 #0.4676614699d+02 -.2256794794d+02 -.4332148227d+01 -.4073073437d+02 2887 #0.0000000000d+00 -.0000000000d+00 -.0000000000d+00 -.0000000000d+00 2888 #0.4780341163d+02 0.0000000000d+00 0.0000000000d+00 0.4780341163d+02 2889 #0.4822581633d+02 0.0000000000d+00 0.0000000000d+00 -.4822581633d+02 2890 #0.4729127470d+02 0.2347155377d+02 0.5153455534d+01 0.4073073437d+02 2891 #0.4627255267d+02 -.2167412893d+02 -.3519736379d+01 -.4073073437d+02 2892 #0.2465400591d+01 -.1797424844d+01 -.1633719155d+01 -.4224046944d+00 2893 #0.473706252575d-01 0.000000000000d+00 0.000000000000d+00 5 -3 3 -11 11 21 0 0.11849903d-02 0.43683926d-01 0.52807978d+03 0.52807978d+03 0.52807978d+03 1 2 1 0.106660059627d+03 2894 #-.101626389492d-02 0.000000000000d+00 -.181915673961d-03 5 -3 3 -11 11 21 2 0.11849903d-02 0.43683926d-01 0.52807978d+03 0.52807978d+03 0.52807978d+03 1 3 1 -.433615206719d+01 2895 #0.219583436285d-02 0.000000000000d+00 0.000000000000d+00 5 -3 3 -11 11 21 2 0.11849903d-02 0.43683926d-01 0.52807978d+03 0.52807978d+03 0.52807978d+03 1 15 1 0.936909375537d+01 2896 #0.290043597283d-03 0.000000000000d+00 0.000000000000d+00 5 -3 3 -11 11 21 2 0.12292838d-02 0.43683926d-01 0.58606724d+03 0.58606724d+03 0.58606724d+03 1 12 1 0.118841547979d+01 2897 #-.856330613460d-01 0.000000000000d+00 0.000000000000d+00 5 -3 3 -11 11 21 2 0.11849903d-02 0.43683926d-01 0.52807978d+03 0.52807978d+03 0.52807978d+03 1 4 1 -.365375546483d+03 2898 #0.854918237609d-01 0.000000000000d+00 0.000000000000d+00 5 -3 3 -11 11 21 2 0.12112732d-02 0.45047393d-01 0.58606724d+03 0.58606724d+03 0.58606724d+03 2 11 1 0.337816057347d+03 2899 #0.359257891118d-05 0.000000000000d+00 0.000000000000d+00 5 21 3 -11 11 3 2 0.12292838d-02 0.43683926d-01 0.58606724d+03 0.58606724d+03 0.58606724d+03 1 12 3 0.334254554762d+00 2900 #0.929944817736d-03 0.000000000000d+00 0.000000000000d+00 5 21 3 -11 11 3 2 0.12112732d-02 0.45047393d-01 0.58606724d+03 0.58606724d+03 0.58606724d+03 2 11 3 0.835109616010d+02 2901 2902 2903 text = text.lower().replace('d','e') 2904 all_line = text.split('\n') 2905 #get global information 2906 first_line ='' 2907 while not first_line.strip(): 2908 first_line = all_line.pop(0) 2909 2910 wgt, nb_wgt, nb_event, _ = first_line.split() 2911 self.total_wgt = float(wgt.replace('d','e')) 2912 nb_wgt, nb_event = int(nb_wgt), int(nb_event) 2913 self.nb_wgt, self.nb_event = nb_wgt, nb_event 2914 2915 momenta = [] 2916 self.momenta = momenta #keep the original list of momenta to be able to rewrite the events 2917 wgts = [] 2918 for line in all_line: 2919 data = line.split() 2920 if len(data) == 4: 2921 p = FourMomentum(data) 2922 momenta.append(p) 2923 elif len(data)>0: 2924 wgt = OneNLOWeight(line, real_type=self.real_type) 2925 wgts.append(wgt) 2926 2927 assert len(wgts) == int(nb_wgt) 2928 2929 get_weights_for_momenta = dict( (i,[]) for i in range(1,nb_event+1) ) 2930 size_momenta = 0 2931 for wgt in wgts: 2932 if wgt.momenta_config in get_weights_for_momenta: 2933 get_weights_for_momenta[wgt.momenta_config].append(wgt) 2934 else: 2935 if size_momenta == 0: size_momenta = wgt.nexternal 2936 assert size_momenta == wgt.nexternal 2937 get_weights_for_momenta[wgt.momenta_config] = [wgt] 2938 2939 assert sum(len(c) for c in get_weights_for_momenta.values()) == int(nb_wgt), "%s != %s" % (sum(len(c) for c in get_weights_for_momenta.values()), nb_wgt) 2940 2941 # check singular behavior 2942 for key in range(1, nb_event+1): 2943 wgts = get_weights_for_momenta[key] 2944 if not wgts: 2945 continue 2946 if size_momenta == 0: size_momenta = wgts[0].nexternal 2947 p = momenta[size_momenta*(key-1):key*size_momenta] 2948 evt = self.BasicEvent(p, wgts, self.event, self.real_type) 2949 if len(evt) == size_momenta: #real type 2950 for wgt in wgts: 2951 if not wgt.type in self.real_type: 2952 continue 2953 if evt.check_fks_singularity(wgt.to_merge_pdg[0]-1, 2954 wgt.to_merge_pdg[1]-1, 2955 nb_init=sum(1 for p in self.event if p.status==-1), 2956 threshold=self.threshold): 2957 get_weights_for_momenta[wgt.momenta_config].remove(wgt) 2958 get_weights_for_momenta[wgt.born_related].append(wgt) 2959 wgt.momenta_config = wgt.born_related 2960 2961 assert sum(len(c) for c in get_weights_for_momenta.values()) == int(nb_wgt), "%s != %s" % (sum(len(c) for c in get_weights_for_momenta.values()), nb_wgt) 2962 2963 self.cevents = [] 2964 for key in range(1, nb_event+1): 2965 if key in get_weights_for_momenta: 2966 wgt = get_weights_for_momenta[key] 2967 if not wgt: 2968 continue 2969 pdg_to_event = {} 2970 for w in wgt: 2971 pdgs = w.pdgs 2972 if w.momenta_config == w.born_related: 2973 pdgs = list(pdgs) 2974 ind1, ind2 = [ind-1 for ind in w.to_merge_pdg] 2975 if ind1> ind2: 2976 ind1, ind2 = ind2, ind1 2977 pdgs.pop(ind1) 2978 pdgs.insert(ind1, w.merge_new_pdg ) 2979 pdgs.pop(ind2) 2980 pdgs = tuple(pdgs) 2981 if pdgs not in pdg_to_event: 2982 p = momenta[size_momenta*(key-1):key*size_momenta] 2983 evt = self.BasicEvent(p, [w], self.event, self.real_type) 2984 self.cevents.append(evt) 2985 pdg_to_event[pdgs] = evt 2986 else: 2987 pdg_to_event[pdgs].wgts.append(w) 2988 2989 if __debug__: 2990 nb_wgt_check = 0 2991 for cevt in self.cevents: 2992 nb_wgt_check += len(cevt.wgts) 2993 assert nb_wgt_check == int(nb_wgt)
2994 2995 2996 2997 if '__main__' == __name__: 2998 2999 if False: 3000 lhe = EventFile('unweighted_events.lhe.gz') 3001 #lhe.parsing = False 3002 start = time.time() 3003 for event in lhe: 3004 event.parse_lo_weight() 3005 print 'old method -> ', time.time()-start 3006 lhe = EventFile('unweighted_events.lhe.gz') 3007 #lhe.parsing = False 3008 start = time.time() 3009 for event in lhe: 3010 event.parse_lo_weight_test() 3011 print 'new method -> ', time.time()-start 3012 3013 3014 # Example 1: adding some missing information to the event (here distance travelled) 3015 if False: 3016 start = time 3017 lhe = EventFile('unweighted_events.lhe.gz') 3018 output = open('output_events.lhe', 'w') 3019 #write the banner to the output file 3020 output.write(lhe.banner) 3021 # Loop over all events 3022 for event in lhe: 3023 for particle in event: 3024 # modify particle attribute: here remove the mass 3025 particle.mass = 0 3026 particle.vtim = 2 # The one associate to distance travelled by the particle. 3027 3028 #write this modify event 3029 output.write(str(event)) 3030 output.write('</LesHouchesEvent>\n') 3031 3032 # Example 3: Plotting some variable 3033 if False: 3034 lhe = EventFile('unweighted_events.lhe.gz') 3035 import matplotlib.pyplot as plt 3036 import matplotlib.gridspec as gridspec 3037 nbins = 100 3038 3039 nb_pass = 0 3040 data = [] 3041 for event in lhe: 3042 etaabs = 0 3043 etafinal = 0 3044 for particle in event: 3045 if particle.status==1: 3046 p = FourMomentum(particle) 3047 eta = p.pseudorapidity 3048 if abs(eta) > etaabs: 3049 etafinal = eta 3050 etaabs = abs(eta) 3051 if etaabs < 4: 3052 data.append(etafinal) 3053 nb_pass +=1 3054 3055 3056 print nb_pass 3057 gs1 = gridspec.GridSpec(2, 1, height_ratios=[5,1]) 3058 gs1.update(wspace=0, hspace=0) # set the spacing between axes. 3059 ax = plt.subplot(gs1[0]) 3060 3061 n, bins, patches = ax.hist(data, nbins, histtype='step', label='original') 3062 ax_c = ax.twinx() 3063 ax_c.set_ylabel('MadGraph5_aMC@NLO') 3064 ax_c.yaxis.set_label_coords(1.01, 0.25) 3065 ax_c.set_yticks(ax.get_yticks()) 3066 ax_c.set_yticklabels([]) 3067 ax.set_xlim([-4,4]) 3068 print "bin value:", n 3069 print "start/end point of bins", bins 3070 plt.axis('on') 3071 plt.xlabel('weight ratio') 3072 plt.show() 3073 3074 3075 # Example 4: More complex plotting example (with ratio plot) 3076 if False: 3077 lhe = EventFile('unweighted_events.lhe') 3078 import matplotlib.pyplot as plt 3079 import matplotlib.gridspec as gridspec 3080 nbins = 100 3081 3082 #mtau, wtau = 45, 5.1785e-06 3083 mtau, wtau = 1.777, 4.027000e-13 3084 nb_pass = 0 3085 data, data2, data3 = [], [], [] 3086 for event in lhe: 3087 nb_pass +=1 3088 if nb_pass > 10000: 3089 break 3090 tau1 = FourMomentum() 3091 tau2 = FourMomentum() 3092 for part in event: 3093 if part.pid in [-12,11,16]: 3094 momenta = FourMomentum(part) 3095 tau1 += momenta 3096 elif part.pid == 15: 3097 tau2 += FourMomentum(part) 3098 3099 if abs((mtau-tau2.mass())/wtau)<1e6 and tau2.mass() >1: 3100 data.append((tau1.mass()-mtau)/wtau) 3101 data2.append((tau2.mass()-mtau)/wtau) 3102 gs1 = gridspec.GridSpec(2, 1, height_ratios=[5,1]) 3103 gs1.update(wspace=0, hspace=0) # set the spacing between axes. 3104 ax = plt.subplot(gs1[0]) 3105 3106 n, bins, patches = ax.hist(data2, nbins, histtype='step', label='original') 3107 n2, bins2, patches2 = ax.hist(data, bins=bins, histtype='step',label='reconstructed') 3108 import cmath 3109 3110 breit = lambda m : math.sqrt(4*math.pi)*1/(((m)**2-mtau**2)**2+(mtau*wtau)**2)*wtau 3111 3112 data3 = [breit(mtau + x*wtau)*wtau*16867622.6624*50 for x in bins] 3113 3114 ax.plot(bins, data3,label='breit-wigner') 3115 # add the legend 3116 ax.legend() 3117 # add on the right program tag 3118 ax_c = ax.twinx() 3119 ax_c.set_ylabel('MadGraph5_aMC@NLO') 3120 ax_c.yaxis.set_label_coords(1.01, 0.25) 3121 ax_c.set_yticks(ax.get_yticks()) 3122 ax_c.set_yticklabels([]) 3123 3124 plt.title('invariant mass of tau LHE/reconstructed') 3125 plt.axis('on') 3126 ax.set_xticklabels([]) 3127 # ratio plot 3128 ax = plt.subplot(gs1[1]) 3129 data4 = [n[i]/(data3[i]) for i in range(nbins)] 3130 ax.plot(bins, data4 + [0] , 'b') 3131 data4 = [n2[i]/(data3[i]) for i in range(nbins)] 3132 ax.plot(bins, data4 + [0] , 'g') 3133 ax.set_ylim([0,2]) 3134 #remove last y tick to avoid overlap with above plot: 3135 tick = ax.get_yticks() 3136 ax.set_yticks(tick[:-1]) 3137 3138 3139 plt.axis('on') 3140 plt.xlabel('(M - Mtau)/Wtau') 3141 plt.show() 3142