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

Source Code for Module madgraph.various.systematics

   1  ################################################################################ 
   2  # 
   3  # Copyright (c) 2016 The MadGraph5_aMC@NLO Development team and Contributors 
   4  # 
   5  # This file is a part of the MadGraph5_aMC@NLO project, an application which  
   6  # automatically generates Feynman diagrams and matrix elements for arbitrary 
   7  # high-energy processes in the Standard Model and beyond. 
   8  # 
   9  # It is subject to the MadGraph5_aMC@NLO license which should accompany this  
  10  # distribution. 
  11  # 
  12  # For more information, visit madgraph.phys.ucl.ac.be and amcatnlo.web.cern.ch 
  13  # 
  14  ################################################################################ 
  15  from __future__ import division 
  16   
  17  if __name__ == "__main__": 
  18      import sys 
  19      import os 
  20      root = os.path.dirname(__file__) 
  21      if os.path.basename(root) == 'internal': 
  22          sys.path.append(os.path.dirname(root)) 
  23      else: 
  24          sys.path.append(os.path.dirname(os.path.dirname(root))) 
  25           
  26  import lhe_parser 
  27  import banner 
  28  import banner as banner_mod 
  29  import itertools 
  30  import misc 
  31  import math 
  32  import os  
  33  import re 
  34  import sys 
  35  import time 
  36  import StringIO 
  37   
  38  pjoin = os.path.join 
  39  root = os.path.dirname(__file__) 
  40   
41 -class SystematicsError(Exception):
42 pass
43
44 -class Systematics(object):
45
46 - def __init__(self, input_file, output_file, 47 start_event=0, stop_event=sys.maxint, write_banner=False, 48 mur=[0.5,1,2], 49 muf=[0.5,1,2], 50 alps=[1], 51 pdf='errorset', #[(id, subset)] 52 dyn=[-1,1,2,3,4], 53 together=[('mur', 'muf', 'dyn')], 54 remove_wgts=[], 55 keep_wgts=[], 56 start_id=None, 57 lhapdf_config=misc.which('lhapdf-config'), 58 log=lambda x: sys.stdout.write(str(x)+'\n'), 59 only_beam=False, 60 ion_scaling=True, 61 weight_format=None, 62 weight_info=None, 63 ):
64 65 66 # INPUT/OUTPUT FILE 67 if isinstance(input_file, str): 68 self.input = lhe_parser.EventFile(input_file) 69 else: 70 self.input = input_file 71 self.output_path = output_file 72 self.weight_format = weight_format 73 self.weight_info_format = weight_info 74 if output_file != None: 75 if isinstance(output_file, str): 76 if output_file == input_file: 77 directory,name = os.path.split(output_file) 78 new_name = pjoin(directory, '.tmp_'+name) 79 self.output = lhe_parser.EventFile(new_name, 'w') 80 else: 81 self.output = lhe_parser.EventFile(output_file, 'w') 82 else: 83 self.output = output_file 84 self.log = log 85 86 #get some information from the run_card. 87 self.banner = banner_mod.Banner(self.input.banner) 88 self.force_write_banner = bool(write_banner) 89 self.orig_dyn = self.banner.get('run_card', 'dynamical_scale_choice') 90 self.orig_pdf = self.banner.run_card.get_lhapdf_id() 91 matching_mode = self.banner.get('run_card', 'ickkw') 92 93 #check for beam 94 beam1, beam2 = self.banner.get_pdg_beam() 95 if abs(beam1) != 2212 and abs(beam2) != 2212: 96 self.b1 = 0 97 self.b2 = 0 98 pdf = 'central' 99 #raise SystematicsError, 'can only reweight proton beam' 100 elif abs(beam1) != 2212: 101 self.b1 = 0 102 self.b2 = beam2//2212 103 elif abs(beam2) != 2212: 104 self.b1 = beam1//2212 105 self.b2 = 0 106 else: 107 self.b1 = beam1//2212 108 self.b2 = beam2//2212 109 110 self.orig_ion_pdf = False 111 self.ion_scaling = ion_scaling 112 self.only_beam = only_beam 113 if isinstance(self.banner.run_card, banner_mod.RunCardLO): 114 self.is_lo = True 115 if not self.banner.run_card['use_syst']: 116 raise SystematicsError, 'The events have not been generated with use_syst=True. Cannot evaluate systematics error on these events.' 117 118 if self.banner.run_card['nb_neutron1'] != 0 or \ 119 self.banner.run_card['nb_neutron2'] != 0 or \ 120 self.banner.run_card['nb_proton1'] != 1 or \ 121 self.banner.run_card['nb_proton2'] != 1: 122 self.orig_ion_pdf = True 123 else: 124 self.is_lo = False 125 if not self.banner.run_card['store_rwgt_info']: 126 raise SystematicsError, 'The events have not been generated with store_rwgt_info=True. Cannot evaluate systematics error on these events.' 127 128 # MUR/MUF/ALPS PARSING 129 if isinstance(mur, str): 130 mur = mur.split(',') 131 self.mur=[float(i) for i in mur] 132 if isinstance(muf, str): 133 muf = muf.split(',') 134 self.muf=[float(i) for i in muf] 135 136 if isinstance(alps, str): 137 alps = alps.split(',') 138 self.alps=[float(i) for i in alps] 139 140 # DYNAMICAL SCALE PARSING + together 141 if isinstance(dyn, str): 142 dyn = dyn.split(',') 143 self.dyn=[int(i) for i in dyn] 144 # For FxFx only mode -1 makes sense 145 if matching_mode == 3: 146 self.dyn = [-1] 147 # avoid sqrts at NLO if ISR is possible 148 if 4 in self.dyn and self.b1 and self.b2 and not self.is_lo: 149 self.dyn.remove(4) 150 151 if isinstance(together, str): 152 self.together = together.split(',') 153 else: 154 self.together = together 155 156 # START/STOP EVENT 157 self.start_event=int(start_event) 158 self.stop_event=int(stop_event) 159 if start_event != 0: 160 self.log( "#starting from event #%s" % start_event) 161 if stop_event != sys.maxint: 162 self.log( "#stopping at event #%s" % stop_event) 163 164 # LHAPDF set 165 if isinstance(lhapdf_config, list): 166 lhapdf_config = lhapdf_config[0] 167 lhapdf = misc.import_python_lhapdf(lhapdf_config) 168 if not lhapdf: 169 log('fail to load lhapdf: doe not perform systematics') 170 return 171 lhapdf.setVerbosity(0) 172 self.pdfsets = {} 173 if isinstance(pdf, str): 174 pdf = pdf.split(',') 175 176 if isinstance(pdf,list) and isinstance(pdf[0],(str,int)): 177 self.pdf = [] 178 for data in pdf: 179 if data == 'errorset': 180 data = '%s' % self.orig_pdf 181 if data == 'central': 182 data = '%s@0' % self.orig_pdf 183 if '@' in data: 184 #one particular dataset 185 name, arg = data.rsplit('@',1) 186 if int(arg) == 0: 187 if name.isdigit(): 188 self.pdf.append(lhapdf.mkPDF(int(name))) 189 else: 190 self.pdf.append(lhapdf.mkPDF(name)) 191 elif name.isdigit(): 192 try: 193 self.pdf.append(lhapdf.mkPDF(int(name)+int(arg))) 194 except: 195 raise Exception, 'Individual error sets need to be called with LHAPDF NAME not with LHAGLUE NUMBER' 196 else: 197 self.pdf.append(lhapdf.mkPDF(name, int(arg))) 198 else: 199 if data.isdigit(): 200 pdfset = lhapdf.mkPDF(int(data)).set() 201 else: 202 pdfset = lhapdf.mkPDF(data).set() 203 self.pdfsets[pdfset.lhapdfID] = pdfset 204 self.pdf += pdfset.mkPDFs() 205 else: 206 self.pdf = pdf 207 208 for p in self.pdf: 209 if p.lhapdfID == self.orig_pdf: 210 self.orig_pdf = p 211 break 212 else: 213 self.orig_pdf = lhapdf.mkPDF(self.orig_pdf) 214 if not self.b1 == 0 == self.b2: 215 self.log( "# events generated with PDF: %s (%s)" %(self.orig_pdf.set().name,self.orig_pdf.lhapdfID )) 216 # create all the function that need to be called 217 self.get_all_fct() # define self.fcts and self.args 218 219 # For e+/e- type of collision initialise the running of alps 220 if self.b1 == 0 == self.b2: 221 try: 222 from models.model_reader import Alphas_Runner 223 except ImportError: 224 root_path = pjoin(root, os.pardir, os.pardir) 225 try: 226 import internal.madevent_interface as me_int 227 cmd = me_int.MadEventCmd(root_path,force_run=True) 228 except ImportError: 229 import internal.amcnlo_run_interface as me_int 230 cmd = me_int.Cmd(root_path,force_run=True) 231 if 'mg5_path' in cmd.options and cmd.options['mg5_path']: 232 sys.path.append(cmd.options['mg5_path']) 233 from models.model_reader import Alphas_Runner 234 235 if not hasattr(self.banner, 'param_card'): 236 param_card = self.banner.charge_card('param_card') 237 else: 238 param_card = self.banner.param_card 239 240 asmz = param_card.get_value('sminputs', 3, 0.13) 241 nloop =2 242 zmass = param_card.get_value('mass', 23, 91.188) 243 cmass = param_card.get_value('mass', 4, 1.4) 244 if cmass == 0: 245 cmass = 1.4 246 bmass = param_card.get_value('mass', 5, 4.7) 247 if bmass == 0: 248 bmass = 4.7 249 self.alpsrunner = Alphas_Runner(asmz, nloop, zmass, cmass, bmass) 250 251 # Store which weight to keep/removed 252 self.remove_wgts = [] 253 for id in remove_wgts: 254 if id == 'all': 255 self.remove_wgts = ['all'] 256 break 257 elif ',' in id: 258 min_value, max_value = [int(v) for v in id.split(',')] 259 self.remove_wgts += [i for i in range(min_value, max_value+1)] 260 else: 261 self.remove_wgts.append(id) 262 self.keep_wgts = [] 263 for id in keep_wgts: 264 if id == 'all': 265 self.keep_wgts = ['all'] 266 break 267 elif ',' in id: 268 min_value, max_value = [int(v) for v in id.split(',')] 269 self.remove_wgts += [i for i in range(min_value, max_value+1)] 270 else: 271 self.remove_wgts.append(id) 272 273 # input to start the id in the weight 274 self.start_wgt_id = int(start_id[0]) if (start_id is not None) else None 275 self.has_wgts_pattern = False # tag to check if the pattern for removing
276 # the weights was computed already 277
278 - def is_wgt_kept(self, name):
279 """ determine if we have to keep/remove such weight """ 280 281 if 'all' in self.keep_wgts or not self.remove_wgts: 282 return True 283 284 #start by checking what we want to keep 285 if name in self.keep_wgts: 286 return True 287 288 # check for regular expression 289 if not self.has_wgts_pattern: 290 pat = r'|'.join(w for w in self.keep_wgts if any(letter in w for letter in '*?.([+\\')) 291 if pat: 292 self.keep_wgts_pattern = re.compile(pat) 293 else: 294 self.keep_wgts_pattern = None 295 pat = r'|'.join(w for w in self.remove_wgts if any(letter in w for letter in '*?.([+\\')) 296 if pat: 297 self.rm_wgts_pattern = re.compile(pat) 298 else: 299 self.rm_wgts_pattern = None 300 self.has_wgts_pattern=True 301 302 if self.keep_wgts_pattern and re.match(self.keep_wgts_pattern,name): 303 return True 304 305 #check what we want to remove 306 if 'all' in self.remove_wgts: 307 return False 308 elif name in self.remove_wgts: 309 return False 310 elif self.rm_wgts_pattern and re.match(self.rm_wgts_pattern, name): 311 return False 312 else: 313 return True
314
315 - def remove_old_wgts(self, event):
316 """remove the weight as requested by the user""" 317 318 rwgt_data = event.parse_reweight() 319 for name in rwgt_data.keys(): 320 if not self.is_wgt_kept(name): 321 del rwgt_data[name] 322 event.reweight_order.remove(name)
323 324
325 - def run(self, stdout=sys.stdout):
326 """ """ 327 start_time = time.time() 328 if self.start_event == 0 or self.force_write_banner: 329 lowest_id = self.write_banner(self.output) 330 else: 331 lowest_id = self.get_id() 332 333 ids = [self.get_wgt_name(*self.args[i][:5], cid=lowest_id+i) for i in range(len(self.args)-1)] 334 #ids = [lowest_id+i for i in range(len(self.args)-1)] 335 all_cross = [0 for i in range(len(self.args))] 336 337 self.input.parsing = False 338 for nb_event,event in enumerate(self.input): 339 if nb_event < self.start_event: 340 continue 341 elif nb_event == self.start_event: 342 self.input.parsing = True 343 event = lhe_parser.Event(event) 344 elif nb_event >= self.stop_event: 345 if self.force_write_banner: 346 self.output.write('</LesHouchesEvents>\n') 347 break 348 349 if self.is_lo: 350 if (nb_event-self.start_event)>=0 and (nb_event-self.start_event) % 2500 ==0: 351 self.log( '# currently at event %s [elapsed time: %.2g s]' % (nb_event, time.time()-start_time)) 352 else: 353 if (nb_event-self.start_event)>=0 and (nb_event-self.start_event) % 1000 ==0: 354 self.log( '# currently at event %i [elapsed time: %.2g s]' % (nb_event, time.time()-start_time)) 355 356 self.new_event() #re-init the caching of alphas/pdf 357 self.remove_old_wgts(event) 358 if self.is_lo: 359 wgts = [self.get_lo_wgt(event, *arg) for arg in self.args] 360 else: 361 wgts = [self.get_nlo_wgt(event, *arg) for arg in self.args] 362 363 if wgts[0] == 0: 364 print wgts 365 print event 366 raise Exception 367 368 wgt = [event.wgt*wgts[i]/wgts[0] for i in range(1,len(wgts))] 369 all_cross = [(all_cross[j] + event.wgt*wgts[j]/wgts[0]) for j in range(len(wgts))] 370 371 rwgt_data = event.parse_reweight() 372 rwgt_data.update(zip(ids, wgt)) 373 event.reweight_order += ids 374 # order the 375 self.output.write(str(event)) 376 else: 377 self.output.write('</LesHouchesEvents>\n') 378 self.output.close() 379 self.print_cross_sections(all_cross, min(nb_event,self.stop_event)-self.start_event+1, stdout) 380 381 if self.output.name != self.output_path: 382 #check problem for .gz missinf in output.name 383 if not os.path.exists(self.output.name) and os.path.exists('%s.gz' % self.output.name): 384 to_check = '%s.gz' % self.output.name 385 else: 386 to_check = self.output.name 387 388 if to_check != self.output_path: 389 if '%s.gz' % to_check == self.output_path: 390 misc.gzip(to_check) 391 else: 392 import shutil 393 shutil.move(to_check, self.output_path) 394 395 return all_cross
396
397 - def print_cross_sections(self, all_cross, nb_event, stdout):
398 """print the cross-section.""" 399 400 norm = self.banner.get('run_card', 'event_norm', default='sum') 401 #print "normalisation is ", norm 402 #print "nb_event is ", nb_event 403 404 max_scale, min_scale = 0,sys.maxint 405 max_alps, min_alps = 0, sys.maxint 406 max_dyn, min_dyn = 0,sys.maxint 407 pdfs = {} 408 dyns = {} # dyn : {'max': , 'min':} 409 410 if norm == 'sum': 411 norm = 1 412 elif norm in ['average', 'bias']: 413 norm = 1./nb_event 414 elif norm == 'unity': 415 norm = 1 416 417 all_cross = [c*norm for c in all_cross] 418 stdout.write("# mur\t\tmuf\t\talpsfact\tdynamical_scale\tpdf\t\tcross-section\n") 419 for i,arg in enumerate(self.args): 420 421 to_print = list(arg) 422 to_print[4] = to_print[4].lhapdfID 423 to_print.append(all_cross[i]) 424 to_report = [] 425 stdout.write('%s\t\t%s\t\t%s\t\t%s\t\t%s\t\t%s\n' % tuple(to_print)) 426 427 mur, muf, alps, dyn, pdf = arg[:5] 428 if pdf == self.orig_pdf and (dyn==self.orig_dyn or dyn==-1)\ 429 and (mur!=1 or muf!=1 or alps!=1): 430 max_scale = max(max_scale,all_cross[i]) 431 min_scale = min(min_scale,all_cross[i]) 432 if pdf == self.orig_pdf and mur==1 and muf==1 and \ 433 (dyn==self.orig_dyn or dyn==-1) and alps!=1: 434 max_alps = max(max_alps,all_cross[i]) 435 min_alps = min(min_alps,all_cross[i]) 436 437 if pdf == self.orig_pdf and mur==1 and muf==1 and alps==1: 438 max_dyn = max(max_dyn,all_cross[i]) 439 min_dyn = min(min_dyn,all_cross[i]) 440 441 if pdf == self.orig_pdf and (alps!=1 or mur!=1 or muf!=1) and \ 442 (dyn!=self.orig_dyn or dyn!=-1): 443 if dyn not in dyns: 444 dyns[dyn] = {'max':0, 'min':sys.maxint,'central':0} 445 curr = dyns[dyn] 446 curr['max'] = max(curr['max'],all_cross[i]) 447 curr['min'] = min(curr['min'],all_cross[i]) 448 if pdf == self.orig_pdf and (alps==1 and mur==1 and muf==1) and \ 449 (dyn!=self.orig_dyn or dyn!=-1): 450 if dyn not in dyns: 451 dyns[dyn] = {'max':0, 'min':sys.maxint,'central':all_cross[i]} 452 else: 453 dyns[dyn]['central'] = all_cross[i] 454 455 if alps==1 and mur==1 and muf==1 and (dyn==self.orig_dyn or dyn==-1): 456 pdfset = pdf.set() 457 if pdfset.lhapdfID in self.pdfsets: 458 if pdfset.lhapdfID not in pdfs : 459 pdfs[pdfset.lhapdfID] = [0] * pdfset.size 460 pdfs[pdfset.lhapdfID][pdf.memberID] = all_cross[i] 461 else: 462 to_report.append('# PDF %s : %s\n' % (pdf.lhapdfID, all_cross[i])) 463 464 stdout.write('\n') 465 466 resume = StringIO.StringIO() 467 468 resume.write( '#***************************************************************************\n') 469 resume.write( "#\n") 470 resume.write( '# original cross-section: %s\n' % all_cross[0]) 471 if max_scale: 472 resume.write( '# scale variation: +%2.3g%% -%2.3g%%\n' % ((max_scale-all_cross[0])/all_cross[0]*100,(all_cross[0]-min_scale)/all_cross[0]*100)) 473 if max_alps: 474 resume.write( '# emission scale variation: +%2.3g%% -%2.3g%%\n' % ((max_alps-all_cross[0])/all_cross[0]*100,(max_alps-min_scale)/all_cross[0]*100)) 475 if max_dyn and (max_dyn!= all_cross[0] or min_dyn != all_cross[0]): 476 resume.write( '# central scheme variation: +%2.3g%% -%2.3g%%\n' % ((max_dyn-all_cross[0])/all_cross[0]*100,(all_cross[0]-min_dyn)/all_cross[0]*100)) 477 if self.orig_pdf.lhapdfID in pdfs: 478 lhapdfid = self.orig_pdf.lhapdfID 479 values = pdfs[lhapdfid] 480 pdfset = self.pdfsets[lhapdfid] 481 try: 482 pdferr = pdfset.uncertainty(values) 483 except RuntimeError: 484 resume.write( '# PDF variation: missing combination\n') 485 else: 486 resume.write( '# PDF variation: +%2.3g%% -%2.3g%%\n' % (pdferr.errplus*100/all_cross[0], pdferr.errminus*100/all_cross[0])) 487 # report error/central not directly linked to the central 488 resume.write( "#\n") 489 for lhapdfid,values in pdfs.items(): 490 if lhapdfid == self.orig_pdf.lhapdfID: 491 continue 492 if len(values) == 1 : 493 continue 494 pdfset = self.pdfsets[lhapdfid] 495 496 if pdfset.errorType == 'unknown' : 497 # Don't know how to determine uncertainty for 'unknown' errorType : 498 # File "lhapdf.pyx", line 329, in lhapdf.PDFSet.uncertainty (lhapdf.cpp:6621) 499 # RuntimeError: "ErrorType: unknown" not supported by LHAPDF::PDFSet::uncertainty. 500 continue 501 try: 502 pdferr = pdfset.uncertainty(values) 503 except RuntimeError: 504 # the same error can happend to some other type of error like custom. 505 pass 506 else: 507 resume.write( '#PDF %s: %g +%2.3g%% -%2.3g%%\n' % (pdfset.name, pdferr.central,pdferr.errplus*100/all_cross[0], pdferr.errminus*100/all_cross[0])) 508 509 dyn_name = {1: '\sum ET', 2:'\sum\sqrt{m^2+pt^2}', 3:'0.5 \sum\sqrt{m^2+pt^2}',4:'\sqrt{\hat s}' } 510 for key, curr in dyns.items(): 511 if key ==-1: 512 continue 513 central, maxvalue, minvalue = curr['central'], curr['max'], curr['min'] 514 if central == 0: 515 continue 516 if maxvalue == 0: 517 resume.write("# dynamical scheme # %s : %g # %s\n" %(key, central, dyn_name[key])) 518 else: 519 resume.write("# dynamical scheme # %s : %g +%2.3g%% -%2.3g%% # %s\n" %(key, central, (maxvalue-central)/central*100,(central-minvalue)/central*100, dyn_name[key])) 520 521 resume.write('\n'.join(to_report)) 522 523 resume.write( '#***************************************************************************\n') 524 525 stdout.write(resume.getvalue()) 526 self.log(resume.getvalue())
527 528
529 - def write_banner(self, fsock):
530 """create the new banner with the information of the weight""" 531 532 cid = self.get_id() 533 lowest_id = cid 534 535 in_scale = False 536 in_pdf = False 537 in_alps = False 538 539 text = '' 540 541 default = self.args[0] 542 for arg in self.args[1:]: 543 mur, muf, alps, dyn, pdf = arg[:5] 544 if pdf == self.orig_pdf and alps ==1 and (mur!=1 or muf!=1 or dyn!=-1): 545 if not in_scale: 546 text += "<weightgroup name=\"Central scale variation\" combine=\"envelope\">\n" 547 in_scale=True 548 elif in_scale: 549 if (pdf == self.orig_pdf and alps ==1) and arg != default: 550 pass 551 else: 552 text += "</weightgroup> # scale\n" 553 in_scale = False 554 555 if pdf == self.orig_pdf and mur == muf == 1 and dyn==-1 and alps!=1: 556 if not in_alps: 557 text += "<weightgroup name=\"Emission scale variation\" combine=\"envelope\">\n" 558 in_alps=True 559 elif in_alps: 560 text += "</weightgroup> # ALPS\n" 561 in_alps=False 562 563 if mur == muf == 1 and dyn==-1 and alps ==1: 564 565 if pdf.lhapdfID in self.pdfsets: 566 if in_pdf: 567 text += "</weightgroup> # PDFSET -> PDFSET\n" 568 pdfset = self.pdfsets[pdf.lhapdfID] 569 descrip = pdfset.description.replace('=>',';').replace('>','.gt.').replace('<','.lt.') 570 text +="<weightgroup name=\"%s\" combine=\"%s\"> # %s: %s\n" %\ 571 (pdfset.name, pdfset.errorType,pdfset.lhapdfID, descrip) 572 in_pdf=pdf.lhapdfID 573 elif pdf.memberID == 0 and (pdf.lhapdfID - pdf.memberID) in self.pdfsets: 574 if in_pdf: 575 text += "</weightgroup> # PDFSET -> PDFSET\n" 576 pdfset = self.pdfsets[pdf.lhapdfID - 1] 577 descrip = pdfset.description.replace('=>',';').replace('>','.gt.').replace('<','.lt.') 578 text +="<weightgroup name=\"%s\" combine=\"%s\"> # %s: %s\n" %\ 579 (pdfset.name, pdfset.errorType,pdfset.lhapdfID, descrip) 580 in_pdf=pdfset.lhapdfID 581 elif in_pdf and pdf.lhapdfID - pdf.memberID != in_pdf: 582 text += "</weightgroup> # PDFSET -> PDF\n" 583 in_pdf = False 584 elif in_pdf: 585 text += "</weightgroup> PDF \n" 586 in_pdf=False 587 588 589 590 tag, info = '','' 591 if mur!=1.: 592 tag += 'MUR="%s" ' % mur 593 info += 'MUR=%s ' % mur 594 else: 595 tag += 'MUR="%s" ' % mur 596 if muf!=1.: 597 tag += 'MUF="%s" ' % muf 598 info += 'MUF=%s ' % muf 599 else: 600 tag += 'MUF="%s" ' % muf 601 602 if alps!=1.: 603 tag += 'ALPSFACT="%s" ' % alps 604 info += 'alpsfact=%s ' % alps 605 if dyn!=-1.: 606 tag += 'DYN_SCALE="%s" ' % dyn 607 info += 'dyn_scale_choice=%s ' % {1:'sum pt', 2:'HT',3:'HT/2',4:'sqrts'}[dyn] 608 609 if pdf != self.orig_pdf: 610 tag += 'PDF="%s" ' % pdf.lhapdfID 611 info += 'PDF=%s MemberID=%s' % (pdf.lhapdfID-pdf.memberID, pdf.memberID) 612 else: 613 tag += 'PDF="%s" ' % pdf.lhapdfID 614 615 wgt_name = self.get_wgt_name(mur, muf, alps, dyn, pdf, cid) 616 tag = self.get_wgt_tag(mur, muf, alps, dyn, pdf, cid) 617 info = self.get_wgt_info(mur, muf, alps, dyn, pdf, cid) 618 text +='<weight id="%s" %s> %s </weight>\n' % (wgt_name, tag, info) 619 cid+=1 620 621 if in_scale or in_alps or in_pdf: 622 text += "</weightgroup>\n" 623 624 if 'initrwgt' in self.banner: 625 if not self.remove_wgts: 626 self.banner['initrwgt'] += text 627 else: 628 # remove the line which correspond to removed weight 629 # removed empty group. 630 wgt_in_group=0 631 tmp_group_txt =[] 632 out =[] 633 keep_last = False 634 for line in self.banner['initrwgt'].split('\n'): 635 sline = line.strip() 636 if sline.startswith('</weightgroup'): 637 if wgt_in_group: 638 out += tmp_group_txt 639 out.append('</weightgroup>') 640 if '<weightgroup' in line: 641 wgt_in_group=0 642 tmp_group_txt = [line[line.index('<weightgroup'):]] 643 elif sline.startswith('<weightgroup'): 644 wgt_in_group=0 645 tmp_group_txt = [line] 646 elif sline.startswith('<weight'): 647 name = re.findall(r'\bid=[\'\"]([^\'\"]*)[\'\"]', sline) 648 if self.is_wgt_kept(name[0]): 649 tmp_group_txt.append(line) 650 keep_last = True 651 wgt_in_group +=1 652 else: 653 keep_last = False 654 elif keep_last: 655 tmp_group_txt.append(line) 656 out.append(text) 657 self.banner['initrwgt'] = '\n'.join(out) 658 else: 659 self.banner['initrwgt'] = text 660 661 662 self.banner.write(self.output, close_tag=False) 663 664 return lowest_id
665
666 - def get_wgt_name(self, mur, muf, alps, dyn, pdf, cid=0):
667 668 if self.weight_format: 669 wgt_name = self.weight_format[0] % {'mur': mur, 'muf':muf, 'alps': alps, 'pdf':pdf.lhapdfID, 'dyn':dyn, 'id': cid} 670 else: 671 wgt_name = cid 672 return wgt_name
673
674 - def get_wgt_info(self, mur, muf, alps, dyn, pdf, cid=0):
675 676 if self.weight_info_format: 677 info = self.weight_info_format[0] % {'mur': mur, 'muf':muf, 'alps': alps, 'pdf':pdf.lhapdfID, 'dyn':dyn, 'id': cid, 's':' ', 'n':'\n'} 678 else: 679 info = '' 680 if mur!=1.: 681 info += 'MUR=%s ' % mur 682 if muf!=1.: 683 info += 'MUF=%s ' % muf 684 if alps!=1.: 685 info += 'alpsfact=%s ' % alps 686 if dyn!=-1.: 687 info += 'dyn_scale_choice=%s ' % {1:'sum pt', 2:'HT',3:'HT/2',4:'sqrts'}[dyn] 688 if pdf != self.orig_pdf: 689 info += 'PDF=%s MemberID=%s' % (pdf.lhapdfID-pdf.memberID, pdf.memberID) 690 691 return info
692
693 - def get_wgt_tag (self, mur, muf, alps, dyn, pdf, cid=0):
694 tags = [] 695 tags.append('MUR="%s" ' % mur) 696 tags.append('MUF="%s" ' % muf) 697 if alps!=1.: 698 tags.append('ALPSFACT="%s" ' % alps) 699 if dyn!=-1.: 700 tags.append('DYN_SCALE="%s" ' % dyn) 701 tags.append('PDF="%s" ' % pdf.lhapdfID) 702 return " ".join(tags)
703 704
705 - def get_id(self):
706 707 if self.start_wgt_id is not None: 708 return int(self.start_wgt_id) 709 710 if 'initrwgt' in self.banner: 711 pattern = re.compile('<weight id=(?:\'|\")([_\w]+)(?:\'|\")', re.S+re.I+re.M) 712 matches = pattern.findall(self.banner['initrwgt']) 713 matches.append('0') #ensure to have a valid entry for the max 714 return max([int(wid) for wid in matches if wid.isdigit()])+1 715 else: 716 return 1
717 718
719 - def get_all_fct(self):
720 721 all_args = [] 722 default = [1.,1.,1.,-1,self.orig_pdf] 723 #all_args.append(default) 724 pos = {'mur':0, 'muf':1, 'alps':2, 'dyn':3, 'pdf':4} 725 done = set() 726 for one_block in self.together: 727 for name in one_block: 728 done.add(name) 729 for together in itertools.product(*[getattr(self,name) for name in one_block]): 730 new_args = list(default) 731 for name,value in zip(one_block, together): 732 new_args[pos[name]] = value 733 all_args.append(new_args) 734 for name in pos: 735 if name in done: 736 continue 737 for value in getattr(self, name): 738 new_args = list(default) 739 new_args[pos[name]] = value 740 all_args.append(new_args) 741 742 self.args = [default] + [arg for arg in all_args if arg!= default] 743 744 # add the default before the pdf scan to have a full grouping 745 pdfplusone = [pdf for pdf in self.pdf if pdf.lhapdfID == self.orig_pdf.lhapdfID+1] 746 if pdfplusone: 747 pdfplusone = default[:-1] + [pdfplusone[0]] 748 index = self.args.index(pdfplusone) 749 self.args.insert(index, default) 750 751 self.log( "#Will Compute %s weights per event." % (len(self.args)-1)) 752 return
753
754 - def new_event(self):
755 self.alphas = {} 756 self.pdfQ2 = {}
757 758
759 - def get_pdfQ(self, pdf, pdg, x, scale, beam=1):
760 761 if pdg in [-21,-22]: 762 pdg = abs(pdg) 763 elif pdg == 0: 764 return 1 765 766 if self.only_beam and self.only_beam!= beam and pdf.lhapdfID != self.orig_pdf: 767 return self.getpdfQ(self.pdfsets[self.orig_pdf], pdg, x, scale, beam) 768 769 if self.orig_ion_pdf and (self.ion_scaling or pdf.lhapdfID == self.orig_pdf): 770 nb_p = self.banner.run_card["nb_proton%s" % beam] 771 nb_n = self.banner.run_card["nb_neutron%s" % beam] 772 773 774 if pdg in [1,2]: 775 pdf1 = pdf.xfxQ(1, x, scale)/x 776 pdf2 = pdf.xfxQ(2, x, scale)/x 777 if pdg == 1: 778 f = nb_p * pdf1 + nb_n * pdf2 779 else: 780 f = nb_p * pdf2 + nb_n * pdf1 781 elif pdg in [-1,-2]: 782 pdf1 = pdf.xfxQ(-1, x, scale)/x 783 pdf2 = pdf.xfxQ(-2, x, scale)/x 784 if pdg == -1: 785 f = nb_p * pdf1 + nb_n * pdf2 786 else: 787 f = nb_p * pdf2 + nb_n * pdf1 788 else: 789 f = (nb_p + nb_n) * pdf.xfxQ(pdg, x, scale)/x 790 791 f = f * (nb_p+nb_n) 792 else: 793 f = pdf.xfxQ(pdg, x, scale)/x 794 # if f == 0 and pdf.memberID ==0: 795 # pdfset = pdf.set() 796 # allnumber= [p.xfxQ(pdg, x, scale) for p in pdfset.mkPDFs()] 797 # f = pdfset.uncertainty(allnumber).central /x 798 return f
799
800 - def get_pdfQ2(self, pdf, pdg, x, scale, beam=1):
801 802 if pdg in [-21,-22]: 803 pdg = abs(pdg) 804 elif pdg == 0: 805 return 1 806 807 if (pdf, pdg,x,scale, beam) in self.pdfQ2: 808 return self.pdfQ2[(pdf, pdg,x,scale,beam)] 809 810 if self.orig_ion_pdf and (self.ion_scaling or pdf.lhapdfID == self.orig_pdf): 811 nb_p = self.banner.run_card["nb_proton%s" % beam] 812 nb_n = self.banner.run_card["nb_neutron%s" % beam] 813 814 815 if pdg in [1,2]: 816 pdf1 = pdf.xfxQ2(1, x, scale)/x 817 pdf2 = pdf.xfxQ2(2, x, scale)/x 818 if pdg == 1: 819 f = nb_p * pdf1 + nb_n * pdf2 820 else: 821 f = nb_p * pdf2 + nb_n * pdf1 822 elif pdg in [-1,-2]: 823 pdf1 = pdf.xfxQ2(-1, x, scale)/x 824 pdf2 = pdf.xfxQ2(-2, x, scale)/x 825 if pdg == -1: 826 f = nb_p * pdf1 + nb_n * pdf2 827 else: 828 f = nb_p * pdf2 + nb_n * pdf1 829 else: 830 f = (nb_p + nb_n) * pdf.xfxQ2(pdg, x, scale)/x 831 832 f = f * (nb_p+nb_n) 833 else: 834 f = pdf.xfxQ2(pdg, x, scale)/x 835 self.pdfQ2[(pdf, pdg,x,scale,beam)] = f 836 return f 837 838 839 840 #one method to handle the nnpd2.3 problem -> now just move to central 841 if f == 0 and pdf.memberID ==0: 842 # to avoid problem with nnpdf2.3 in lhapdf6.1.6 843 #print 'central pdf returns 0', pdg, x, scale 844 #print self.pdfsets 845 pdfset = pdf.set() 846 allnumber= [0] + [self.get_pdfQ2(p, pdg, x, scale) for p in pdfset.mkPDFs()[1:]] 847 f = pdfset.uncertainty(allnumber).central 848 self.pdfQ2[(pdf, pdg,x,scale)] = f 849 return f
850
851 - def get_lo_wgt(self,event, Dmur, Dmuf, Dalps, dyn, pdf):
852 """ 853 pdf is a lhapdf object!""" 854 855 loinfo = event.parse_lo_weight() 856 857 if dyn == -1: 858 mur = loinfo['ren_scale'] 859 muf1 = loinfo['pdf_q1'][-1] 860 muf2 = loinfo['pdf_q2'][-1] 861 else: 862 if dyn == 1: 863 mur = event.get_et_scale(1.) 864 elif dyn == 2: 865 mur = event.get_ht_scale(1.) 866 elif dyn == 3: 867 mur = event.get_ht_scale(0.5) 868 elif dyn == 4: 869 mur = event.get_sqrts_scale(1.) 870 if math.isnan(mur): 871 return mur 872 muf1 = mur 873 muf2 = mur 874 loinfo = dict(loinfo) 875 loinfo['pdf_q1'] = loinfo['pdf_q1'] [:-1] + [mur] 876 loinfo['pdf_q2'] = loinfo['pdf_q2'] [:-1] + [mur] 877 878 879 # MUR part 880 if self.b1 == 0 == self.b2: 881 if loinfo['n_qcd'] != 0: 882 wgt = self.alpsrunner(Dmur*mur)**loinfo['n_qcd'] 883 else: 884 wgt = 1.0 885 else: 886 wgt = pdf.alphasQ(Dmur*mur)**loinfo['n_qcd'] 887 # MUF/PDF part 888 wgt *= self.get_pdfQ(pdf, self.b1*loinfo['pdf_pdg_code1'][-1], loinfo['pdf_x1'][-1], Dmuf*muf1, beam=1) 889 wgt *= self.get_pdfQ(pdf, self.b2*loinfo['pdf_pdg_code2'][-1], loinfo['pdf_x2'][-1], Dmuf*muf2, beam=2) 890 891 for scale in loinfo['asrwt']: 892 if self.b1 == 0 == self.b2: 893 wgt = self.alpsrunner(Dalps*scale) 894 else: 895 wgt *= pdf.alphasQ(Dalps*scale) 896 897 # ALS part 898 for i in range(loinfo['n_pdfrw1']-1): 899 scale = min(Dalps*loinfo['pdf_q1'][i], Dmuf*muf1) 900 wgt *= self.get_pdfQ(pdf, self.b1*loinfo['pdf_pdg_code1'][i], loinfo['pdf_x1'][i], scale, beam=1) 901 wgt /= self.get_pdfQ(pdf, self.b1*loinfo['pdf_pdg_code1'][i], loinfo['pdf_x1'][i+1], scale, beam=1) 902 903 for i in range(loinfo['n_pdfrw2']-1): 904 scale = min(Dalps*loinfo['pdf_q2'][i], Dmuf*muf2) 905 wgt *= self.get_pdfQ(pdf, self.b2*loinfo['pdf_pdg_code2'][i], loinfo['pdf_x2'][i], scale, beam=2) 906 wgt /= self.get_pdfQ(pdf, self.b2*loinfo['pdf_pdg_code2'][i], loinfo['pdf_x2'][i+1], scale, beam=2) 907 908 return wgt
909
910 - def get_nlo_wgt(self,event, Dmur, Dmuf, Dalps, dyn, pdf):
911 """return the new weight for NLO event --with weight information-- """ 912 913 wgt = 0 914 nloinfo = event.parse_nlo_weight(real_type=(1,11,12,13)) 915 for cevent in nloinfo.cevents: 916 if dyn == 1: 917 mur2 = max(1.0, cevent.get_et_scale(1.)**2) 918 elif dyn == 2: 919 mur2 = max(1.0, cevent.get_ht_scale(1.)**2) 920 elif dyn == 3: 921 mur2 = max(1.0, cevent.get_ht_scale(0.5)**2) 922 elif dyn == 4: 923 mur2 = cevent.get_sqrts_scale(event,1)**2 924 else: 925 mur2 = 0 926 muf2 = mur2 927 928 for onewgt in cevent.wgts: 929 if not __debug__ and (dyn== -1 and Dmur==1 and Dmuf==1 and pdf==self.orig_pdf): 930 wgt += onewgt.ref_wgt 931 continue 932 933 if dyn == -1: 934 mur2 = onewgt.scales2[1] 935 muf2 = onewgt.scales2[2] 936 Q2 = onewgt.scales2[0] # Ellis-Sexton scale 937 938 wgtpdf = self.get_pdfQ2(pdf, self.b1*onewgt.pdgs[0], onewgt.bjks[0], 939 Dmuf**2 * muf2) 940 wgtpdf *= self.get_pdfQ2(pdf, self.b2*onewgt.pdgs[1], onewgt.bjks[1], 941 Dmuf**2 * muf2) 942 943 tmp = onewgt.pwgt[0] 944 tmp += onewgt.pwgt[1] * math.log(Dmur**2 * mur2/ Q2) 945 tmp += onewgt.pwgt[2] * math.log(Dmuf**2 * muf2/ Q2) 946 947 if self.b1 == 0 == self.b2: 948 alps = self.alpsrunner(Dmur*math.sqrt(mur2)) 949 else: 950 alps = pdf.alphasQ2(Dmur**2*mur2) 951 952 tmp *= math.sqrt(4*math.pi*alps)**onewgt.qcdpower 953 954 if wgtpdf == 0: #happens for nn23pdf due to wrong set in lhapdf 955 key = (self.b1*onewgt.pdgs[0], self.b2*onewgt.pdgs[1], onewgt.bjks[0],onewgt.bjks[1], muf2) 956 if dyn== -1 and Dmuf==1 and Dmur==1 and pdf==self.orig_pdf: 957 wgtpdf = onewgt.ref_wgt / tmp 958 self.pdfQ2[key] = wgtpdf 959 elif key in self.pdfQ2: 960 wgtpdf = self.pdfQ2[key] 961 else: 962 # real zero! 963 wgtpdf = 0 964 965 tmp *= wgtpdf 966 wgt += tmp 967 968 969 if __debug__ and dyn== -1 and Dmur==1 and Dmuf==1 and pdf==self.orig_pdf: 970 if not misc.equal(tmp, onewgt.ref_wgt, sig_fig=2): 971 misc.sprint(tmp, onewgt.ref_wgt, (tmp-onewgt.ref_wgt)/tmp) 972 misc.sprint(onewgt) 973 misc.sprint(cevent) 974 misc.sprint(mur2,muf2) 975 raise Exception, 'not enough agreement between stored value and computed one' 976 977 return wgt
978 979
980 -def call_systematics(args, result=sys.stdout, running=True, 981 log=lambda x:sys.stdout.write(str(x)+'\n')):
982 """calling systematics from a list of arguments""" 983 984 985 input, output = args[0:2] 986 987 start_opts = 2 988 if output and output.startswith('-'): 989 start_opts = 1 990 output = input 991 992 opts = {} 993 for arg in args[start_opts:]: 994 if '=' in arg: 995 key,values= arg.split('=') 996 key = key.replace('-','') 997 values = values.strip() 998 if values[0] in ["'",'"'] and values[-1]==values[0]: 999 values = values[1:-1] 1000 values = values.split(',') 1001 if key == 'together': 1002 if key in opts: 1003 opts[key].append(tuple(values)) 1004 else: 1005 opts[key]=[tuple(values)] 1006 elif key == 'result': 1007 result = open(values[0],'w') 1008 elif key in ['start_event', 'stop_event', 'only_beam']: 1009 opts[key] = banner_mod.ConfigFile.format_variable(values[0], int, key) 1010 elif key in ['write_banner', 'ion_scalling']: 1011 opts[key] = banner_mod.ConfigFile.format_variable(values[0], bool, key) 1012 else: 1013 if key in opts: 1014 opts[key] += values 1015 else: 1016 opts[key] = values 1017 else: 1018 raise SystematicsError, "unknow argument %s" % arg 1019 1020 #load run_card and extract parameter if needed. 1021 if 'from_card' in opts: 1022 if opts['from_card'] != ['internal']: 1023 card = banner.RunCard(opts['from_card'][0]) 1024 else: 1025 for i in range(10): 1026 try: 1027 lhe = lhe_parser.EventFile(input) 1028 break 1029 except OSError,error: 1030 print error 1031 time.sleep(15*(i+1)) 1032 else: 1033 raise 1034 1035 card = banner.RunCard(banner.Banner(lhe.banner)['mgruncard']) 1036 lhe.close() 1037 1038 if isinstance(card, banner.RunCardLO): 1039 # LO case 1040 if 'systematics_arguments' in card.user_set: 1041 return call_systematics([input, output] + card['systematics_arguments'] 1042 , result=result, running=running, 1043 log=log) 1044 1045 else: 1046 opts['mur'] = [float(x) for x in card['sys_scalefact'].split()] 1047 opts['muf'] = opts['mur'] 1048 if card['sys_alpsfact'] != 'None': 1049 opts['alps'] = [float(x) for x in card['sys_alpsfact'].split()] 1050 else: 1051 opts['alps'] = [1.0] 1052 opts['together'] = [('mur','muf','alps','dyn')] 1053 if '&&' in card['sys_pdf']: 1054 pdfs = card['sys_pdf'].split('&&') 1055 else: 1056 data = card['sys_pdf'].split() 1057 pdfs = [] 1058 for d in data: 1059 if not d.isdigit(): 1060 pdfs.append(d) 1061 elif int(d) > 500: 1062 pdfs.append(d) 1063 else: 1064 pdfs[-1] = '%s %s' % (pdfs[-1], d) 1065 1066 opts['dyn'] = [-1,1,2,3,4] 1067 opts['pdf'] = [] 1068 for pdf in pdfs: 1069 split = pdf.split() 1070 if len(split)==1: 1071 opts['pdf'].append('%s' %pdf) 1072 else: 1073 pdf,nb = split 1074 for i in range(int(nb)): 1075 opts['pdf'].append('%s@%s' % (pdf, i)) 1076 if not opts['pdf']: 1077 opts['pdf'] = 'central' 1078 else: 1079 #NLO case 1080 if 'systematics_arguments' in card.user_set: 1081 return call_systematics([input, output] + card['systematics_arguments'] 1082 , result=result, running=running, 1083 log=log) 1084 else: 1085 raise Exception 1086 del opts['from_card'] 1087 1088 1089 obj = Systematics(input, output, log=log, **opts) 1090 if running and obj: 1091 obj.run(result) 1092 return obj
1093 1094 if __name__ == "__main__": 1095 1096 sys_args = sys.argv[1:] 1097 for i, arg in enumerate(sys_args) : 1098 if arg.startswith('--lhapdf_config=') : 1099 lhapdf = misc.import_python_lhapdf(arg[len('--lhapdf_config='):]) 1100 del sys_args[i] 1101 break 1102 1103 if 'lhapdf' not in globals(): 1104 lhapdf = misc.import_python_lhapdf('lhapdf-config') 1105 1106 if not lhapdf: 1107 sys.exit('Can not run systematics since can not link python to lhapdf, specify --lhapdf_config=') 1108 call_systematics(sys_args) 1109