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