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']: 776 splitline = line.split() 777 if len(splitline)==4: 778 cross, error, wgt, group = splitline 779 grouped_cross[int(group)] += cross 780 grouped_error[int(group)] += error**2 781 782 783 nb_group = len(grouped_cross) 784 785 # compute the information for the first line 786 try: 787 run_card = self.banner.run_card 788 except: 789 run_card = self.banner.charge_card("run_card") 790 791 init_information = run_card.get_banner_init_information() 792 if init_information["idbmup1"] == 0: 793 event = self.next() 794 init_information["idbmup1"]= event[0].pdg 795 if init_information["idbmup2"] == 0: 796 init_information["idbmup2"]= event[1].pdg 797 self.seek(0) 798 if init_information["idbmup2"] == 0: 799 event = self.next() 800 init_information["idbmup2"] = event[1].pdg 801 self.seek(0) 802 803 init_information["nprup"] = nb_group 804 805 if run_card["lhe_version"] < 3: 806 init_information["generator_info"] = "" 807 else: 808 init_information["generator_info"] = "<generator name='MadGraph5_aMC@NLO' version='%s'>please cite 1405.0301 </generator>\n" \ 809 % misc.get_pkg_info()['version'] 810 811 # cross_information: 812 cross_info = "%(cross)e %(error)e %(wgt)e %(id)i" 813 init_information["cross_info"] = [] 814 for id in grouped_cross: 815 conv = {"id": int(id), "cross": grouped_cross[id], "error": math.sqrt(grouped_error[id]), 816 "wgt": wgt} 817 init_information["cross_info"].append( cross_info % conv) 818 init_information["cross_info"] = '\n'.join(init_information["cross_info"]) 819 init_information['lha_stra'] = -1 * abs(lha_strategy) 820 821 822 template_init =\ 823 """ %(idbmup1)i %(idbmup2)i %(ebmup1)e %(ebmup2)e %(pdfgup1)i %(pdfgup2)i %(pdfsup1)i %(pdfsup2)i %(lha_stra)i %(nprup)i 824 %(cross_info)s 825 %(generator_info)s 826 """ 827 828 self.banner["init"] = template_init % init_information
829 830 831
832 - def initialize_unweighting(self, getwgt, trunc_error):
833 """ scan once the file to return 834 - the list of the hightest weight (of size trunc_error*NB_EVENT 835 - the cross-section by type of process 836 - the total number of events in the files 837 In top of that it initialise the information for the next routine 838 to determine how to choose which file to read 839 """ 840 self.seek(0) 841 all_wgt = [] 842 total_event = 0 843 sum_cross = collections.defaultdict(int) 844 for i,f in enumerate(self.files): 845 nb_event = 0 846 # We need to loop over the event file to get some information about the 847 # new cross-section/ wgt of event. 848 cross = collections.defaultdict(int) 849 new_wgt =[] 850 for event in f: 851 nb_event += 1 852 total_event += 1 853 event.sample_scale = 1 854 wgt = getwgt(event) 855 cross['all'] += wgt 856 cross['abs'] += abs(wgt) 857 cross[event.ievent] += wgt 858 new_wgt.append(abs(wgt)) 859 # avoid all_wgt to be too large 860 if nb_event % 20000 == 0: 861 new_wgt.sort() 862 # drop the lowest weight 863 nb_keep = max(20, int(nb_event*trunc_error*15)) 864 new_wgt = new_wgt[-nb_keep:] 865 if nb_event == 0: 866 raise Exception 867 # store the information 868 self.initial_nb_events[i] = nb_event 869 self.scales[i] = self.across[i]/cross['abs'] if self.across[i] else 1 870 #misc.sprint("sum of wgt in event %s is %s. Should be %s => scale %s (nb_event: %s)" 871 # % (i, cross['all'], self.allcross[i], self.scales[i], nb_event)) 872 for key in cross: 873 sum_cross[key] += cross[key]* self.scales[i] 874 all_wgt +=[self.scales[i] * w for w in new_wgt] 875 all_wgt.sort() 876 nb_keep = max(20, int(total_event*trunc_error*10)) 877 all_wgt = all_wgt[-nb_keep:] 878 879 self.total_event_in_files = total_event 880 #final selection of the interesting weight to keep 881 all_wgt.sort() 882 # drop the lowest weight 883 nb_keep = max(20, int(total_event*trunc_error*10)) 884 all_wgt = all_wgt[-nb_keep:] 885 self.seek(0) 886 self._configure = True 887 return all_wgt, sum_cross, total_event
888
889 - def configure(self):
890 891 self._configure = True 892 for i,f in enumerate(self.files): 893 self.initial_nb_events[i] = len(f) 894 self.total_event_in_files = sum(self.initial_nb_events)
895
896 - def __len__(self):
897 898 return len(self.files)
899
900 - def seek(self, pos):
901 """ """ 902 903 if pos !=0: 904 raise Exception 905 for i in range(len(self)): 906 self.curr_nb_events[i] = 0 907 for f in self.files: 908 f.seek(pos)
909
910 - def unweight(self, outputpath, get_wgt, **opts):
911 """unweight the current file according to wgt information wgt. 912 which can either be a fct of the event or a tag in the rwgt list. 913 max_wgt allow to do partial unweighting. 914 trunc_error allow for dynamical partial unweighting 915 event_target reweight for that many event with maximal trunc_error. 916 (stop to write event when target is reached) 917 """ 918 919 if isinstance(get_wgt, str): 920 unwgt_name =get_wgt 921 def get_wgt_multi(event): 922 event.parse_reweight() 923 return event.reweight_data[unwgt_name] * event.sample_scale
924 else: 925 unwgt_name = get_wgt.func_name 926 get_wgt_multi = lambda event: get_wgt(event) * event.sample_scale 927 #define the weighting such that we have built-in the scaling 928 929 930 if 'event_target' in opts and opts['event_target']: 931 if 'normalization' in opts: 932 if opts['normalization'] == 'sum': 933 new_wgt = sum(self.across)/opts['event_target'] 934 strategy = 3 935 elif opts['normalization'] == 'average': 936 strategy = 4 937 new_wgt = sum(self.across) 938 elif opts['normalization'] == 'unit': 939 strategy =3 940 new_wgt = 1. 941 else: 942 strategy = 4 943 new_wgt = sum(self.across) 944 self.define_init_banner(new_wgt, strategy) 945 self.written_weight = new_wgt 946 elif 'write_init' in opts and opts['write_init']: 947 self.define_init_banner(0,0) 948 del opts['write_init'] 949 return super(MultiEventFile, self).unweight(outputpath, get_wgt_multi, **opts)
950
951 - def write(self, path, random=False, banner=None, get_info=False):
952 """ """ 953 954 if isinstance(path, str): 955 out = EventFile(path, 'w') 956 if self.parsefile and not banner: 957 banner = self.files[0].banner 958 elif not banner: 959 firstlhe = EventFile(self.files[0]) 960 banner = firstlhe.banner 961 else: 962 out = path 963 if banner: 964 out.write(banner) 965 nb_event = 0 966 info = collections.defaultdict(float) 967 if random and self.open: 968 for event in self: 969 nb_event +=1 970 out.write(event) 971 if get_info: 972 event.parse_reweight() 973 for key, value in event.reweight_data.items(): 974 info[key] += value 975 info['central'] += event.wgt 976 elif not random: 977 for i,f in enumerate(self.files): 978 #check if we need to parse the file or not 979 if not self.parsefile: 980 if i==0: 981 try: 982 lhe = firstlhe 983 except: 984 lhe = EventFile(f) 985 else: 986 lhe = EventFile(f) 987 else: 988 lhe = f 989 for event in lhe: 990 nb_event +=1 991 if get_info: 992 event.parse_reweight() 993 for key, value in event.reweight_data.items(): 994 info[key] += value 995 info['central'] += event.wgt 996 out.write(str(event)) 997 lhe.close() 998 out.write("</LesHouchesEvents>\n") 999 return nb_event, info
1000
1001 - def remove(self):
1002 """ """ 1003 if self.parsefile: 1004 for f in self.files: 1005 os.remove(f.name) 1006 else: 1007 for f in self.files: 1008 os.remove(f)
1009
1010 1011 1012 -class Event(list):
1013 """Class storing a single event information (list of particles + global information)""" 1014 1015 warning_order = True # raise a warning if the order of the particle are not in accordance of child/mother 1016
1017 - def __init__(self, text=None):
1018 """The initialization of an empty Event (or one associate to a text file)""" 1019 list.__init__(self) 1020 1021 # First line information 1022 self.nexternal = 0 1023 self.ievent = 0 1024 self.wgt = 0 1025 self.aqcd = 0 1026 self.scale = 0 1027 self.aqed = 0 1028 self.aqcd = 0 1029 # Weight information 1030 self.tag = '' 1031 self.eventflag = {} # for information in <event > 1032 self.comment = '' 1033 self.reweight_data = {} 1034 self.matched_scale_data = None 1035 self.syscalc_data = {} 1036 if text: 1037 self.parse(text)
1038 1039 1040
1041 - def parse(self, text):
1042 """Take the input file and create the structured information""" 1043 #text = re.sub(r'</?event>', '', text) # remove pointless tag 1044 status = 'first' 1045 for line in text.split('\n'): 1046 line = line.strip() 1047 if not line: 1048 continue 1049 elif line[0] == '#': 1050 self.comment += '%s\n' % line 1051 continue 1052 elif line.startswith('<event'): 1053 if '=' in line: 1054 found = re.findall(r"""(\w*)=(?:(?:['"])([^'"]*)(?=['"])|(\S*))""",line) 1055 #for '<event line=4 value=\'3\' error="5" test=" 1 and 2">\n' 1056 #return [('line', '', '4'), ('value', '3', ''), ('error', '5', ''), ('test', ' 1 and 2', '')] 1057 self.eventflag = dict((n, a1) if a1 else (n,a2) for n,a1,a2 in found) 1058 # return {'test': ' 1 and 2', 'line': '4', 'value': '3', 'error': '5'} 1059 continue 1060 1061 elif 'first' == status: 1062 if '<rwgt>' in line: 1063 status = 'tag' 1064 else: 1065 self.assign_scale_line(line) 1066 status = 'part' 1067 continue 1068 if '<' in line: 1069 status = 'tag' 1070 1071 if 'part' == status: 1072 self.append(Particle(line, event=self)) 1073 else: 1074 if '</event>' in line: 1075 line = line.replace('</event>','',1) 1076 self.tag += '%s\n' % line 1077 1078 self.assign_mother()
1079
1080 - def assign_mother(self):
1081 # assign the mother: 1082 for i,particle in enumerate(self): 1083 if i < particle.mother1 or i < particle.mother2: 1084 if self.warning_order: 1085 logger.warning("Order of particle in the event did not agree with parent/child order. This might be problematic for some code.") 1086 Event.warning_order = False 1087 self.reorder_mother_child() 1088 return self.assign_mother() 1089 1090 if particle.mother1: 1091 try: 1092 particle.mother1 = self[int(particle.mother1) -1] 1093 except Exception: 1094 logger.warning("WRONG MOTHER INFO %s", self) 1095 particle.mother1 = 0 1096 if particle.mother2: 1097 try: 1098 particle.mother2 = self[int(particle.mother2) -1] 1099 except Exception: 1100 logger.warning("WRONG MOTHER INFO %s", self) 1101 particle.mother2 = 0
1102 1103
1104 - def reorder_mother_child(self):
1105 """check and correct the mother/child position. 1106 only correct one order by call (but this is a recursive call)""" 1107 1108 tomove, position = None, None 1109 for i,particle in enumerate(self): 1110 if i < particle.mother1: 1111 # move i after particle.mother1 1112 tomove, position = i, particle.mother1-1 1113 break 1114 if i < particle.mother2: 1115 tomove, position = i, particle.mother2-1 1116 1117 # nothing to change -> we are done 1118 if not tomove: 1119 return 1120 1121 # move the particles: 1122 particle = self.pop(tomove) 1123 self.insert(int(position), particle) 1124 1125 #change the mother id/ event_id in the event. 1126 for i, particle in enumerate(self): 1127 particle.event_id = i 1128 #misc.sprint( i, particle.event_id) 1129 m1, m2 = particle.mother1, particle.mother2 1130 if m1 == tomove +1: 1131 particle.mother1 = position+1 1132 elif tomove < m1 <= position +1: 1133 particle.mother1 -= 1 1134 if m2 == tomove +1: 1135 particle.mother2 = position+1 1136 elif tomove < m2 <= position +1: 1137 particle.mother2 -= 1 1138 # re-call the function for the next potential change 1139 return self.reorder_mother_child()
1140 1141 1142 1143 1144 1145
1146 - def parse_reweight(self):
1147 """Parse the re-weight information in order to return a dictionary 1148 {key: value}. If no group is define group should be '' """ 1149 if self.reweight_data: 1150 return self.reweight_data 1151 self.reweight_data = {} 1152 self.reweight_order = [] 1153 start, stop = self.tag.find('<rwgt>'), self.tag.find('</rwgt>') 1154 if start != -1 != stop : 1155 pattern = re.compile(r'''<\s*wgt id=(?:\'|\")(?P<id>[^\'\"]+)(?:\'|\")\s*>\s*(?P<val>[\ded+-.]*)\s*</wgt>''',re.I) 1156 data = pattern.findall(self.tag[start:stop]) 1157 try: 1158 self.reweight_data = dict([(pid, float(value)) for (pid, value) in data 1159 if not self.reweight_order.append(pid)]) 1160 # the if is to create the order file on the flight 1161 except ValueError, error: 1162 raise Exception, 'Event File has unvalid weight. %s' % error 1163 self.tag = self.tag[:start] + self.tag[stop+7:] 1164 return self.reweight_data
1165
1166 - def parse_nlo_weight(self, real_type=(1,11)):
1167 """ """ 1168 if hasattr(self, 'nloweight'): 1169 return self.nloweight 1170 1171 start, stop = self.tag.find('<mgrwgt>'), self.tag.find('</mgrwgt>') 1172 if start != -1 != stop : 1173 1174 text = self.tag[start+8:stop] 1175 self.nloweight = NLO_PARTIALWEIGHT(text, self, real_type=real_type) 1176 return self.nloweight
1177
1178 - def parse_lo_weight(self):
1179 """ """ 1180 if hasattr(self, 'loweight'): 1181 return self.loweight 1182 1183 start, stop = self.tag.find('<mgrwt>'), self.tag.find('</mgrwt>') 1184 1185 if start != -1 != stop : 1186 text = self.tag[start+8:stop] 1187 #<rscale> 3 0.29765919e+03</rscale> 1188 #<asrwt>0</asrwt> 1189 #<pdfrwt beam="1"> 1 21 0.15134321e+00 0.29765919e+03</pdfrwt> 1190 #<pdfrwt beam="2"> 1 21 0.38683649e-01 0.29765919e+03</pdfrwt> 1191 #<totfact> 0.17315115e+03</totfact> 1192 self.loweight={} 1193 for line in text.split('\n'): 1194 line = line.replace('<', ' <').replace("'",'"') 1195 if 'rscale' in line: 1196 _, nqcd, scale, _ = line.split() 1197 self.loweight['n_qcd'] = int(nqcd) 1198 self.loweight['ren_scale'] = float(scale) 1199 elif '<pdfrwt beam="1"' in line: 1200 args = line.split() 1201 self.loweight['n_pdfrw1'] = int(args[2]) 1202 npdf = self.loweight['n_pdfrw1'] 1203 self.loweight['pdf_pdg_code1'] = [int(i) for i in args[3:3+npdf]] 1204 self.loweight['pdf_x1'] = [float(i) for i in args[3+npdf:3+2*npdf]] 1205 self.loweight['pdf_q1'] = [float(i) for i in args[3+2*npdf:3+3*npdf]] 1206 elif '<pdfrwt beam="2"' in line: 1207 args = line.split() 1208 self.loweight['n_pdfrw2'] = int(args[2]) 1209 npdf = self.loweight['n_pdfrw2'] 1210 self.loweight['pdf_pdg_code2'] = [int(i) for i in args[3:3+npdf]] 1211 self.loweight['pdf_x2'] = [float(i) for i in args[3+npdf:3+2*npdf]] 1212 self.loweight['pdf_q2'] = [float(i) for i in args[3+2*npdf:3+3*npdf]] 1213 elif '<asrwt>' in line: 1214 args = line.replace('>','> ').split() 1215 nalps = int(args[1]) 1216 self.loweight['asrwt'] = [float(a) for a in args[2:2+nalps]] 1217 1218 elif 'totfact' in line: 1219 args = line.split() 1220 self.loweight['tot_fact'] = float(args[1]) 1221 else: 1222 return None 1223 return self.loweight
1224
1225 - def parse_matching_scale(self):
1226 """Parse the line containing the starting scale for the shower""" 1227 1228 if self.matched_scale_data is not None: 1229 return self.matched_scale_data 1230 1231 self.matched_scale_data = [] 1232 1233 1234 pattern = re.compile("<scales\s|</scales>") 1235 data = re.split(pattern,self.tag) 1236 if len(data) == 1: 1237 return [] 1238 else: 1239 tmp = {} 1240 start,content, end = data 1241 self.tag = "%s%s" % (start, end) 1242 pattern = re.compile("pt_clust_(\d*)=\"([\de+-.]*)\"") 1243 for id,value in pattern.findall(content): 1244 tmp[int(id)] = float(value) 1245 1246 for i in range(1, len(tmp)+1): 1247 self.matched_scale_data.append(tmp[i]) 1248 1249 return self.matched_scale_data
1250
1251 - def parse_syscalc_info(self):
1252 """ parse the flag for syscalc between <mgrwt></mgrwt> 1253 <mgrwt> 1254 <rscale> 3 0.26552898E+03</rscale> 1255 <asrwt>0</asrwt> 1256 <pdfrwt beam="1"> 1 21 0.14527945E+00 0.26552898E+03</pdfrwt> 1257 <pdfrwt beam="2"> 1 21 0.15249110E-01 0.26552898E+03</pdfrwt> 1258 <totfact> 0.10344054E+04</totfact> 1259 </mgrwt> 1260 """ 1261 if self.syscalc_data: 1262 return self.syscalc_data 1263 1264 pattern = re.compile("<mgrwt>|</mgrwt>") 1265 pattern2 = re.compile("<(?P<tag>[\w]*)(?:\s*(\w*)=[\"'](.*)[\"']\s*|\s*)>(.*)</(?P=tag)>") 1266 data = re.split(pattern,self.tag) 1267 if len(data) == 1: 1268 return [] 1269 else: 1270 tmp = {} 1271 start,content, end = data 1272 self.tag = "%s%s" % (start, end) 1273 for tag, key, keyval, tagval in pattern2.findall(content): 1274 if key: 1275 self.syscalc_data[(tag, key, keyval)] = tagval 1276 else: 1277 self.syscalc_data[tag] = tagval 1278 return self.syscalc_data
1279 1280
1281 - def add_decay_to_particle(self, position, decay_event):
1282 """define the decay of the particle id by the event pass in argument""" 1283 1284 this_particle = self[position] 1285 #change the status to internal particle 1286 this_particle.status = 2 1287 this_particle.helicity = 0 1288 1289 # some usefull information 1290 decay_particle = decay_event[0] 1291 this_4mom = FourMomentum(this_particle) 1292 nb_part = len(self) #original number of particle 1293 1294 thres = decay_particle.E*1e-10 1295 assert max(decay_particle.px, decay_particle.py, decay_particle.pz) < thres,\ 1296 "not on rest particle %s %s %s %s" % (decay_particle.E, decay_particle.px,decay_particle.py,decay_particle.pz) 1297 1298 self.nexternal += decay_event.nexternal -1 1299 old_scales = list(self.parse_matching_scale()) 1300 if old_scales: 1301 jet_position = sum(1 for i in range(position) if self[i].status==1) 1302 self.matched_scale_data.pop(jet_position) 1303 # add the particle with only handling the 4-momenta/mother 1304 # color information will be corrected later. 1305 for particle in decay_event[1:]: 1306 # duplicate particle to avoid border effect 1307 new_particle = Particle(particle, self) 1308 new_particle.event_id = len(self) 1309 self.append(new_particle) 1310 if old_scales: 1311 self.matched_scale_data.append(old_scales[jet_position]) 1312 # compute and assign the new four_momenta 1313 new_momentum = this_4mom.boost(FourMomentum(new_particle)) 1314 new_particle.set_momentum(new_momentum) 1315 # compute the new mother 1316 for tag in ['mother1', 'mother2']: 1317 mother = getattr(particle, tag) 1318 if isinstance(mother, Particle): 1319 mother_id = getattr(particle, tag).event_id 1320 if mother_id == 0: 1321 setattr(new_particle, tag, this_particle) 1322 else: 1323 try: 1324 setattr(new_particle, tag, self[nb_part + mother_id -1]) 1325 except Exception, error: 1326 print error 1327 misc.sprint( self) 1328 misc.sprint(nb_part + mother_id -1) 1329 misc.sprint(tag) 1330 misc.sprint(position, decay_event) 1331 misc.sprint(particle) 1332 misc.sprint(len(self), nb_part + mother_id -1) 1333 raise 1334 elif tag == "mother2" and isinstance(particle.mother1, Particle): 1335 new_particle.mother2 = this_particle 1336 else: 1337 raise Exception, "Something weird happens. Please report it for investigation" 1338 # Need to correct the color information of the particle 1339 # first find the first available color index 1340 max_color=501 1341 for particle in self[:nb_part]: 1342 max_color=max(max_color, particle.color1, particle.color2) 1343 1344 # define a color mapping and assign it: 1345 color_mapping = {} 1346 color_mapping[decay_particle.color1] = this_particle.color1 1347 color_mapping[decay_particle.color2] = this_particle.color2 1348 for particle in self[nb_part:]: 1349 if particle.color1: 1350 if particle.color1 not in color_mapping: 1351 max_color +=1 1352 color_mapping[particle.color1] = max_color 1353 particle.color1 = max_color 1354 else: 1355 particle.color1 = color_mapping[particle.color1] 1356 if particle.color2: 1357 if particle.color2 not in color_mapping: 1358 max_color +=1 1359 color_mapping[particle.color2] = max_color 1360 particle.color2 = max_color 1361 else: 1362 particle.color2 = color_mapping[particle.color2]
1363 1364 1365
1366 - def remove_decay(self, pdg_code=0, event_id=None):
1367 1368 to_remove = [] 1369 if event_id is not None: 1370 to_remove.append(self[event_id]) 1371 1372 if pdg_code: 1373 for particle in self: 1374 if particle.pid == pdg_code: 1375 to_remove.append(particle) 1376 1377 new_event = Event() 1378 # copy first line information + ... 1379 for tag in ['nexternal', 'ievent', 'wgt', 'aqcd', 'scale', 'aqed','tag','comment']: 1380 setattr(new_event, tag, getattr(self, tag)) 1381 1382 for particle in self: 1383 if isinstance(particle.mother1, Particle) and particle.mother1 in to_remove: 1384 to_remove.append(particle) 1385 if particle.status == 1: 1386 new_event.nexternal -= 1 1387 continue 1388 elif isinstance(particle.mother2, Particle) and particle.mother2 in to_remove: 1389 to_remove.append(particle) 1390 if particle.status == 1: 1391 new_event.nexternal -= 1 1392 continue 1393 else: 1394 new_event.append(Particle(particle)) 1395 1396 #ensure that the event_id is correct for all_particle 1397 # and put the status to 1 for removed particle 1398 for pos, particle in enumerate(new_event): 1399 particle.event_id = pos 1400 if particle in to_remove: 1401 particle.status = 1 1402 return new_event
1403
1404 - def get_decay(self, pdg_code=0, event_id=None):
1405 1406 to_start = [] 1407 if event_id is not None: 1408 to_start.append(self[event_id]) 1409 1410 elif pdg_code: 1411 for particle in self: 1412 if particle.pid == pdg_code: 1413 to_start.append(particle) 1414 break 1415 1416 new_event = Event() 1417 # copy first line information + ... 1418 for tag in ['ievent', 'wgt', 'aqcd', 'scale', 'aqed','tag','comment']: 1419 setattr(new_event, tag, getattr(self, tag)) 1420 1421 # Add the decaying particle 1422 old2new = {} 1423 new_decay_part = Particle(to_start[0]) 1424 new_decay_part.mother1 = None 1425 new_decay_part.mother2 = None 1426 new_decay_part.status = -1 1427 old2new[new_decay_part.event_id] = len(old2new) 1428 new_event.append(new_decay_part) 1429 1430 1431 # add the other particle 1432 for particle in self: 1433 if isinstance(particle.mother1, Particle) and particle.mother1.event_id in old2new\ 1434 or isinstance(particle.mother2, Particle) and particle.mother2.event_id in old2new: 1435 old2new[particle.event_id] = len(old2new) 1436 new_event.append(Particle(particle)) 1437 1438 #ensure that the event_id is correct for all_particle 1439 # and correct the mother1/mother2 by the new reference 1440 nexternal = 0 1441 for pos, particle in enumerate(new_event): 1442 particle.event_id = pos 1443 if particle.mother1: 1444 particle.mother1 = new_event[old2new[particle.mother1.event_id]] 1445 if particle.mother2: 1446 particle.mother2 = new_event[old2new[particle.mother2.event_id]] 1447 if particle.status in [-1,1]: 1448 nexternal +=1 1449 new_event.nexternal = nexternal 1450 1451 return new_event
1452 1453
1454 - def check(self):
1455 """check various property of the events""" 1456 1457 #1. Check that the 4-momenta are conserved 1458 E, px, py, pz = 0,0,0,0 1459 absE, abspx, abspy, abspz = 0,0,0,0 1460 for particle in self: 1461 coeff = 1 1462 if particle.status == -1: 1463 coeff = -1 1464 elif particle.status != 1: 1465 continue 1466 E += coeff * particle.E 1467 absE += abs(particle.E) 1468 px += coeff * particle.px 1469 py += coeff * particle.py 1470 pz += coeff * particle.pz 1471 abspx += abs(particle.px) 1472 abspy += abs(particle.py) 1473 abspz += abs(particle.pz) 1474 # check that relative error is under control 1475 threshold = 5e-7 1476 if E/absE > threshold: 1477 logger.critical(self) 1478 raise Exception, "Do not conserve Energy %s, %s" % (E/absE, E) 1479 if px/abspx > threshold: 1480 logger.critical(self) 1481 raise Exception, "Do not conserve Px %s, %s" % (px/abspx, px) 1482 if py/abspy > threshold: 1483 logger.critical(self) 1484 raise Exception, "Do not conserve Py %s, %s" % (py/abspy, py) 1485 if pz/abspz > threshold: 1486 logger.critical(self) 1487 raise Exception, "Do not conserve Pz %s, %s" % (pz/abspz, pz) 1488 1489 #2. check the color of the event 1490 self.check_color_structure()
1491
1492 - def assign_scale_line(self, line):
1493 """read the line corresponding to global event line 1494 format of the line is: 1495 Nexternal IEVENT WEIGHT SCALE AEW AS 1496 """ 1497 inputs = line.split() 1498 assert len(inputs) == 6 1499 self.nexternal=int(inputs[0]) 1500 self.ievent=int(inputs[1]) 1501 self.wgt=float(inputs[2]) 1502 self.scale=float(inputs[3]) 1503 self.aqed=float(inputs[4]) 1504 self.aqcd=float(inputs[5])
1505
1506 - def get_tag_and_order(self):
1507 """Return the unique tag identifying the SubProcesses for the generation. 1508 Usefull for program like MadSpin and Reweight module.""" 1509 1510 initial, final, order = [], [], [[], []] 1511 for particle in self: 1512 if particle.status == -1: 1513 initial.append(particle.pid) 1514 order[0].append(particle.pid) 1515 elif particle.status == 1: 1516 final.append(particle.pid) 1517 order[1].append(particle.pid) 1518 initial.sort(), final.sort() 1519 tag = (tuple(initial), tuple(final)) 1520 return tag, order
1521
1522 - def get_helicity(self, get_order, allow_reversed=True):
1523 """return a list with the helicities in the order asked for""" 1524 1525 #avoid to modify the input 1526 order = [list(get_order[0]), list(get_order[1])] 1527 out = [9] *(len(order[0])+len(order[1])) 1528 for i, part in enumerate(self): 1529 if part.status == 1: #final 1530 try: 1531 ind = order[1].index(part.pid) 1532 except ValueError, error: 1533 if not allow_reversed: 1534 raise error 1535 else: 1536 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 1537 try: 1538 return self.get_helicity(order, False) 1539 except ValueError: 1540 raise error 1541 position = len(order[0]) + ind 1542 order[1][ind] = 0 1543 elif part.status == -1: 1544 try: 1545 ind = order[0].index(part.pid) 1546 except ValueError, error: 1547 if not allow_reversed: 1548 raise error 1549 else: 1550 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 1551 try: 1552 return self.get_helicity(order, False) 1553 except ValueError: 1554 raise error 1555 1556 position = ind 1557 order[0][ind] = 0 1558 else: #intermediate 1559 continue 1560 out[position] = int(part.helicity) 1561 return out
1562 1563
1564 - def check_color_structure(self):
1565 """check the validity of the color structure""" 1566 1567 #1. check that each color is raised only once. 1568 color_index = collections.defaultdict(int) 1569 for particle in self: 1570 if particle.status in [-1,1]: 1571 if particle.color1: 1572 color_index[particle.color1] +=1 1573 if -7 < particle.pdg < 0: 1574 raise Exception, "anti-quark with color tag" 1575 if particle.color2: 1576 color_index[particle.color2] +=1 1577 if 7 > particle.pdg > 0: 1578 raise Exception, "quark with anti-color tag" 1579 1580 1581 for key,value in color_index.items(): 1582 if value > 2: 1583 print self 1584 print key, value 1585 raise Exception, 'Wrong color_flow' 1586 1587 1588 #2. check that each parent present have coherent color-structure 1589 check = [] 1590 popup_index = [] #check that the popup index are created in a unique way 1591 for particle in self: 1592 mothers = [] 1593 childs = [] 1594 if particle.mother1: 1595 mothers.append(particle.mother1) 1596 if particle.mother2 and particle.mother2 is not particle.mother1: 1597 mothers.append(particle.mother2) 1598 if not mothers: 1599 continue 1600 if (particle.mother1.event_id, particle.mother2.event_id) in check: 1601 continue 1602 check.append((particle.mother1.event_id, particle.mother2.event_id)) 1603 1604 childs = [p for p in self if p.mother1 is particle.mother1 and \ 1605 p.mother2 is particle.mother2] 1606 1607 mcolors = [] 1608 manticolors = [] 1609 for m in mothers: 1610 if m.color1: 1611 if m.color1 in manticolors: 1612 manticolors.remove(m.color1) 1613 else: 1614 mcolors.append(m.color1) 1615 if m.color2: 1616 if m.color2 in mcolors: 1617 mcolors.remove(m.color2) 1618 else: 1619 manticolors.append(m.color2) 1620 ccolors = [] 1621 canticolors = [] 1622 for m in childs: 1623 if m.color1: 1624 if m.color1 in canticolors: 1625 canticolors.remove(m.color1) 1626 else: 1627 ccolors.append(m.color1) 1628 if m.color2: 1629 if m.color2 in ccolors: 1630 ccolors.remove(m.color2) 1631 else: 1632 canticolors.append(m.color2) 1633 for index in mcolors[:]: 1634 if index in ccolors: 1635 mcolors.remove(index) 1636 ccolors.remove(index) 1637 for index in manticolors[:]: 1638 if index in canticolors: 1639 manticolors.remove(index) 1640 canticolors.remove(index) 1641 1642 if mcolors != []: 1643 #only case is a epsilon_ijk structure. 1644 if len(canticolors) + len(mcolors) != 3: 1645 logger.critical(str(self)) 1646 raise Exception, "Wrong color flow for %s -> %s" ([m.pid for m in mothers], [c.pid for c in childs]) 1647 else: 1648 popup_index += canticolors 1649 elif manticolors != []: 1650 #only case is a epsilon_ijk structure. 1651 if len(ccolors) + len(manticolors) != 3: 1652 logger.critical(str(self)) 1653 raise Exception, "Wrong color flow for %s -> %s" ([m.pid for m in mothers], [c.pid for c in childs]) 1654 else: 1655 popup_index += ccolors 1656 1657 # Check that color popup (from epsilon_ijk) are raised only once 1658 if len(popup_index) != len(set(popup_index)): 1659 logger.critical(self) 1660 raise Exception, "Wrong color flow: identical poping-up index, %s" % (popup_index)
1661
1662 - def __str__(self, event_id=''):
1663 """return a correctly formatted LHE event""" 1664 1665 out="""<event%(event_flag)s> 1666 %(scale)s 1667 %(particles)s 1668 %(comments)s 1669 %(tag)s 1670 %(reweight)s 1671 </event> 1672 """ 1673 if event_id not in ['', None]: 1674 self.eventflag['event'] = str(event_id) 1675 1676 if self.eventflag: 1677 event_flag = ' %s' % ' '.join('%s="%s"' % (k,v) for (k,v) in self.eventflag.items()) 1678 else: 1679 event_flag = '' 1680 1681 if self.nexternal: 1682 scale_str = "%2d %6d %+13.7e %14.8e %14.8e %14.8e" % \ 1683 (self.nexternal,self.ievent,self.wgt,self.scale,self.aqed,self.aqcd) 1684 else: 1685 scale_str = '' 1686 1687 if self.reweight_data: 1688 # check that all key have an order if not add them at the end 1689 if set(self.reweight_data.keys()) != set(self.reweight_order): 1690 self.reweight_order += [k for k in self.reweight_data.keys() \ 1691 if k not in self.reweight_order] 1692 1693 reweight_str = '<rwgt>\n%s\n</rwgt>' % '\n'.join( 1694 '<wgt id=\'%s\'> %+13.7e </wgt>' % (i, float(self.reweight_data[i])) 1695 for i in self.reweight_order) 1696 else: 1697 reweight_str = '' 1698 1699 tag_str = self.tag 1700 if self.matched_scale_data: 1701 tag_str = "<scales %s></scales>%s" % ( 1702 ' '.join(['pt_clust_%i=\"%s\"' % (i,v) 1703 for i,v in enumerate(self.matched_scale_data)]), 1704 self.tag) 1705 if self.syscalc_data: 1706 keys= ['rscale', 'asrwt', ('pdfrwt', 'beam', '1'), ('pdfrwt', 'beam', '2'), 1707 'matchscale', 'totfact'] 1708 sys_str = "<mgrwt>\n" 1709 template = """<%(key)s%(opts)s>%(values)s</%(key)s>\n""" 1710 for k in keys: 1711 if k not in self.syscalc_data: 1712 continue 1713 replace = {} 1714 replace['values'] = self.syscalc_data[k] 1715 if isinstance(k, str): 1716 replace['key'] = k 1717 replace['opts'] = '' 1718 else: 1719 replace['key'] = k[0] 1720 replace['opts'] = ' %s=\"%s\"' % (k[1],k[2]) 1721 sys_str += template % replace 1722 sys_str += "</mgrwt>\n" 1723 reweight_str = sys_str + reweight_str 1724 1725 out = out % {'event_flag': event_flag, 1726 'scale': scale_str, 1727 'particles': '\n'.join([str(p) for p in self]), 1728 'tag': tag_str, 1729 'comments': self.comment, 1730 'reweight': reweight_str} 1731 1732 return re.sub('[\n]+', '\n', out)
1733
1734 - def get_momenta(self, get_order, allow_reversed=True):
1735 """return the momenta vector in the order asked for""" 1736 1737 #avoid to modify the input 1738 order = [list(get_order[0]), list(get_order[1])] 1739 out = [''] *(len(order[0])+len(order[1])) 1740 for i, part in enumerate(self): 1741 if part.status == 1: #final 1742 try: 1743 ind = order[1].index(part.pid) 1744 except ValueError, error: 1745 if not allow_reversed: 1746 raise error 1747 else: 1748 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 1749 try: 1750 return self.get_momenta_str(order, False) 1751 except ValueError: 1752 raise error 1753 position = len(order[0]) + ind 1754 order[1][ind] = 0 1755 elif part.status == -1: 1756 try: 1757 ind = order[0].index(part.pid) 1758 except ValueError, error: 1759 if not allow_reversed: 1760 raise error 1761 else: 1762 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 1763 try: 1764 return self.get_momenta_str(order, False) 1765 except ValueError: 1766 raise error 1767 1768 position = ind 1769 order[0][ind] = 0 1770 else: #intermediate 1771 continue 1772 1773 out[position] = (part.E, part.px, part.py, part.pz) 1774 1775 return out
1776 1777 1778
1779 - def get_ht_scale(self, prefactor=1):
1780 1781 scale = 0 1782 for particle in self: 1783 if particle.status != 1: 1784 continue 1785 p=FourMomentum(particle) 1786 scale += math.sqrt(p.mass_sqr + p.pt**2) 1787 1788 return prefactor * scale
1789 1790
1791 - def get_et_scale(self, prefactor=1):
1792 1793 scale = 0 1794 for particle in self: 1795 if particle.status != 1: 1796 continue 1797 p = FourMomentum(particle) 1798 pt = p.pt 1799 if (pt>0): 1800 scale += p.E*pt/math.sqrt(pt**2+p.pz**2) 1801 1802 return prefactor * scale
1803
1804 - def get_sqrts_scale(self, prefactor=1):
1805 1806 scale = 0 1807 init = [] 1808 for particle in self: 1809 if particle.status == -1: 1810 init.append(FourMomentum(particle)) 1811 if len(init) == 1: 1812 return init[0].mass 1813 elif len(init)==2: 1814 return math.sqrt((init[0]+init[1])**2)
1815 1816 1817 1818
1819 - def get_momenta_str(self, get_order, allow_reversed=True):
1820 """return the momenta str in the order asked for""" 1821 1822 out = self.get_momenta(get_order, allow_reversed) 1823 #format 1824 format = '%.12f' 1825 format_line = ' '.join([format]*4) + ' \n' 1826 out = [format_line % one for one in out] 1827 out = ''.join(out).replace('e','d') 1828 return out
1829
1830 -class WeightFile(EventFile):
1831 """A class to allow to read both gzip and not gzip file. 1832 containing only weight from pythia --generated by SysCalc""" 1833
1834 - def __new__(self, path, mode='r', *args, **opt):
1835 if path.endswith(".gz"): 1836 try: 1837 return gzip.GzipFile.__new__(WeightFileGzip, path, mode, *args, **opt) 1838 except IOError, error: 1839 raise 1840 except Exception, error: 1841 if mode == 'r': 1842 misc.gunzip(path) 1843 return file.__new__(WeightFileNoGzip, path[:-3], mode, *args, **opt) 1844 else: 1845 return file.__new__(WeightFileNoGzip, path, mode, *args, **opt)
1846 1847
1848 - def __init__(self, path, mode='r', *args, **opt):
1849 """open file and read the banner [if in read mode]""" 1850 1851 super(EventFile, self).__init__(path, mode, *args, **opt) 1852 self.banner = '' 1853 if mode == 'r': 1854 line = '' 1855 while '</header>' not in line.lower(): 1856 try: 1857 line = super(EventFile, self).next() 1858 except StopIteration: 1859 self.seek(0) 1860 self.banner = '' 1861 break 1862 if "<event" in line.lower(): 1863 self.seek(0) 1864 self.banner = '' 1865 break 1866 1867 self.banner += line
1868
1869 1870 -class WeightFileGzip(WeightFile, EventFileGzip):
1871 pass
1872
1873 -class WeightFileNoGzip(WeightFile, EventFileNoGzip):
1874 pass
1875
1876 1877 -class FourMomentum(object):
1878 """a convenient object for 4-momenta operation""" 1879
1880 - def __init__(self, obj=0, px=0, py=0, pz=0, E=0):
1881 """initialize the four momenta""" 1882 1883 if obj is 0 and E: 1884 obj = E 1885 1886 if isinstance(obj, (FourMomentum, Particle)): 1887 px = obj.px 1888 py = obj.py 1889 pz = obj.pz 1890 E = obj.E 1891 elif isinstance(obj, (list, tuple)): 1892 assert len(obj) ==4 1893 E = obj[0] 1894 px = obj[1] 1895 py = obj[2] 1896 pz = obj[3] 1897 elif isinstance(obj, str): 1898 obj = [float(i) for i in obj.split()] 1899 assert len(obj) ==4 1900 E = obj[0] 1901 px = obj[1] 1902 py = obj[2] 1903 pz = obj[3] 1904 else: 1905 E =obj 1906 1907 1908 self.E = float(E) 1909 self.px = float(px) 1910 self.py = float(py) 1911 self.pz = float(pz)
1912 1913 @property
1914 - def mass(self):
1915 """return the mass""" 1916 return math.sqrt(self.E**2 - self.px**2 - self.py**2 - self.pz**2)
1917 1918 @property
1919 - def mass_sqr(self):
1920 """return the mass square""" 1921 return self.E**2 - self.px**2 - self.py**2 - self.pz**2
1922 1923 @property
1924 - def pt(self):
1925 return math.sqrt(max(0, self.pt2()))
1926 1927 @property
1928 - def pseudorapidity(self):
1929 norm = math.sqrt(self.px**2 + self.py**2+self.pz**2) 1930 return 0.5* math.log((norm - self.pz) / (norm + self.pz))
1931 1932 @property
1933 - def rapidity(self):
1934 return 0.5* math.log((self.E +self.pz) / (self.E - self.pz))
1935 1936 1937
1938 - def pt2(self):
1939 """ return the pt square """ 1940 1941 return self.px**2 + self.py**2
1942
1943 - def __add__(self, obj):
1944 1945 assert isinstance(obj, FourMomentum) 1946 new = FourMomentum(self.E+obj.E, 1947 self.px + obj.px, 1948 self.py + obj.py, 1949 self.pz + obj.pz) 1950 return new
1951
1952 - def __iadd__(self, obj):
1953 """update the object with the sum""" 1954 self.E += obj.E 1955 self.px += obj.px 1956 self.py += obj.py 1957 self.pz += obj.pz 1958 return self
1959
1960 - def __sub__(self, obj):
1961 1962 assert isinstance(obj, FourMomentum) 1963 new = FourMomentum(self.E-obj.E, 1964 self.px - obj.px, 1965 self.py - obj.py, 1966 self.pz - obj.pz) 1967 return new
1968
1969 - def __isub__(self, obj):
1970 """update the object with the sum""" 1971 self.E -= obj.E 1972 self.px -= obj.px 1973 self.py -= obj.py 1974 self.pz -= obj.pz 1975 return self
1976
1977 - def __mul__(self, obj):
1978 if isinstance(obj, FourMomentum): 1979 return self.E*obj.E - self.px *obj.px - self.py * obj.py - self.pz * obj.pz 1980 elif isinstance(obj, (float, int)): 1981 return FourMomentum(obj*self.E,obj*self.px,obj*self.py,obj*self.pz ) 1982 else: 1983 raise NotImplemented
1984 __rmul__ = __mul__ 1985
1986 - def __pow__(self, power):
1987 assert power in [1,2] 1988 1989 if power == 1: 1990 return FourMomentum(self) 1991 elif power == 2: 1992 return self.mass_sqr
1993
1994 - def __repr__(self):
1995 return 'FourMomentum(%s,%s,%s,%s)' % (self.E, self.px, self.py,self.pz)
1996
1997 - def get_tuple(self):
1998 return (self.E, self.px, self.py,self.pz)
1999
2000 - def boost(self, mom):
2001 """mom 4-momenta is suppose to be given in the rest frame of this 4-momenta. 2002 the output is the 4-momenta in the frame of this 4-momenta 2003 function copied from HELAS routine.""" 2004 2005 2006 pt = self.px**2 + self.py**2 + self.pz**2 2007 if pt: 2008 s3product = self.px * mom.px + self.py * mom.py + self.pz * mom.pz 2009 mass = self.mass 2010 lf = (mom.E + (self.E - mass) * s3product / pt ) / mass 2011 return FourMomentum(E=(self.E*mom.E+s3product)/mass, 2012 px=mom.px + self.px * lf, 2013 py=mom.py + self.py * lf, 2014 pz=mom.pz + self.pz * lf) 2015 else: 2016 return FourMomentum(mom)
2017
2018 - def zboost(self, pboost=None, E=0, pz=0):
2019 """Both momenta should be in the same frame. 2020 The boost perform correspond to the boost required to set pboost at 2021 rest (only z boost applied). 2022 """ 2023 if isinstance(pboost, FourMomentum): 2024 E = pboost.E 2025 pz = pboost.pz 2026 2027 #beta = pz/E 2028 gamma = E / math.sqrt(E**2-pz**2) 2029 gammabeta = pz / math.sqrt(E**2-pz**2) 2030 2031 out = FourMomentum([gamma*self.E - gammabeta*self.pz, 2032 self.px, 2033 self.py, 2034 gamma*self.pz - gammabeta*self.E]) 2035 2036 if abs(out.pz) < 1e-6 * out.E: 2037 out.pz = 0 2038 return out
2039
2040 2041 2042 2043 -class OneNLOWeight(object):
2044
2045 - def __init__(self, input, real_type=(1,11)):
2046 """ """ 2047 2048 self.real_type = real_type 2049 if isinstance(input, str): 2050 self.parse(input)
2051
2052 - def __str__(self):
2053 2054 out = """ pwgt: %(pwgt)s 2055 born, real : %(born)s %(real)s 2056 pdgs : %(pdgs)s 2057 bjks : %(bjks)s 2058 scales**2, gs: %(scales2)s %(gs)s 2059 born/real related : %(born_related)s %(real_related)s 2060 type / nfks : %(type)s %(nfks)s 2061 to merge : %(to_merge_pdg)s in %(merge_new_pdg)s 2062 ref_wgt : %(ref_wgt)s""" % self.__dict__ 2063 return out
2064 2065
2066 - def parse(self, text):
2067 """parse the line and create the related object""" 2068 #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 2069 # below comment are from Rik description email 2070 2071 data = text.split() 2072 # 1. The first three doubles are, as before, the 'wgt', i.e., the overall event of this 2073 # contribution, and the ones multiplying the log[mu_R/QES] and the log[mu_F/QES] 2074 # stripped of alpha_s and the PDFs. 2075 self.pwgt = [float(f) for f in data[:3]] 2076 # 2. The next two doubles are the values of the (corresponding) Born and 2077 # real-emission matrix elements. You can either use these values to check 2078 # that the newly computed original matrix element weights are correct, 2079 # or directly use these so that you don't have to recompute the original weights. 2080 # For contributions for which the real-emission matrix elements were 2081 # not computed, the 2nd of these numbers is zero. The opposite is not true, 2082 # because each real-emission phase-space configuration has an underlying Born one 2083 # (this is not unique, but on our code we made a specific choice here). 2084 # This latter information is useful if the real-emission matrix elements 2085 # are unstable; you can then reweight with the Born instead. 2086 # (see also point 9 below, where the momentum configurations are assigned). 2087 # I don't think this instability is real problem when reweighting the real-emission 2088 # with tree-level matrix elements (as we generally would do), but is important 2089 # when reweighting with loop-squared contributions as we have been doing for gg->H. 2090 # (I'm not sure that reweighting tree-level with loop^2 is something that 2091 # we can do in general, because we don't really know what to do with the 2092 # virtual matrix elements because we cannot generate 2-loop diagrams.) 2093 self.born = float(data[3]) 2094 self.real = float(data[4]) 2095 # 3. integer: number of external particles of the real-emission configuration (as before) 2096 self.nexternal = int(data[5]) 2097 # 4. PDG codes corresponding to the real-emission configuration (as before) 2098 self.pdgs = [int(i) for i in data[6:6+self.nexternal]] 2099 flag = 6+self.nexternal # new starting point for the position 2100 # 5. next integer is the power of g_strong in the matrix elements (as before) 2101 self.qcdpower = int(data[flag]) 2102 # 6. 2 doubles: The bjorken x's used for this contribution (as before) 2103 self.bjks = [float(f) for f in data[flag+1:flag+3]] 2104 # 7. 3 doubles: The Ellis-sexton scale, the renormalisation scale and the factorisation scale, all squared, used for this contribution (as before) 2105 self.scales2 = [float(f) for f in data[flag+3:flag+6]] 2106 # 8.the value of g_strong 2107 self.gs = float(data[flag+6]) 2108 # 9. 2 integers: the corresponding Born and real-emission type kinematics. (in the list of momenta) 2109 # Note that also the Born-kinematics has n+1 particles, with, in general, 2110 # one particle with zero momentum (this is not ALWAYS the case, 2111 # there could also be 2 particles with perfectly collinear momentum). 2112 # To convert this from n+1 to a n particles, you have to sum the momenta 2113 # of the two particles that 'merge', see point 12 below. 2114 self.born_related = int(data[flag+7]) 2115 self.real_related = int(data[flag+8]) 2116 # 10. 1 integer: the 'type'. This is the information you should use to determine 2117 # if to reweight with Born, virtual or real-emission matrix elements. 2118 # (Apart from the possible problems with complicated real-emission matrix elements 2119 # that need to be computed very close to the soft/collinear limits, see point 2 above. 2120 # I guess that for tree-level this is always okay, but when reweighting 2121 # a tree-level contribution with a one-loop squared one, as we do 2122 # for gg->Higgs, this is important). 2123 # type=1 : real-emission: 2124 # type=2 : Born: 2125 # type=3 : integrated counter terms: 2126 # type=4 : soft counter-term : 2127 # type=5 : collinear counter-term : 2128 # type=6 : soft-collinear counter-term: 2129 # type=7 : O(alphaS) expansion of Sudakov factor for NNLL+NLO : 2130 # type=8 : soft counter-term (with n+1-body kin.): 2131 # type=9 : collinear counter-term (with n+1-body kin.): 2132 # type=10: soft-collinear counter-term (with n+1-body kin.): 2133 # type=11: real-emission (with n-body kin.): 2134 # type=12: MC subtraction with n-body kin.: 2135 # type=13: MC subtraction with n+1-body kin.: 2136 # type=14: virtual corrections minus approximate virtual 2137 # type=15: approximate virtual corrections: 2138 self.type = int(data[flag+9]) 2139 # 11. 1 integer: The FKS configuration for this contribution (not really 2140 # relevant for anything, but is used in checking the reweighting to 2141 # get scale & PDF uncertainties). 2142 self.nfks = int(data[flag+10]) 2143 # 12. 2 integers: the two particles that should be merged to form the 2144 # born contribution from the real-emission one. Remove these two particles 2145 # from the (ordered) list of PDG codes, and insert a newly created particle 2146 # at the location of the minimum of the two particles removed. 2147 # I.e., if you merge particles 2 and 4, you have to insert the new particle 2148 # as the 2nd particle. And particle 5 and above will be shifted down by one. 2149 self.to_merge_pdg = [int (f) for f in data[flag+11:flag+13]] 2150 # 13. 1 integer: the PDG code of the particle that is created after merging the two particles at point 12. 2151 self.merge_new_pdg = int(data[flag+13]) 2152 # 14. 1 double: the reference number that one should be able to reconstruct 2153 # form the weights (point 1 above) and the rest of the information of this line. 2154 # This is really the contribution to this event as computed by the code 2155 # (and is passed to the integrator). It contains everything. 2156 self.ref_wgt = float(data[flag+14]) 2157 2158 #check the momenta configuration linked to the event 2159 if self.type in self.real_type: 2160 self.momenta_config = self.real_related 2161 else: 2162 self.momenta_config = self.born_related
2163
2164 2165 -class NLO_PARTIALWEIGHT(object):
2166
2167 - class BasicEvent(list):
2168
2169 - def __init__(self, momenta, wgts, event, real_type=(1,11)):
2170 list.__init__(self, momenta) 2171 2172 assert self 2173 self.wgts = wgts 2174 self.pdgs = list(wgts[0].pdgs) 2175 self.event = event 2176 self.real_type = real_type 2177 2178 if wgts[0].momenta_config == wgts[0].born_related: 2179 # need to remove one momenta. 2180 ind1, ind2 = [ind-1 for ind in wgts[0].to_merge_pdg] 2181 if ind1> ind2: 2182 ind1, ind2 = ind2, ind1 2183 if ind1 >= sum(1 for p in event if p.status==-1): 2184 new_p = self[ind1] + self[ind2] 2185 else: 2186 new_p = self[ind1] - self[ind2] 2187 self.pop(ind1) 2188 self.insert(ind1, new_p) 2189 self.pop(ind2) 2190 self.pdgs.pop(ind1) 2191 self.pdgs.insert(ind1, wgts[0].merge_new_pdg ) 2192 self.pdgs.pop(ind2) 2193 # DO NOT update the pdgs of the partial weight! 2194 elif any(w.type in self.real_type for w in wgts): 2195 if any(w.type not in self.real_type for w in wgts): 2196 raise Exception 2197 # check if this is too soft/colinear if so use the born 2198 ind1, ind2 = [ind-1 for ind in wgts[0].to_merge_pdg] 2199 if ind1> ind2: 2200 ind1, ind2 = ind2, ind1 2201 if ind1 >= sum(1 for p in event if p.status==-1): 2202 new_p = self[ind1] + self[ind2] 2203 else: 2204 new_p = self[ind1] - self[ind2] 2205 2206 if __debug__: 2207 ptot = FourMomentum() 2208 for i in xrange(len(self)): 2209 if i <2: 2210 ptot += self[i] 2211 else: 2212 ptot -= self[i] 2213 if ptot.mass_sqr > 1e-16: 2214 misc.sprint(ptot, ptot.mass_sqr) 2215 2216 inv_mass = new_p.mass_sqr 2217 shat = (self[0]+self[1]).mass_sqr 2218 if (abs(inv_mass)/shat < 1e-6): 2219 # misc.sprint(abs(inv_mass)/shat) 2220 self.pop(ind1) 2221 self.insert(ind1, new_p) 2222 self.pop(ind2) 2223 self.pdgs.pop(ind1) 2224 self.pdgs.insert(ind1, wgts[0].merge_new_pdg ) 2225 self.pdgs.pop(ind2) 2226 # DO NOT update the pdgs of the partial weight! 2227 2228 if __debug__: 2229 ptot = FourMomentum() 2230 for i in xrange(len(self)): 2231 if i <2: 2232 ptot += self[i] 2233 else: 2234 ptot -= self[i] 2235 if ptot.mass_sqr > 1e-16: 2236 misc.sprint(ptot, ptot.mass_sqr) 2237 # raise Exception 2238 else: 2239 raise Exception
2240
2241 - def get_pdg_code(self):
2242 return self.pdgs
2243
2244 - def get_tag_and_order(self):
2245 """ return the tag and order for this basic event""" 2246 (initial, _), _ = self.event.get_tag_and_order() 2247 order = self.get_pdg_code() 2248 2249 2250 initial, out = order[:len(initial)], order[len(initial):] 2251 initial.sort() 2252 out.sort() 2253 return (tuple(initial), tuple(out)), order
2254
2255 - def get_momenta(self, get_order, allow_reversed=True):
2256 """return the momenta vector in the order asked for""" 2257 2258 #avoid to modify the input 2259 order = [list(get_order[0]), list(get_order[1])] 2260 out = [''] *(len(order[0])+len(order[1])) 2261 pdgs = self.get_pdg_code() 2262 for pos, part in enumerate(self): 2263 if pos < len(get_order[0]): #initial 2264 try: 2265 ind = order[0].index(pdgs[pos]) 2266 except ValueError, error: 2267 if not allow_reversed: 2268 raise error 2269 else: 2270 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 2271 try: 2272 return self.get_momenta(order, False) 2273 except ValueError: 2274 raise error 2275 2276 2277 position = ind 2278 order[0][ind] = 0 2279 else: #final 2280 try: 2281 ind = order[1].index(pdgs[pos]) 2282 except ValueError, error: 2283 if not allow_reversed: 2284 raise error 2285 else: 2286 order = [[-i for i in get_order[0]],[-i for i in get_order[1]]] 2287 try: 2288 return self.get_momenta(order, False) 2289 except ValueError: 2290 raise error 2291 position = len(order[0]) + ind 2292 order[1][ind] = 0 2293 2294 out[position] = (part.E, part.px, part.py, part.pz) 2295 2296 return out
2297 2298
2299 - def get_helicity(self, *args):
2300 return [9] * len(self)
2301 2302 @property
2303 - def aqcd(self):
2304 return self.event.aqcd
2305
2306 - def get_ht_scale(self, prefactor=1):
2307 2308 scale = 0 2309 for particle in self: 2310 p = particle 2311 scale += math.sqrt(max(0, p.mass_sqr + p.pt**2)) 2312 2313 return prefactor * scale
2314
2315 - def get_et_scale(self, prefactor=1):
2316 2317 scale = 0 2318 for particle in self: 2319 p = particle 2320 pt = p.pt 2321 if (pt>0): 2322 scale += p.E*pt/math.sqrt(pt**2+p.pz**2) 2323 2324 return prefactor * scale
2325 2326
2327 - def get_sqrts_scale(self, event,prefactor=1):
2328 2329 scale = 0 2330 nb_init = 0 2331 for particle in event: 2332 if particle.status == -1: 2333 nb_init+=1 2334 if nb_init == 1: 2335 return self[0].mass 2336 elif nb_init==2: 2337 return math.sqrt((self[0]+self[1])**2)
2338 2339 2340 2341
2342 - def __init__(self, input, event, real_type=(1,11)):
2343 2344 self.real_type = real_type 2345 self.event = event 2346 if isinstance(input, str): 2347 self.parse(input)
2348 2349
2350 - def parse(self, text):
2351 """create the object from the string information (see example below)""" 2352 #0.2344688900d+00 8 2 0 2353 #0.4676614699d+02 0.0000000000d+00 0.0000000000d+00 0.4676614699d+02 2354 #0.4676614699d+02 0.0000000000d+00 0.0000000000d+00 -.4676614699d+02 2355 #0.4676614699d+02 0.2256794794d+02 0.4332148227d+01 0.4073073437d+02 2356 #0.4676614699d+02 -.2256794794d+02 -.4332148227d+01 -.4073073437d+02 2357 #0.0000000000d+00 -.0000000000d+00 -.0000000000d+00 -.0000000000d+00 2358 #0.4780341163d+02 0.0000000000d+00 0.0000000000d+00 0.4780341163d+02 2359 #0.4822581633d+02 0.0000000000d+00 0.0000000000d+00 -.4822581633d+02 2360 #0.4729127470d+02 0.2347155377d+02 0.5153455534d+01 0.4073073437d+02 2361 #0.4627255267d+02 -.2167412893d+02 -.3519736379d+01 -.4073073437d+02 2362 #0.2465400591d+01 -.1797424844d+01 -.1633719155d+01 -.4224046944d+00 2363 #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 2364 #-.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 2365 #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 2366 #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 2367 #-.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 2368 #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 2369 #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 2370 #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 2371 2372 text = text.lower().replace('d','e') 2373 all_line = text.split('\n') 2374 #get global information 2375 first_line ='' 2376 while not first_line.strip(): 2377 first_line = all_line.pop(0) 2378 2379 wgt, nb_wgt, nb_event, _ = first_line.split() 2380 nb_wgt, nb_event = int(nb_wgt), int(nb_event) 2381 2382 momenta = [] 2383 wgts = [] 2384 for line in all_line: 2385 data = line.split() 2386 if len(data) == 4: 2387 p = FourMomentum(data) 2388 momenta.append(p) 2389 elif len(data)>0: 2390 wgt = OneNLOWeight(line, real_type=self.real_type) 2391 wgts.append(wgt) 2392 2393 assert len(wgts) == int(nb_wgt) 2394 2395 get_weights_for_momenta = {} 2396 size_momenta = 0 2397 for wgt in wgts: 2398 if wgt.momenta_config in get_weights_for_momenta: 2399 get_weights_for_momenta[wgt.momenta_config].append(wgt) 2400 else: 2401 if size_momenta == 0: size_momenta = wgt.nexternal 2402 assert size_momenta == wgt.nexternal 2403 get_weights_for_momenta[wgt.momenta_config] = [wgt] 2404 2405 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) 2406 2407 2408 2409 self.cevents = [] 2410 for key in range(1, nb_event+1): 2411 if key in get_weights_for_momenta: 2412 wgt = get_weights_for_momenta[key] 2413 evt = self.BasicEvent(momenta[:size_momenta], get_weights_for_momenta[key], self.event, self.real_type) 2414 self.cevents.append(evt) 2415 momenta = momenta[size_momenta:] 2416 2417 nb_wgt_check = 0 2418 for cevt in self.cevents: 2419 nb_wgt_check += len(cevt.wgts) 2420 assert nb_wgt_check == int(nb_wgt)
2421 2422 2423 2424 if '__main__' == __name__: 2425 2426 # Example 1: adding some missing information to the event (here distance travelled) 2427 if False: 2428 lhe = EventFile('unweighted_events.lhe.gz') 2429 output = open('output_events.lhe', 'w') 2430 #write the banner to the output file 2431 output.write(lhe.banner) 2432 # Loop over all events 2433 for event in lhe: 2434 for particle in event: 2435 # modify particle attribute: here remove the mass 2436 particle.mass = 0 2437 particle.vtim = 2 # The one associate to distance travelled by the particle. 2438 2439 #write this modify event 2440 output.write(str(event)) 2441 output.write('</LesHouchesEvent>\n') 2442 2443 # Example 3: Plotting some variable 2444 if False: 2445 lhe = EventFile('unweighted_events.lhe.gz') 2446 import matplotlib.pyplot as plt 2447 import matplotlib.gridspec as gridspec 2448 nbins = 100 2449 2450 nb_pass = 0 2451 data = [] 2452 for event in lhe: 2453 etaabs = 0 2454 etafinal = 0 2455 for particle in event: 2456 if particle.status==1: 2457 p = FourMomentum(particle) 2458 eta = p.pseudorapidity 2459 if abs(eta) > etaabs: 2460 etafinal = eta 2461 etaabs = abs(eta) 2462 if etaabs < 4: 2463 data.append(etafinal) 2464 nb_pass +=1 2465 2466 2467 print nb_pass 2468 gs1 = gridspec.GridSpec(2, 1, height_ratios=[5,1]) 2469 gs1.update(wspace=0, hspace=0) # set the spacing between axes. 2470 ax = plt.subplot(gs1[0]) 2471 2472 n, bins, patches = ax.hist(data, nbins, histtype='step', label='original') 2473 ax_c = ax.twinx() 2474 ax_c.set_ylabel('MadGraph5_aMC@NLO') 2475 ax_c.yaxis.set_label_coords(1.01, 0.25) 2476 ax_c.set_yticks(ax.get_yticks()) 2477 ax_c.set_yticklabels([]) 2478 ax.set_xlim([-4,4]) 2479 print "bin value:", n 2480 print "start/end point of bins", bins 2481 plt.axis('on') 2482 plt.xlabel('weight ratio') 2483 plt.show() 2484 2485 2486 # Example 4: More complex plotting example (with ratio plot) 2487 if False: 2488 lhe = EventFile('unweighted_events.lhe') 2489 import matplotlib.pyplot as plt 2490 import matplotlib.gridspec as gridspec 2491 nbins = 100 2492 2493 #mtau, wtau = 45, 5.1785e-06 2494 mtau, wtau = 1.777, 4.027000e-13 2495 nb_pass = 0 2496 data, data2, data3 = [], [], [] 2497 for event in lhe: 2498 nb_pass +=1 2499 if nb_pass > 10000: 2500 break 2501 tau1 = FourMomentum() 2502 tau2 = FourMomentum() 2503 for part in event: 2504 if part.pid in [-12,11,16]: 2505 momenta = FourMomentum(part) 2506 tau1 += momenta 2507 elif part.pid == 15: 2508 tau2 += FourMomentum(part) 2509 2510 if abs((mtau-tau2.mass())/wtau)<1e6 and tau2.mass() >1: 2511 data.append((tau1.mass()-mtau)/wtau) 2512 data2.append((tau2.mass()-mtau)/wtau) 2513 gs1 = gridspec.GridSpec(2, 1, height_ratios=[5,1]) 2514 gs1.update(wspace=0, hspace=0) # set the spacing between axes. 2515 ax = plt.subplot(gs1[0]) 2516 2517 n, bins, patches = ax.hist(data2, nbins, histtype='step', label='original') 2518 n2, bins2, patches2 = ax.hist(data, bins=bins, histtype='step',label='reconstructed') 2519 import cmath 2520 2521 breit = lambda m : math.sqrt(4*math.pi)*1/(((m)**2-mtau**2)**2+(mtau*wtau)**2)*wtau 2522 2523 data3 = [breit(mtau + x*wtau)*wtau*16867622.6624*50 for x in bins] 2524 2525 ax.plot(bins, data3,label='breit-wigner') 2526 # add the legend 2527 ax.legend() 2528 # add on the right program tag 2529 ax_c = ax.twinx() 2530 ax_c.set_ylabel('MadGraph5_aMC@NLO') 2531 ax_c.yaxis.set_label_coords(1.01, 0.25) 2532 ax_c.set_yticks(ax.get_yticks()) 2533 ax_c.set_yticklabels([]) 2534 2535 plt.title('invariant mass of tau LHE/reconstructed') 2536 plt.axis('on') 2537 ax.set_xticklabels([]) 2538 # ratio plot 2539 ax = plt.subplot(gs1[1]) 2540 data4 = [n[i]/(data3[i]) for i in range(nbins)] 2541 ax.plot(bins, data4 + [0] , 'b') 2542 data4 = [n2[i]/(data3[i]) for i in range(nbins)] 2543 ax.plot(bins, data4 + [0] , 'g') 2544 ax.set_ylim([0,2]) 2545 #remove last y tick to avoid overlap with above plot: 2546 tick = ax.get_yticks() 2547 ax.set_yticks(tick[:-1]) 2548 2549 2550 plt.axis('on') 2551 plt.xlabel('(M - Mtau)/Wtau') 2552 plt.show() 2553