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