Package madgraph :: Package interface :: Module reweight_interface
[hide private]
[frames] | no frames]

Source Code for Module madgraph.interface.reweight_interface

   1  ################################################################################ 
   2  # 
   3  # Copyright (c) 2009 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  """ Command interface for Re-Weighting """ 
  16  from __future__ import division 
  17  import difflib 
  18  import logging 
  19  import math 
  20  import os 
  21  import re 
  22  import shutil 
  23  import sys 
  24  import tempfile 
  25  import time 
  26  import subprocess 
  27  from subprocess import Popen, PIPE, STDOUT 
  28   
  29   
  30  pjoin = os.path.join 
  31   
  32  import madgraph.interface.extended_cmd as extended_cmd 
  33  import madgraph.interface.madgraph_interface as mg_interface 
  34  import madgraph.interface.master_interface as master_interface 
  35  import madgraph.interface.common_run_interface as common_run_interface 
  36  import madgraph.interface.madevent_interface as madevent_interface 
  37  import madgraph.iolibs.files as files 
  38  #import MadSpin.interface_madspin as madspin_interface 
  39  import madgraph.various.misc as misc 
  40  import madgraph.various.banner as banner 
  41  import madgraph.various.lhe_parser as lhe_parser 
  42  import madgraph.various.combine_plots as combine_plots 
  43  import madgraph.various.cluster as cluster 
  44  import madgraph.fks.fks_common as fks_common 
  45  import madgraph.core.diagram_generation as diagram_generation 
  46   
  47  import models.import_ufo as import_ufo 
  48  import models.check_param_card as check_param_card  
  49  #import MadSpin.decay as madspin 
  50   
  51   
  52  logger = logging.getLogger('decay.stdout') # -> stdout 
  53  logger_stderr = logging.getLogger('decay.stderr') # ->stderr 
  54  cmd_logger = logging.getLogger('cmdprint2') # -> print 
  55   
  56  # global to check which f2py module have been already loaded. (to avoid border effect) 
  57  dir_to_f2py_free_mod = {} 
  58  nb_f2py_module = 0 # each time the process/model is changed this number is modified to  
  59                     # forced the python module to re-create an executable 
  60   
  61  lhapdf = None 
62 63 64 -class ReweightInterface(extended_cmd.Cmd):
65 """Basic interface for reweighting operation""" 66 67 prompt = 'Reweight>' 68 debug_output = 'Reweight_debug' 69 rwgt_dir_possibility = ['rw_me','rw_me_second','rw_mevirt','rw_mevirt_second'] 70 71 @misc.mute_logger()
72 - def __init__(self, event_path=None, allow_madspin=False, mother=None, *completekey, **stdin):
73 """initialize the interface with potentially an event_path""" 74 75 76 self.me_dir = os.getcwd() 77 if not event_path: 78 cmd_logger.info('************************************************************') 79 cmd_logger.info('* *') 80 cmd_logger.info('* Welcome to Reweight Module *') 81 cmd_logger.info('* *') 82 cmd_logger.info('************************************************************') 83 extended_cmd.Cmd.__init__(self, *completekey, **stdin) 84 85 self.model = None 86 self.has_standalone_dir = False 87 self.mother= mother # calling interface 88 self.multicore=False 89 90 self.options = {'curr_dir': os.path.realpath(os.getcwd()), 91 'rwgt_name':None} 92 93 self.events_file = None 94 self.processes = {} 95 self.f2pylib = {} 96 self.second_model = None 97 self.second_process = None 98 self.dedicated_path = {} 99 self.soft_threshold = None 100 self.systematics = False # allow to run systematics in ouput2.0 mode 101 self.mg5cmd = master_interface.MasterCmd() 102 if mother: 103 self.mg5cmd.options.update(mother.options) 104 self.seed = None 105 self.output_type = "default" 106 self.helicity_reweighting = True 107 self.rwgt_mode = '' # can be LO, NLO, NLO_tree, '' is default 108 self.has_nlo = False 109 self.rwgt_dir = None 110 self.exitted = False # Flag to know if do_quit was already called. 111 112 if event_path: 113 logger.info("Extracting the banner ...") 114 self.do_import(event_path, allow_madspin=allow_madspin) 115 116 # dictionary to fortan evaluator 117 self.calculator = {} 118 self.calculator_nbcall = {} 119 120 #all the cross-section for convenience 121 self.all_cross_section = {}
122
123 - def do_import(self, inputfile, allow_madspin=False):
124 """import the event file""" 125 126 args = self.split_arg(inputfile) 127 if not args: 128 return self.InvalidCmd, 'import requires arguments' 129 130 # change directory where to write the output 131 self.options['curr_dir'] = os.path.realpath(os.path.dirname(inputfile)) 132 if os.path.basename(os.path.dirname(os.path.dirname(inputfile))) == 'Events': 133 self.options['curr_dir'] = pjoin(self.options['curr_dir'], 134 os.path.pardir, os.pardir) 135 136 137 if not os.path.exists(inputfile): 138 if inputfile.endswith('.gz'): 139 if not os.path.exists(inputfile[:-3]): 140 raise self.InvalidCmd('No such file or directory : %s' % inputfile) 141 else: 142 inputfile = inputfile[:-3] 143 elif os.path.exists(inputfile + '.gz'): 144 inputfile = inputfile + '.gz' 145 else: 146 raise self.InvalidCmd('No such file or directory : %s' % inputfile) 147 148 if inputfile.endswith('.gz'): 149 misc.gunzip(inputfile) 150 inputfile = inputfile[:-3] 151 152 # Read the banner of the inputfile 153 self.lhe_input = lhe_parser.EventFile(os.path.realpath(inputfile)) 154 if not self.lhe_input.banner: 155 value = self.ask("What is the path to banner", 0, [0], "please enter a path", timeout=0) 156 self.lhe_input.banner = open(value).read() 157 self.banner = self.lhe_input.get_banner() 158 159 #get original cross-section/error 160 if 'init' not in self.banner: 161 self.orig_cross = (0,0) 162 #raise self.InvalidCmd('Event file does not contain init information') 163 else: 164 for line in self.banner['init'].split('\n'): 165 split = line.split() 166 if len(split) == 4: 167 cross, error = float(split[0]), float(split[1]) 168 self.orig_cross = (cross, error) 169 170 171 172 # Check the validity of the banner: 173 if 'slha' not in self.banner: 174 self.events_file = None 175 raise self.InvalidCmd('Event file does not contain model information') 176 elif 'mg5proccard' not in self.banner: 177 self.events_file = None 178 raise self.InvalidCmd('Event file does not contain generation information') 179 180 if 'madspin' in self.banner and not allow_madspin: 181 raise self.InvalidCmd('Reweight should be done before running MadSpin') 182 183 184 # load information 185 process = self.banner.get_detail('proc_card', 'generate') 186 if '[' in process and isinstance(self.banner.get('run_card'), banner.RunCardNLO): 187 if not self.banner.get_detail('run_card', 'store_rwgt_info'): 188 logger.warning("The information to perform a proper NLO reweighting is not present in the event file.") 189 logger.warning(" We will perform a LO reweighting instead. This does not guarantee NLO precision.") 190 self.rwgt_mode = 'LO' 191 192 if 'OLP' in self.mother.options: 193 if self.mother.options['OLP'].lower() != 'madloop': 194 logger.warning("Accurate NLO mode only works for OLP=MadLoop not for OLP=%s. An approximate (LO) reweighting will be performed instead") 195 self.rwgt_mode = 'LO' 196 197 if 'lhapdf' in self.mother.options and not self.mother.options['lhapdf']: 198 logger.warning('NLO accurate reweighting requires lhapdf to be installed. Pass in approximate LO mode.') 199 self.rwgt_mode = 'LO' 200 else: 201 self.rwgt_mode = 'LO' 202 203 if not process: 204 msg = 'Invalid proc_card information in the file (no generate line):\n %s' % self.banner['mg5proccard'] 205 raise Exception, msg 206 process, option = mg_interface.MadGraphCmd.split_process_line(process) 207 self.proc_option = option 208 self.is_decay = len(process.split('>',1)[0].split()) == 1 209 210 logger.info("process: %s" % process) 211 logger.info("options: %s" % option)
212 213 214 @staticmethod
215 - def get_LO_definition_from_NLO(proc, model, real_only=False):
216 """return the LO definitions of the process corresponding to the born/real""" 217 218 # split the line definition with the part before and after the NLO tag 219 process, order, final = re.split('\[\s*(.*)\s*\]', proc) 220 # add the part without any additional jet. 221 commandline="add process %s %s --no_warning=duplicate;" % (process, final) 222 if not order: 223 #NO NLO tag => nothing to do actually return input 224 return proc 225 elif not order.startswith(('virt','loonly','noborn')): 226 # OK this a standard NLO process 227 if real_only: 228 commandline= '' 229 230 if '=' in order: 231 # get the type NLO QCD/QED/... 232 order = order.split('=',1)[1] 233 234 # define the list of particles that are needed for the radiation 235 pert = fks_common.find_pert_particles_interactions(model, 236 pert_order = order)['soft_particles'] 237 commandline += "define pert_%s = %s;" % (order.replace(' ',''), ' '.join(map(str,pert)) ) 238 239 # check if we have to increase by one the born order 240 241 if '%s=' % order in process or '%s<=' % order in process: 242 result=re.split(' ',process) 243 process='' 244 for r in result: 245 if '%s=' % order in r: 246 ior=re.split('=',r) 247 r='QCD=%i' % (int(ior[1])+1) 248 elif '%s<=' % order in r: 249 ior=re.split('=',r) 250 r='QCD<=%i' % (int(ior[1])+1) 251 process=process+r+' ' 252 #handle special tag $ | / @ 253 result = re.split('([/$@]|\w+(?:^2)?(?:=|<=|>)+\w+)', process, 1) 254 if len(result) ==3: 255 process, split, rest = result 256 commandline+="add process %s pert_%s %s%s %s --no_warning=duplicate;" % (process, order.replace(' ','') ,split, rest, final) 257 else: 258 commandline +='add process %s pert_%s %s --no_warning=duplicate;' % (process,order.replace(' ',''), final) 259 elif order.startswith(('noborn=')): 260 # pass in sqrvirt= 261 return "add process %s ;" % proc.replace('noborn=', 'sqrvirt=') 262 263 else: 264 #just return the input. since this Madloop. 265 return "add process %s ;" % proc 266 return commandline
267 268
269 - def check_events(self):
270 """Check some basic property of the events file""" 271 272 sum_of_weight = 0 273 sum_of_abs_weight = 0 274 negative_event = 0 275 positive_event = 0 276 277 start = time.time() 278 for event_nb,event in enumerate(self.lhe_input): 279 #control logger 280 if (event_nb % max(int(10**int(math.log10(float(event_nb)+1))),10)==0): 281 running_time = misc.format_timer(time.time()-start) 282 logger.info('Event nb %s %s' % (event_nb, running_time)) 283 if (event_nb==10001): logger.info('reducing number of print status. Next status update in 10000 events') 284 285 try: 286 event.check() #check 4 momenta/... 287 except Exception, error: 288 print event 289 raise error 290 sum_of_weight += event.wgt 291 sum_of_abs_weight += abs(event.wgt) 292 if event.wgt < 0 : 293 negative_event +=1 294 else: 295 positive_event +=1 296 297 logger.info("total cross-section: %s" % sum_of_weight) 298 logger.info("total abs cross-section: %s" % sum_of_abs_weight) 299 logger.info("fraction of negative event %s", negative_event/(negative_event+positive_event)) 300 logger.info("total number of events %s", (negative_event+positive_event)) 301 logger.info("negative event %s", negative_event)
302 303 304 305 306 @extended_cmd.debug()
307 - def complete_import(self, text, line, begidx, endidx):
308 "Complete the import command" 309 310 args=self.split_arg(line[0:begidx]) 311 312 if len(args) == 1: 313 base_dir = '.' 314 else: 315 base_dir = args[1] 316 317 return self.path_completion(text, base_dir) 318 319 # Directory continuation 320 if os.path.sep in args[-1] + text: 321 return self.path_completion(text, 322 pjoin(*[a for a in args if \ 323 a.endswith(os.path.sep)]))
324
325 - def help_change(self):
326 """help for change command""" 327 328 print "change model X :use model X for the reweighting" 329 print "change process p p > e+ e-: use a new process for the reweighting" 330 print "change process p p > mu+ mu- --add : add one new process to existing ones" 331 print "change output [default|2.0|unweight]:" 332 print " default: add weight(s) to the current file"
333
334 - def do_change(self, line):
335 """allow to define a second model/processes""" 336 337 global nb_f2py_module 338 339 args = self.split_arg(line) 340 if len(args)<2: 341 logger.critical("not enough argument (need at least two). Discard line") 342 if args[0] == "model": 343 nb_f2py_module += 1 # tag to force the f2py to reload 344 self.second_model = " ".join(args[1:]) 345 if self.has_standalone_dir: 346 self.terminate_fortran_executables() 347 self.has_standalone_dir = False 348 elif args[0] == "process": 349 nb_f2py_module += 1 350 if self.has_standalone_dir: 351 self.terminate_fortran_executables() 352 self.has_standalone_dir = False 353 if args[-1] == "--add": 354 self.second_process.append(" ".join(args[1:-1])) 355 else: 356 self.second_process = [" ".join(args[1:])] 357 elif args[0] in ['virtual_path', 'tree_path']: 358 self.dedicated_path[args[0]] = os.path.abspath(args[1]) 359 elif args[0] == "output": 360 if args[1] in ['default', '2.0', 'unweight']: 361 self.output_type = args[1] 362 elif args[0] == "helicity": 363 self.helicity_reweighting = banner.ConfigFile.format_variable(args[1], bool, "helicity") 364 elif args[0] == "mode": 365 if args[1] != 'LO': 366 if 'OLP' in self.mother.options and self.mother.options['OLP'].lower() != 'madloop': 367 logger.warning("Only LO reweighting is allowed for OLP!=MadLoop. Keeping the mode to LO.") 368 self.rwgt_mode = 'LO' 369 elif not self.banner.get_detail('run_card','store_rwgt_info', default=False): 370 logger.warning("Missing information for NLO type of reweighting. Keeping the mode to LO.") 371 self.rwgt_mode = 'LO' 372 elif 'lhapdf' in self.mother.options and not self.mother.options['lhapdf']: 373 logger.warning('NLO accurate reweighting requires lhapdf to be installed. Pass in approximate LO mode.') 374 self.rwgt_mode = 'LO' 375 else: 376 self.rwgt_mode = args[1] 377 else: 378 self.rwgt_mode = args[1] 379 elif args[0] == "rwgt_dir": 380 self.rwgt_dir = args[1] 381 if not os.path.exists(self.rwgt_dir): 382 os.mkdir(self.rwgt_dir) 383 self.rwgt_dir = os.path.abspath(self.rwgt_dir) 384 elif args[0] == 'systematics': 385 if self.output_type == 'default': 386 logger.warning('systematics can only be computed for non default output type. pass to output mode \'2.0\'') 387 self.output_type = '2.0' 388 if len(args) == 2: 389 try: 390 self.systematics = banner.ConfigFile.format_variable(args[1], bool) 391 except Exception, error: 392 self.systematics = args[1:] 393 else: 394 self.systematics = args[1:] 395 elif args[0] == 'soft_threshold': 396 self.soft_threshold = banner.ConfigFile.format_variable(args[1], float, 'soft_threshold') 397 elif args[0] == 'multicore': 398 pass 399 # this line is meant to be parsed by common_run_interface and change the way this class is called. 400 #It has no direct impact on this class. 401 else: 402 logger.critical("unknown option! %s. Discard line." % args[0])
403 404
405 - def check_launch(self, args):
406 """check the validity of the launch command""" 407 408 if not self.lhe_input: 409 if isinstance(self.lhe_input, lhe_parser.EventFile): 410 self.lhe_input = lhe_parser.EventFile(self.lhe_input.name) 411 else: 412 raise self.InvalidCmd("No events files defined.") 413 414 opts = {'rwgt_name':None} 415 if any(a.startswith('--') for a in args): 416 for a in args[:]: 417 if a.startswith('--') and '=' in a: 418 key,value = a[2:].split('=') 419 opts[key] = value .replace("'","") .replace('"','') 420 return opts
421
422 - def help_launch(self):
423 """help for the launch command""" 424 425 logger.info('''Add to the loaded events a weight associated to a 426 new param_card (to be define). The weight returned is the ratio of the 427 square matrix element by the squared matrix element of production. 428 All scale are kept fix for this re-weighting.''')
429 430
431 - def get_weight_names(self):
432 """ return the various name for the computed weights """ 433 434 if self.rwgt_mode == 'LO': 435 return [''] 436 elif self.rwgt_mode == 'NLO': 437 return ['_nlo'] 438 elif self.rwgt_mode == 'LO+NLO': 439 return ['_lo', '_nlo'] 440 elif self.rwgt_mode == 'NLO_tree': 441 return ['_tree'] 442 elif not self.rwgt_mode and self.has_nlo : 443 return ['_nlo'] 444 else: 445 return ['']
446 447 @misc.mute_logger()
448 - def do_launch(self, line):
449 """end of the configuration launched the code""" 450 451 args = self.split_arg(line) 452 opts = self.check_launch(args) 453 if opts['rwgt_name']: 454 self.options['rwgt_name'] = opts['rwgt_name'] 455 456 model_line = self.banner.get('proc_card', 'full_model_line') 457 458 if not self.has_standalone_dir: 459 if self.rwgt_dir and os.path.exists(pjoin(self.rwgt_dir,'rw_me','rwgt.pkl')): 460 self.load_from_pickle() 461 if not self.rwgt_dir: 462 self.me_dir = self.rwgt_dir 463 self.load_module() # load the fortran information from the f2py module 464 elif self.multicore == 'wait': 465 i=0 466 while not os.path.exists(pjoin(self.me_dir,'rw_me','rwgt.pkl')): 467 time.sleep(10+i) 468 i+=5 469 print 'wait for pickle' 470 print "loading from pickle" 471 if not self.rwgt_dir: 472 self.rwgt_dir = self.me_dir 473 self.load_from_pickle(keep_name=True) 474 self.load_module() 475 else: 476 self.create_standalone_directory() 477 self.compile() 478 self.load_module() 479 if self.multicore == 'create': 480 self.load_module() 481 if not self.rwgt_dir: 482 self.rwgt_dir = self.me_dir 483 self.save_to_pickle() 484 485 # get the mode of reweighting #LO/NLO/NLO_tree/... 486 type_rwgt = self.get_weight_names() 487 # get iterator over param_card and the name associated to the current reweighting. 488 param_card_iterator, tag_name = self.handle_param_card(model_line, args, type_rwgt) 489 490 if self.rwgt_dir: 491 path_me =self.rwgt_dir 492 else: 493 path_me = self.me_dir 494 495 if self.second_model or self.second_process or self.dedicated_path: 496 rw_dir = pjoin(path_me, 'rw_me_second') 497 else: 498 rw_dir = pjoin(path_me, 'rw_me') 499 500 start = time.time() 501 # initialize the collector for the various re-weighting 502 cross, ratio, ratio_square,error = {},{},{}, {} 503 for name in type_rwgt + ['orig']: 504 cross[name], error[name] = 0.,0. 505 ratio[name],ratio_square[name] = 0., 0.# to compute the variance and associate error 506 507 if self.output_type == "default": 508 output = open( self.lhe_input.name +'rw', 'w') 509 #write the banner to the output file 510 self.banner.write(output, close_tag=False) 511 else: 512 output = {} 513 if tag_name.isdigit(): 514 name_tag= 'rwgt_%s' % tag_name 515 else: 516 name_tag = tag_name 517 base = os.path.dirname(self.lhe_input.name) 518 for rwgttype in type_rwgt: 519 output[(name_tag,rwgttype)] = lhe_parser.EventFile(pjoin(base,'rwgt_events%s_%s.lhe.gz' %(rwgttype,tag_name)), 'w') 520 #write the banner to the output file 521 self.banner.write(output[(name_tag,rwgttype)], close_tag=False) 522 523 if self.lhe_input.closed: 524 self.lhe_input = lhe_parser.EventFile(self.lhe_input.name) 525 526 self.lhe_input.seek(0) 527 for event_nb,event in enumerate(self.lhe_input): 528 #control logger 529 if (event_nb % max(int(10**int(math.log10(float(event_nb)+1))),10)==0): 530 running_time = misc.format_timer(time.time()-start) 531 logger.info('Event nb %s %s' % (event_nb, running_time)) 532 if (event_nb==10001): logger.info('reducing number of print status. Next status update in 10000 events') 533 if (event_nb==100001): logger.info('reducing number of print status. Next status update in 100000 events') 534 535 weight = self.calculate_weight(event) 536 if not isinstance(weight, dict): 537 weight = {'':weight} 538 539 for name in weight: 540 cross[name] += weight[name] 541 ratio[name] += weight[name]/event.wgt 542 ratio_square[name] += (weight[name]/event.wgt)**2 543 544 # ensure to have a consistent order of the weights. new one are put 545 # at the back, remove old position if already defines 546 for tag in type_rwgt: 547 try: 548 event.reweight_order.remove('%s%s' % (tag_name,tag)) 549 except ValueError: 550 continue 551 552 event.reweight_order += ['%s%s' % (tag_name,name) for name in type_rwgt] 553 if self.output_type == "default": 554 for name in weight: 555 if 'orig' in name: 556 continue 557 event.reweight_data['%s%s' % (tag_name,name)] = weight[name] 558 #write this event with weight 559 output.write(str(event)) 560 else: 561 for i,name in enumerate(weight): 562 if 'orig' in name: 563 continue 564 if weight[name] == 0: 565 continue 566 new_evt = lhe_parser.Event(str(event)) 567 new_evt.wgt = weight[name] 568 new_evt.parse_reweight() 569 new_evt.reweight_data = {} 570 output[(tag_name,name)].write(str(new_evt)) 571 572 # check normalisation of the events: 573 if 'event_norm' in self.run_card: 574 if self.run_card['event_norm'] in ['average','bias']: 575 for key, value in cross.items(): 576 cross[key] = value / (event_nb+1) 577 578 running_time = misc.format_timer(time.time()-start) 579 logger.info('All event done (nb_event: %s) %s' % (event_nb+1, running_time)) 580 581 582 if self.output_type == "default": 583 output.write('</LesHouchesEvents>\n') 584 output.close() 585 else: 586 for key in output: 587 output[key].write('</LesHouchesEvents>\n') 588 output[key].close() 589 if self.systematics and len(output) ==1: 590 try: 591 logger.info('running systematics computation') 592 import madgraph.various.systematics as syst 593 594 if not isinstance(self.systematics, bool): 595 args = [output[key].name, output[key].name] + self.systematics 596 else: 597 args = [output[key].name, output[key].name] 598 if self.mother and self.mother.options['lhapdf']: 599 args.append('--lhapdf_config=%s' % self.mother.options['lhapdf']) 600 syst.call_systematics(args, result=open('rwg_syst_%s.result' % key[0],'w'), 601 log=logger.info) 602 except Exception: 603 logger.error('fail to add systematics') 604 raise 605 # add output information 606 if self.mother and hasattr(self.mother, 'results'): 607 run_name = self.mother.run_name 608 results = self.mother.results 609 results.add_run(run_name, self.run_card, current=True) 610 results.add_detail('nb_event', event_nb+1) 611 name = type_rwgt[0] 612 results.add_detail('cross', cross[name]) 613 event_nb +=1 614 for name in type_rwgt: 615 variance = ratio_square[name]/event_nb - (ratio[name]/event_nb)**2 616 orig_cross, orig_error = self.orig_cross 617 error[name] = variance/math.sqrt(event_nb) * orig_cross + ratio[name]/event_nb * orig_error 618 results.add_detail('error', error[type_rwgt[0]]) 619 import madgraph.interface.madevent_interface as ME_interface 620 621 self.lhe_input.close() 622 if not self.mother: 623 name, ext = self.lhe_input.name.rsplit('.',1) 624 target = '%s_out.%s' % (name, ext) 625 elif self.output_type != "default" : 626 target = pjoin(self.mother.me_dir, 'Events', run_name, 'events.lhe') 627 else: 628 target = self.lhe_input.name 629 630 if self.output_type == "default": 631 files.mv(output.name, target) 632 logger.info('Event %s have now the additional weight' % self.lhe_input.name) 633 elif self.output_type == "unweight": 634 for key in output: 635 output[key].write('</LesHouchesEvents>\n') 636 output.close() 637 lhe = lhe_parser.EventFile(output[key].name) 638 nb_event = lhe.unweight(target) 639 if self.mother and hasattr(self.mother, 'results'): 640 results = self.mother.results 641 results.add_detail('nb_event', nb_event) 642 results.current.parton.append('lhe') 643 logger.info('Event %s is now unweighted under the new theory' % lhe.name) 644 else: 645 if self.mother and hasattr(self.mother, 'results'): 646 results = self.mother.results 647 results.current.parton.append('lhe') 648 logger.info('Eventfiles is/are now created with new central weight') 649 650 if self.multicore != 'create': 651 for name in cross: 652 if name == 'orig': 653 continue 654 logger.info('new cross-section is %s: %g pb (indicative error: %g pb)' %\ 655 ('(%s)' %name if name else '',cross[name], error[name])) 656 657 self.terminate_fortran_executables(new_card_only=True) 658 #store result 659 for name in cross: 660 if name == 'orig': 661 self.all_cross_section[name] = (cross[name], error[name]) 662 else: 663 self.all_cross_section[(tag_name,name)] = (cross[name], error[name]) 664 665 # perform the scanning 666 if param_card_iterator: 667 for i,card in enumerate(param_card_iterator): 668 if self.options['rwgt_name']: 669 self.options['rwgt_name'] = '%s_%s' % (self.options['rwgt_name'].rsplit('_',1)[0], i+1) 670 card.write(pjoin(rw_dir, 'Cards', 'param_card.dat')) 671 self.exec_cmd("launch --keep_card", printcmd=False, precmd=True) 672 673 self.options['rwgt_name'] = None
674 675
676 - def handle_param_card(self, model_line, args, type_rwgt):
677 678 if self.rwgt_dir: 679 path_me =self.rwgt_dir 680 else: 681 path_me = self.me_dir 682 683 if self.second_model or self.second_process or self.dedicated_path: 684 rw_dir = pjoin(path_me, 'rw_me_second') 685 else: 686 rw_dir = pjoin(path_me, 'rw_me') 687 688 689 if not '--keep_card' in args: 690 ff = open(pjoin(rw_dir,'Cards', 'param_card.dat'), 'w') 691 ff.write(self.banner['slha']) 692 ff.close() 693 if self.has_nlo and self.rwgt_mode != "LO": 694 rwdir_virt = rw_dir.replace('rw_me', 'rw_mevirt') 695 files.ln(ff.name, starting_dir=pjoin(rwdir_virt, 'Cards')) 696 ff = open(pjoin(path_me, 'rw_me','Cards', 'param_card_orig.dat'), 'w') 697 ff.write(self.banner['slha']) 698 ff.close() 699 if self.has_nlo and self.rwgt_mode != "LO": 700 files.ln(ff.name, starting_dir=pjoin(path_me, 'rw_mevirt', 'Cards')) 701 cmd = common_run_interface.CommonRunCmd.ask_edit_card_static(cards=['param_card.dat'], 702 ask=self.ask, pwd=rw_dir, first_cmd=self.stored_line) 703 self.stored_line = None 704 705 # check for potential scan in the new card 706 new_card = open(pjoin(rw_dir, 'Cards', 'param_card.dat')).read() 707 pattern_scan = re.compile(r'''^[\s\d]*scan''', re.I+re.M) 708 param_card_iterator = [] 709 if pattern_scan.search(new_card): 710 try: 711 import internal.extended_cmd as extended_internal 712 Shell_internal = extended_internal.CmdShell 713 except: 714 Shell_internal = None 715 import madgraph.interface.extended_cmd as extended_cmd 716 if not isinstance(self.mother, (extended_cmd.CmdShell, Shell_internal)): 717 raise Exception, "scan are not allowed on the Web" 718 # at least one scan parameter found. create an iterator to go trough the cards 719 main_card = check_param_card.ParamCardIterator(new_card) 720 if self.options['rwgt_name']: 721 self.options['rwgt_name'] = '%s_0' % self.options['rwgt_name'] 722 723 param_card_iterator = main_card 724 first_card = param_card_iterator.next(autostart=True) 725 new_card = first_card.write() 726 first_card.write(pjoin(rw_dir, 'Cards', 'param_card.dat')) 727 728 # check if "Auto" is present for a width parameter 729 tmp_card = new_card.lower().split('block',1)[1] 730 if "auto" in tmp_card: 731 self.mother.check_param_card(pjoin(rw_dir, 'Cards', 'param_card.dat')) 732 new_card = open(pjoin(rw_dir, 'Cards', 'param_card.dat')).read() 733 734 735 # Find new tag in the banner and add information if needed 736 if 'initrwgt' in self.banner and self.output_type == 'default': 737 if 'name=\'mg_reweighting\'' in self.banner['initrwgt']: 738 blockpat = re.compile(r'''<weightgroup name=\'mg_reweighting\'\s*>(?P<text>.*?)</weightgroup>''', re.I+re.M+re.S) 739 before, content, after = blockpat.split(self.banner['initrwgt']) 740 header_rwgt_other = before + after 741 pattern = re.compile('<weight id=\'(?:rwgt_(?P<id>\d+)|(?P<id2>[_\w]+))(?P<rwgttype>\s*|_\w+)\'>(?P<info>.*?)</weight>', re.S+re.I+re.M) 742 mg_rwgt_info = pattern.findall(content) 743 744 maxid = 0 745 for k,(i, fulltag, nlotype, diff) in enumerate(mg_rwgt_info): 746 if i: 747 if int(i) > maxid: 748 maxid = int(i) 749 mg_rwgt_info[k] = (i, nlotype, diff) # remove the pointless fulltag tag 750 else: 751 mg_rwgt_info[k] = (fulltag, nlotype, diff) # remove the pointless id tag 752 753 maxid += 1 754 rewgtid = maxid 755 if self.options['rwgt_name']: 756 #ensure that the entry is not already define if so overwrites it 757 for (i, nlotype, diff) in mg_rwgt_info[:]: 758 for flag in type_rwgt: 759 if 'rwgt_%s' % i == '%s%s' %(self.options['rwgt_name'],flag) or \ 760 i == '%s%s' % (self.options['rwgt_name'], flag): 761 logger.warning("tag %s%s already defines, will replace it", self.options['rwgt_name'],flag) 762 mg_rwgt_info.remove((i, nlotype, diff)) 763 764 else: 765 header_rwgt_other = self.banner['initrwgt'] 766 mg_rwgt_info = [] 767 rewgtid = 1 768 else: 769 self.banner['initrwgt'] = '' 770 header_rwgt_other = '' 771 mg_rwgt_info = [] 772 rewgtid = 1 773 774 # add the reweighting in the banner information: 775 #starts by computing the difference in the cards. 776 s_orig = self.banner['slha'] 777 s_new = new_card 778 self.new_param_card = check_param_card.ParamCard(s_new.splitlines()) 779 780 #define tag for the run 781 if self.options['rwgt_name']: 782 tag = self.options['rwgt_name'] 783 else: 784 tag = str(rewgtid) 785 786 if not self.second_model and not self.dedicated_path: 787 old_param = check_param_card.ParamCard(s_orig.splitlines()) 788 new_param = self.new_param_card 789 card_diff = old_param.create_diff(new_param) 790 if card_diff == '' and not self.second_process: 791 logger.warning(' REWEIGHTING: original card and new card are identical.') 792 try: 793 if old_param['sminputs'].get(3)- new_param['sminputs'].get(3) > 1e-3 * new_param['sminputs'].get(3): 794 logger.warning("We found different value of alpha_s. Note that the value of alpha_s used is the one associate with the event and not the one from the cards.") 795 except Exception, error: 796 logger.debug("error in check of alphas: %s" % str(error)) 797 pass #this is a security 798 if not self.second_process: 799 for name in type_rwgt: 800 mg_rwgt_info.append((tag, name, card_diff)) 801 else: 802 str_proc = "\n change process ".join([""]+self.second_process) 803 for name in type_rwgt: 804 mg_rwgt_info.append((tag, name, str_proc + '\n'+ card_diff)) 805 else: 806 if self.second_model: 807 str_info = "change model %s" % self.second_model 808 else: 809 str_info ='' 810 if self.second_process: 811 str_info += "\n change process ".join([""]+self.second_process) 812 if self.dedicated_path: 813 for k,v in self.dedicated_path.items(): 814 str_info += "\n change %s %s" % (k,v) 815 card_diff = str_info 816 str_info += '\n' + s_new 817 for name in type_rwgt: 818 mg_rwgt_info.append((tag, name, str_info)) 819 # re-create the banner. 820 self.banner['initrwgt'] = header_rwgt_other 821 if self.output_type == 'default': 822 self.banner['initrwgt'] += '\n<weightgroup name=\'mg_reweighting\'>\n' 823 else: 824 self.banner['initrwgt'] += '\n<weightgroup name=\'main\'>\n' 825 for tag, rwgttype, diff in mg_rwgt_info: 826 if tag.isdigit(): 827 self.banner['initrwgt'] += '<weight id=\'rwgt_%s%s\'>%s</weight>\n' % \ 828 (tag, rwgttype, diff) 829 else: 830 self.banner['initrwgt'] += '<weight id=\'%s%s\'>%s</weight>\n' % \ 831 (tag, rwgttype, diff) 832 self.banner['initrwgt'] += '\n</weightgroup>\n' 833 self.banner['initrwgt'] = self.banner['initrwgt'].replace('\n\n', '\n') 834 835 logger.info('starts to compute weight for events with the following modification to the param_card:') 836 logger.info(card_diff.replace('\n','\nKEEP:')) 837 self.run_card = banner.Banner(self.banner).charge_card('run_card') 838 839 if self.options['rwgt_name']: 840 tag_name = self.options['rwgt_name'] 841 else: 842 tag_name = 'rwgt_%s' % rewgtid 843 844 #initialise module. 845 for (path,tag), module in self.f2pylib.items(): 846 with misc.chdir(pjoin(os.path.dirname(rw_dir), path)): 847 with misc.stdchannel_redirected(sys.stdout, os.devnull): 848 if 'second' in path or tag == 3: 849 module.initialise(pjoin(rw_dir, 'Cards', 'param_card.dat')) 850 else: 851 module.initialise(pjoin(path_me, 'rw_me', 'Cards', 'param_card_orig.dat')) 852 853 return param_card_iterator, tag_name
854 855
856 - def do_set(self, line):
857 "Not in help" 858 859 logger.warning("Invalid Syntax. The command 'set' should be placed after the 'launch' one. Continuing by adding automatically 'launch'") 860 self.stored_line = "set %s" % line 861 return self.exec_cmd("launch")
862
863 - def default(self, line, log=True):
864 """Default action if line is not recognized""" 865 866 if os.path.isfile(line): 867 if log: 868 logger.warning("Invalid Syntax. The path to a param_card' should be placed after the 'launch' command. Continuing by adding automatically 'launch'") 869 self.stored_line = line 870 return self.exec_cmd("launch") 871 else: 872 return super(ReweightInterface,self).default(line, log=log)
873
874 - def write_reweighted_event(self, event, tag_name, **opt):
875 """a function for running in multicore""" 876 877 if not hasattr(opt['thread_space'], "calculator"): 878 opt['thread_space'].calculator = {} 879 opt['thread_space'].calculator_nbcall = {} 880 opt['thread_space'].cross = 0 881 opt['thread_space'].output = open( self.lhe_input.name +'rw.%s' % opt['thread_id'], 'w') 882 if self.mother: 883 out_path = pjoin(self.mother.me_dir, 'Events', 'reweight.lhe.%s' % opt['thread_id']) 884 opt['thread_space'].output2 = open(out_path, 'w') 885 886 weight = self.calculate_weight(event, space=opt['thread_space']) 887 opt['thread_space'].cross += weight 888 if self.output_type == "default": 889 event.reweight_data[tag_name] = weight 890 #write this event with weight 891 opt['thread_space'].output.write(str(event)) 892 if self.mother: 893 event.wgt = weight 894 event.reweight_data = {} 895 opt['thread_space'].output2.write(str(event)) 896 else: 897 event.wgt = weight 898 event.reweight_data = {} 899 if self.mother: 900 opt['thread_space'].output2.write(str(event)) 901 else: 902 opt['thread_space'].output.write(str(event)) 903 904 return 0
905
906 - def do_compute_widths(self, line):
907 return self.mother.do_compute_widths(line)
908
909 - def calculate_weight(self, event):
910 """space defines where to find the calculator (in multicore)""" 911 912 global lhapdf 913 914 if self.has_nlo and self.rwgt_mode != "LO": 915 return self.calculate_nlo_weight(event) 916 917 event.parse_reweight() 918 orig_wgt = event.wgt 919 # LO reweighting 920 w_orig = self.calculate_matrix_element(event, 0) 921 # reshuffle event for mass effect # external mass only 922 if isinstance(self.run_card, banner.RunCardLO): 923 jac = event.change_ext_mass(self.new_param_card) 924 else: 925 jac =1 926 927 if jac != 1: 928 if self.output_type == 'default': 929 logger.critical('mass reweighting requires dedicated lhe output!. Please include "change output 2.0" in your reweight_card') 930 raise Exception 931 mode = self.run_card['dynamical_scale_choice'] 932 if mode == -1: 933 logger.warning('dynamical_scale is set to -1. New sample will be with HT/2 dynamical scale for renormalisation scale') 934 mode = 3 935 event.scale = event.get_scale(mode) 936 event.aqcd = self.lhe_input.get_alphas(event.scale, lhapdf_config=self.mother.options['lhapdf']) 937 938 939 940 if event.wgt != 0: # impossible reshuffling 941 w_new = self.calculate_matrix_element(event, 1) 942 else: 943 w_new = 0 944 945 if w_orig == 0: 946 tag, order = event.get_tag_and_order() 947 orig_order, Pdir, hel_dict = self.id_to_path[tag] 948 misc.sprint(w_orig, w_new) 949 misc.sprint(event) 950 misc.sprint(self.invert_momenta(event.get_momenta(orig_order))) 951 misc.sprint(event.get_momenta(orig_order)) 952 misc.sprint(event.aqcd) 953 hel_order = event.get_helicity(orig_order) 954 if self.helicity_reweighting and 9 not in hel_order: 955 nhel = hel_dict[tuple(hel_order)] 956 else: 957 nhel = 0 958 misc.sprint(nhel, Pdir, hel_dict) 959 raise Exception, "Invalid matrix element for original computation (weight=0)" 960 961 return {'orig': orig_wgt, '': w_new/w_orig*orig_wgt*jac}
962
963 - def calculate_nlo_weight(self, event):
964 965 966 type_nlo = self.get_weight_names() 967 final_weight = {'orig': event.wgt} 968 969 event.parse_reweight() 970 event.parse_nlo_weight(threshold=self.soft_threshold) 971 if self.output_type != 'default': 972 event.nloweight.modified = True # the internal info will be changed 973 # so set this flage to True to change 974 # the writting of those data 975 976 #initialise the input to the function which recompute the weight 977 scales2 = [] 978 pdg = [] 979 bjx = [] 980 wgt_tree = [] # reweight for loop-improved type 981 wgt_virt = [] #reweight b+v together 982 base_wgt = [] 983 gs=[] 984 qcdpower = [] 985 ref_wgts = [] #for debugging 986 987 orig_wgt = 0 988 for cevent in event.nloweight.cevents: 989 #check if we need to compute the virtual for that cevent 990 need_V = False # the real is nothing else than the born for a N+1 config 991 all_ctype = [w.type for w in cevent.wgts] 992 if '_nlo' in type_nlo and any(c in all_ctype for c in [2,14,15]): 993 need_V =True 994 995 w_orig = self.calculate_matrix_element(cevent, 0) 996 w_new = self.calculate_matrix_element(cevent, 1) 997 ratio_T = w_new/w_orig 998 if need_V: 999 scale2 = cevent.wgts[0].scales2[0] 1000 #for scale2 in set(c.scales2[1] for c in cevent.wgts): 1001 w_origV = self.calculate_matrix_element(cevent, 'V0', scale2=scale2) 1002 w_newV = self.calculate_matrix_element(cevent, 'V1', scale2=scale2) 1003 ratio_BV = (w_newV + w_new) / (w_origV + w_orig) 1004 ratio_V = w_newV/w_origV 1005 else: 1006 ratio_V = "should not be used" 1007 ratio_BV = "should not be used" 1008 for c_wgt in cevent.wgts: 1009 orig_wgt += c_wgt.ref_wgt 1010 #add the information to the input 1011 scales2.append(c_wgt.scales2) 1012 pdg.append(c_wgt.pdgs[:2]) 1013 1014 bjx.append(c_wgt.bjks) 1015 qcdpower.append(c_wgt.qcdpower) 1016 gs.append(c_wgt.gs) 1017 ref_wgts.append(c_wgt.ref_wgt) 1018 1019 if '_nlo' in type_nlo: 1020 if c_wgt.type in [2,14,15]: 1021 R = ratio_BV 1022 else: 1023 R = ratio_T 1024 1025 new_wgt = [c_wgt.pwgt[0] * R, 1026 c_wgt.pwgt[1] * ratio_T, 1027 c_wgt.pwgt[2] * ratio_T] 1028 wgt_virt.append(new_wgt) 1029 1030 if '_tree' in type_nlo: 1031 new_wgt = [c_wgt.pwgt[0] * ratio_T, 1032 c_wgt.pwgt[1] * ratio_T, 1033 c_wgt.pwgt[2] * ratio_T] 1034 wgt_tree.append(new_wgt) 1035 1036 base_wgt.append(c_wgt.pwgt[:3]) 1037 1038 #change the ordering to the fortran one: 1039 scales2 = self.invert_momenta(scales2) 1040 pdg = self.invert_momenta(pdg) 1041 bjx = self.invert_momenta(bjx) 1042 # re-compute original weight to reduce numerical inacurracy 1043 base_wgt = self.invert_momenta(base_wgt) 1044 1045 orig_wgt_check, partial_check = self.combine_wgt(scales2, pdg, bjx, base_wgt, gs, qcdpower, 1., 1.) 1046 1047 if '_nlo' in type_nlo: 1048 wgt = self.invert_momenta(wgt_virt) 1049 with misc.stdchannel_redirected(sys.stdout, os.devnull): 1050 new_out, partial = self.combine_wgt(scales2, pdg, bjx, wgt, gs, qcdpower, 1., 1.) 1051 # try to correct for precision issue 1052 avg = [partial_check[i]/ref_wgts[i] for i in range(len(ref_wgts))] 1053 out = sum(partial[i]/avg[i] if 0.85<avg[i]<1.15 else 0 \ 1054 for i in range(len(avg))) 1055 final_weight['_nlo'] = out/orig_wgt*event.wgt 1056 1057 1058 if '_tree' in type_nlo: 1059 wgt = self.invert_momenta(wgt_tree) 1060 with misc.stdchannel_redirected(sys.stdout, os.devnull): 1061 out, partial = self.combine_wgt(scales2, pdg, bjx, wgt, gs, qcdpower, 1., 1.) 1062 # try to correct for precision issue 1063 avg = [partial_check[i]/ref_wgts[i] for i in range(len(ref_wgts))] 1064 new_out = sum(partial[i]/avg[i] if 0.85<avg[i]<1.15 else partial[i] \ 1065 for i in range(len(avg))) 1066 final_weight['_tree'] = new_out/orig_wgt*event.wgt 1067 1068 1069 if '_lo' in type_nlo: 1070 w_orig = self.calculate_matrix_element(event, 0) 1071 w_new = self.calculate_matrix_element(event, 1) 1072 final_weight['_lo'] = w_new/w_orig*event.wgt 1073 1074 1075 if self.output_type != 'default' and len(type_nlo)==1 and '_lo' not in type_nlo: 1076 to_write = [partial[i]/ref_wgts[i]*partial_check[i] 1077 if 0.85<avg[i]<1.15 else 0 1078 for i in range(len(ref_wgts))] 1079 for cevent in event.nloweight.cevents: 1080 for c_wgt in cevent.wgts: 1081 c_wgt.ref_wgt = to_write.pop(0) 1082 if '_tree' in type_nlo: 1083 c_wgt.pwgt = wgt_tree.pop(0) 1084 else: 1085 c_wgt.pwgt = wgt_virt.pop(0) 1086 assert not to_write 1087 assert not wgt_tree 1088 return final_weight
1089 1090 1091 @staticmethod
1092 - def invert_momenta(p):
1093 """ fortran/C-python do not order table in the same order""" 1094 new_p = [] 1095 for i in range(len(p[0])): new_p.append([0]*len(p)) 1096 for i, onep in enumerate(p): 1097 for j, x in enumerate(onep): 1098 new_p[j][i] = x 1099 return new_p
1100 1101 @staticmethod
1102 - def rename_f2py_lib(Pdir, tag):
1103 if tag == 2: 1104 return 1105 if os.path.exists(pjoin(Pdir, 'matrix%spy.so' % tag)): 1106 return 1107 else: 1108 open(pjoin(Pdir, 'matrix%spy.so' % tag),'w').write(open(pjoin(Pdir, 'matrix2py.so') 1109 ).read().replace('matrix2py', 'matrix%spy' % tag))
1110
1111 - def calculate_matrix_element(self, event, hypp_id, scale2=0):
1112 """routine to return the matrix element""" 1113 1114 if self.has_nlo: 1115 nb_retry, sleep = 10, 60 1116 else: 1117 nb_retry, sleep = 5, 20 1118 1119 tag, order = event.get_tag_and_order() 1120 if isinstance(hypp_id, str) and hypp_id.startswith('V'): 1121 tag = (tag,'V') 1122 hypp_id = int(hypp_id[1:]) 1123 # base = "rw_mevirt" 1124 #else: 1125 # base = "rw_me" 1126 1127 if (not self.second_model and not self.second_process and not self.dedicated_path) or hypp_id==0: 1128 orig_order, Pdir, hel_dict = self.id_to_path[tag] 1129 else: 1130 orig_order, Pdir, hel_dict = self.id_to_path_second[tag] 1131 1132 base = os.path.basename(os.path.dirname(Pdir)) 1133 if '_second' in base: 1134 moduletag = (base, 2) 1135 else: 1136 moduletag = (base, 2+hypp_id) 1137 1138 module = self.f2pylib[moduletag] 1139 1140 p = event.get_momenta(orig_order) 1141 # add helicity information 1142 1143 hel_order = event.get_helicity(orig_order) 1144 if self.helicity_reweighting and 9 not in hel_order: 1145 nhel = hel_dict[tuple(hel_order)] 1146 else: 1147 nhel = -1 1148 1149 # For 2>N pass in the center of mass frame 1150 # - required for helicity by helicity re-weighitng 1151 # - Speed-up loop computation 1152 if (hasattr(event[1], 'status') and event[1].status == -1) or \ 1153 (event[1].px == event[1].py == 0.): 1154 pboost = lhe_parser.FourMomentum(p[0]) + lhe_parser.FourMomentum(p[1]) 1155 for i,thisp in enumerate(p): 1156 p[i] = lhe_parser.FourMomentum(thisp).zboost(pboost).get_tuple() 1157 assert p[0][1] == p[0][2] == 0 == p[1][2] == p[1][2] == 0 1158 1159 pold = list(p) 1160 p = self.invert_momenta(p) 1161 pdg = list(orig_order[0])+list(orig_order[1]) 1162 1163 with misc.chdir(Pdir): 1164 with misc.stdchannel_redirected(sys.stdout, os.devnull): 1165 me_value = module.smatrixhel(pdg,p, event.aqcd, scale2, nhel) 1166 1167 # for loop we have also the stability status code 1168 if isinstance(me_value, tuple): 1169 me_value, code = me_value 1170 #if code points unstability -> returns 0 1171 hundred_value = (code % 1000) //100 1172 if hundred_value in [4]: 1173 me_value = 0. 1174 1175 return me_value
1176
1177 - def terminate_fortran_executables(self, new_card_only=False):
1178 """routine to terminate all fortran executables""" 1179 1180 for (mode, production) in dict(self.calculator): 1181 1182 if new_card_only and production == 0: 1183 continue 1184 del self.calculator[(mode, production)]
1185
1186 - def do_quit(self, line):
1187 if self.exitted: 1188 return 1189 self.exitted = True 1190 1191 if 'init' in self.banner: 1192 cross = 0 1193 error = 0 1194 for line in self.banner['init'].split('\n'): 1195 split = line.split() 1196 if len(split) == 4: 1197 cross, error = float(split[0]), float(split[1]) 1198 1199 if not self.multicore == 'create': 1200 # No print of results for the multicore mode for the one printed on screen 1201 if 'orig' not in self.all_cross_section: 1202 logger.info('Original cross-section: %s +- %s pb' % (cross, error)) 1203 else: 1204 logger.info('Original cross-section: %s +- %s pb (cross-section from sum of weights: %s)' % (cross, error, self.all_cross_section['orig'][0])) 1205 logger.info('Computed cross-section:') 1206 keys = self.all_cross_section.keys() 1207 keys.sort() 1208 for key in keys: 1209 if key == 'orig': 1210 continue 1211 logger.info('%s : %s +- %s pb' % (key[0] if not key[1] else '%s%s' % key, 1212 self.all_cross_section[key][0],self.all_cross_section[key][1] )) 1213 self.terminate_fortran_executables() 1214 1215 if self.rwgt_dir and self.multicore == False: 1216 self.save_to_pickle() 1217 1218 with misc.stdchannel_redirected(sys.stdout, os.devnull): 1219 for run_id in self.calculator: 1220 del self.calculator[run_id] 1221 del self.calculator
1222 1223
1224 - def __del__(self):
1225 self.do_quit('')
1226 1227
1228 - def adding_me(self, matrix_elements, path):
1229 """Adding one element to the list based on the matrix element"""
1230 1231 1232 @misc.mute_logger()
1233 - def create_standalone_tree_directory(self, data ,second=False):
1234 """generate the various directory for the weight evaluation""" 1235 1236 mgcmd = self.mg5cmd 1237 path_me = data['path'] 1238 # 2. compute the production matrix element ----------------------------- 1239 has_nlo = False 1240 mgcmd.exec_cmd("set group_subprocesses False") 1241 1242 if not second: 1243 logger.info('generating the square matrix element for reweighting') 1244 else: 1245 logger.info('generating the square matrix element for reweighting (second model and/or processes)') 1246 start = time.time() 1247 commandline='' 1248 for i,proc in enumerate(data['processes']): 1249 if '[' not in proc: 1250 commandline += "add process %s ;" % proc 1251 else: 1252 has_nlo = True 1253 if self.banner.get('run_card','ickkw') == 3: 1254 if len(proc) == min([len(p.strip()) for p in data['processes']]): 1255 commandline += self.get_LO_definition_from_NLO(proc, self.model) 1256 else: 1257 commandline += self.get_LO_definition_from_NLO(proc, 1258 self.model, real_only=True) 1259 else: 1260 commandline += self.get_LO_definition_from_NLO(proc, self.model) 1261 1262 commandline = commandline.replace('add process', 'generate',1) 1263 logger.info(commandline) 1264 try: 1265 mgcmd.exec_cmd(commandline, precmd=True, errorhandling=False) 1266 except diagram_generation.NoDiagramException: 1267 commandline='' 1268 for proc in data['processes']: 1269 if '[' not in proc: 1270 raise 1271 # pass to virtsq= 1272 base, post = proc.split('[',1) 1273 nlo_order, post = post.split(']',1) 1274 if '=' not in nlo_order: 1275 nlo_order = 'virt=%s' % nlo_order 1276 elif 'noborn' in nlo_order: 1277 nlo_order = nlo_order.replace('noborn', 'virt') 1278 commandline += "add process %s [%s] %s;" % (base,nlo_order,post) 1279 commandline = commandline.replace('add process', 'generate',1) 1280 logger.info("RETRY with %s", commandline) 1281 mgcmd.exec_cmd(commandline, precmd=True) 1282 has_nlo = False 1283 except Exception, error: 1284 raise 1285 1286 commandline = 'output standalone_rw %s --prefix=int' % pjoin(path_me,data['paths'][0]) 1287 mgcmd.exec_cmd(commandline, precmd=True) 1288 logger.info('Done %.4g' % (time.time()-start)) 1289 self.has_standalone_dir = True 1290 1291 1292 # 3. Store id to directory information --------------------------------- 1293 if False: 1294 # keep this for debugging 1295 matrix_elements = mgcmd._curr_matrix_elements.get_matrix_elements() 1296 1297 to_check = [] # list of tag that do not have a Pdir at creation time. 1298 for me in matrix_elements: 1299 for proc in me.get('processes'): 1300 initial = [] #filled in the next line 1301 final = [l.get('id') for l in proc.get('legs')\ 1302 if l.get('state') or initial.append(l.get('id'))] 1303 order = (initial, final) 1304 tag = proc.get_initial_final_ids() 1305 decay_finals = proc.get_final_ids_after_decay() 1306 1307 if tag[1] != decay_finals: 1308 order = (initial, list(decay_finals)) 1309 decay_finals.sort() 1310 tag = (tag[0], tuple(decay_finals)) 1311 Pdir = pjoin(path_me, data['paths'][0], 'SubProcesses', 1312 'P%s' % me.get('processes')[0].shell_string()) 1313 1314 if not os.path.exists(Pdir): 1315 to_check.append(tag) 1316 continue 1317 if tag in data['id2path']: 1318 if not Pdir == data['id2path'][tag][1]: 1319 misc.sprint(tag, Pdir, data['id2path'][tag][1]) 1320 raise self.InvalidCmd, '2 different process have the same final states. This module can not handle such situation' 1321 else: 1322 continue 1323 # build the helicity dictionary 1324 hel_nb = 0 1325 hel_dict = {9:0} # unknown helicity -> use full ME 1326 for helicities in me.get_helicity_matrix(): 1327 hel_nb +=1 #fortran starts at 1 1328 hel_dict[tuple(helicities)] = hel_nb 1329 1330 data['id2path'][tag] = [order, Pdir, hel_dict] 1331 1332 for tag in to_check: 1333 if tag not in self.id_to_path: 1334 logger.warning("no valid path for %s" % (tag,)) 1335 #raise self.InvalidCmd, "no valid path for %s" % (tag,) 1336 1337 # 4. Check MadLoopParam for Loop induced 1338 if os.path.exists(pjoin(path_me, data['paths'][0], 'Cards', 'MadLoopParams.dat')): 1339 MLCard = banner.MadLoopParam(pjoin(path_me, data['paths'][0], 'Cards', 'MadLoopParams.dat')) 1340 MLCard.set('WriteOutFilters', False) 1341 MLCard.set('UseLoopFilter', False) 1342 MLCard.set("DoubleCheckHelicityFilter", False) 1343 MLCard.set("HelicityFilterLevel", 0) 1344 MLCard.write(pjoin(path_me, data['paths'][0], 'SubProcesses', 'MadLoopParams.dat'), 1345 pjoin(path_me, data['paths'][0], 'Cards', 'MadLoopParams.dat'), 1346 commentdefault=False) 1347 1348 #if self.multicore == 'create': 1349 # print "compile OLP", data['paths'][0] 1350 # misc.compile(['OLP_static'], cwd=pjoin(path_me, data['paths'][0],'SubProcesses'), 1351 # nb_core=self.mother.options['nb_core']) 1352 1353 if os.path.exists(pjoin(path_me, data['paths'][1], 'Cards', 'MadLoopParams.dat')): 1354 if self.multicore == 'create': 1355 print "compile OLP", data['paths'][1] 1356 # It is potentially unsafe to use several cores, We limit ourself to one for now 1357 # n_cores = self.mother.options['nb_core'] 1358 n_cores = 1 1359 misc.compile(['OLP_static'], cwd=pjoin(path_me, data['paths'][1],'SubProcesses'), 1360 nb_core=self.mother.options['nb_core']) 1361 1362 return has_nlo
1363 1364 1365 @misc.mute_logger()
1366 - def create_standalone_virt_directory(self, data ,second=False):
1367 """generate the various directory for the weight evaluation""" 1368 1369 mgcmd = self.mg5cmd 1370 path_me = data['path'] 1371 # Do not pass here for LO/NLO_tree 1372 start = time.time() 1373 commandline='' 1374 for proc in data['processes']: 1375 if '[' not in proc: 1376 pass 1377 else: 1378 proc = proc.replace('[', '[ virt=') 1379 commandline += "add process %s ;" % proc 1380 # deactivate golem since it creates troubles 1381 old_options = dict(mgcmd.options) 1382 if mgcmd.options['golem'] or mgcmd.options['pjfry']: 1383 logger.info(" When doing NLO reweighting, MG5aMC cannot use the loop reduction algorithms Golem and/or PJFry++") 1384 mgcmd.options['golem'] = None 1385 mgcmd.options['pjfry'] = None 1386 commandline = commandline.replace('add process', 'generate',1) 1387 logger.info(commandline) 1388 mgcmd.exec_cmd(commandline, precmd=True) 1389 commandline = 'output standalone_rw %s --prefix=int -f' % pjoin(path_me, data['paths'][1]) 1390 mgcmd.exec_cmd(commandline, precmd=True) 1391 1392 #put back golem to original value 1393 mgcmd.options['golem'] = old_options['golem'] 1394 mgcmd.options['pjfry'] = old_options['pjfry'] 1395 # update make_opts 1396 m_opts = {} 1397 if mgcmd.options['lhapdf']: 1398 #lhapdfversion = subprocess.Popen([mgcmd.options['lhapdf'], '--version'], 1399 # stdout = subprocess.PIPE).stdout.read().strip()[0] 1400 m_opts['lhapdf'] = True 1401 m_opts['f2pymode'] = True 1402 m_opts['lhapdfversion'] = 5 # 6 always fail on my computer since 5 is compatible but slower always use 5 1403 m_opts['llhapdf'] = self.mother.get_lhapdf_libdir() 1404 else: 1405 raise Exception, "NLO reweighting requires LHAPDF to work correctly" 1406 1407 path = pjoin(path_me,data['paths'][1], 'Source', 'make_opts') 1408 common_run_interface.CommonRunCmd.update_make_opts_full(path, m_opts) 1409 logger.info('Done %.4g' % (time.time()-start)) 1410 1411 1412 # Download LHAPDF SET 1413 common_run_interface.CommonRunCmd.install_lhapdf_pdfset_static(\ 1414 mgcmd.options['lhapdf'], None, self.banner.run_card.get_lhapdf_id()) 1415 1416 # now store the id information 1417 if False: 1418 # keep it for debugging purposes 1419 matrix_elements = mgcmd._curr_matrix_elements.get_matrix_elements() 1420 for me in matrix_elements: 1421 for proc in me.get('processes'): 1422 initial = [] #filled in the next line 1423 final = [l.get('id') for l in proc.get('legs')\ 1424 if l.get('state') or initial.append(l.get('id'))] 1425 order = (initial, final) 1426 tag = proc.get_initial_final_ids() 1427 decay_finals = proc.get_final_ids_after_decay() 1428 1429 if tag[1] != decay_finals: 1430 order = (initial, list(decay_finals)) 1431 decay_finals.sort() 1432 tag = (tag[0], tuple(decay_finals)) 1433 Pdir = pjoin(path_me, data['paths'][1], 'SubProcesses', 1434 'P%s' % me.get('processes')[0].shell_string()) 1435 assert os.path.exists(Pdir), "Pdir %s do not exists" % Pdir 1436 if (tag,'V') in data['id2path']: 1437 if not Pdir == data['id2path'][(tag,'V')][1]: 1438 misc.sprint(tag, Pdir, self.id_to_path[(tag,'V')][1]) 1439 raise self.InvalidCmd, '2 different process have the same final states. This module can not handle such situation' 1440 else: 1441 continue 1442 # build the helicity dictionary 1443 hel_nb = 0 1444 hel_dict = {9:0} # unknown helicity -> use full ME 1445 for helicities in me.get_helicity_matrix(): 1446 hel_nb +=1 #fortran starts at 1 1447 hel_dict[tuple(helicities)] = hel_nb 1448 1449 data['id2path'][(tag,'V')] = [order, Pdir, hel_dict]
1450 1451 1452 @misc.mute_logger()
1453 - def create_standalone_directory(self, second=False):
1454 """generate the various directory for the weight evaluation""" 1455 1456 data={} 1457 if not second: 1458 data['paths'] = ['rw_me', 'rw_mevirt'] 1459 # model 1460 info = self.banner.get('proc_card', 'full_model_line') 1461 if '-modelname' in info: 1462 data['mg_names'] = False 1463 else: 1464 data['mg_names'] = True 1465 data['model_name'] = self.banner.get('proc_card', 'model') 1466 #processes 1467 data['processes'] = [line[9:].strip() for line in self.banner.proc_card 1468 if line.startswith('generate')] 1469 data['processes'] += [' '.join(line.split()[2:]) for line in self.banner.proc_card 1470 if re.search('^\s*add\s+process', line)] 1471 #object_collector 1472 #self.id_to_path = {} 1473 #data['id2path'] = self.id_to_path 1474 else: 1475 data['paths'] = ['rw_me_second', 'rw_mevirt_second'] 1476 # model 1477 if self.second_model: 1478 data['mg_names'] = True 1479 if ' ' in self.second_model: 1480 args = self.second_model.split() 1481 if '--modelname' in args: 1482 data['mg_names'] = False 1483 data['model_name'] = args[0] 1484 else: 1485 data['model_name'] = self.second_model 1486 else: 1487 data['model_name'] = None 1488 #processes 1489 if self.second_process: 1490 data['processes'] = self.second_process 1491 else: 1492 data['processes'] = [line[9:].strip() for line in self.banner.proc_card 1493 if line.startswith('generate')] 1494 data['processes'] += [' '.join(line.split()[2:]) 1495 for line in self.banner.proc_card 1496 if re.search('^\s*add\s+process', line)] 1497 #object_collector 1498 #self.id_to_path_second = {} 1499 #data['id2path'] = self.id_to_path_second 1500 1501 # 0. clean previous run ------------------------------------------------ 1502 if not self.rwgt_dir: 1503 path_me = self.me_dir 1504 else: 1505 path_me = self.rwgt_dir 1506 data['path'] = path_me 1507 try: 1508 shutil.rmtree(pjoin(path_me,data['paths'][0])) 1509 except Exception: 1510 pass 1511 try: 1512 shutil.rmtree(pjoin(path_me, data['paths'][1])) 1513 except Exception: 1514 pass 1515 1516 # 1. prepare the interface---------------------------------------------- 1517 mgcmd = self.mg5cmd 1518 complex_mass = False 1519 has_cms = re.compile(r'''set\s+complex_mass_scheme\s*(True|T|1|true|$|;)''') 1520 for line in self.banner.proc_card: 1521 if line.startswith('set'): 1522 mgcmd.exec_cmd(line, printcmd=False, precmd=False, postcmd=False) 1523 if has_cms.search(line): 1524 complex_mass = True 1525 elif line.startswith('define'): 1526 try: 1527 mgcmd.exec_cmd(line, printcmd=False, precmd=False, postcmd=False) 1528 except Exception: 1529 pass 1530 1531 # 1. Load model--------------------------------------------------------- 1532 if not data['model_name'] and not second: 1533 raise self.InvalidCmd('Only UFO model can be loaded in this module.') 1534 elif data['model_name']: 1535 self.load_model(data['model_name'], data['mg_names'], complex_mass) 1536 modelpath = self.model.get('modelpath') 1537 if os.path.basename(modelpath) != mgcmd._curr_model['name']: 1538 name, restrict = mgcmd._curr_model['name'].rsplit('-',1) 1539 if os.path.exists(pjoin(os.path.dirname(modelpath),name, 'restrict_%s.dat' % restrict)): 1540 modelpath = pjoin(os.path.dirname(modelpath), mgcmd._curr_model['name']) 1541 1542 commandline="import model %s " % modelpath 1543 if not data['mg_names']: 1544 commandline += ' -modelname ' 1545 mgcmd.exec_cmd(commandline) 1546 1547 #multiparticles 1548 for name, content in self.banner.get('proc_card', 'multiparticles'): 1549 mgcmd.exec_cmd("define %s = %s" % (name, content)) 1550 1551 if second and 'tree_path' in self.dedicated_path: 1552 files.ln(self.dedicated_path['tree_path'], path_me,name=data['paths'][0]) 1553 if 'virtual_path' in self.dedicated_path: 1554 has_nlo=True 1555 else: 1556 has_nlo=False 1557 else: 1558 has_nlo = self.create_standalone_tree_directory(data, second) 1559 1560 1561 # 5. create the virtual for NLO reweighting --------------------------- 1562 if second and 'virtual_path' in self.dedicated_path: 1563 files.ln(self.dedicated_path['virtual_path'], path_me, name=data['paths'][1]) 1564 elif has_nlo and 'NLO' in self.rwgt_mode: 1565 self.create_standalone_virt_directory(data, second) 1566 1567 if not second: 1568 #compile the module to combine the weight 1569 misc.compile(cwd=pjoin(path_me, data['paths'][1], 'Source')) 1570 #link it 1571 if path_me not in sys.path: 1572 sys.path.insert(0, os.path.realpath(path_me)) 1573 with misc.chdir(pjoin(path_me)): 1574 mymod = __import__('%s.Source.rwgt2py' % data['paths'][1], globals(), locals(), [],-1) 1575 mymod = mymod.Source.rwgt2py 1576 with misc.stdchannel_redirected(sys.stdout, os.devnull): 1577 mymod.initialise([self.banner.run_card['lpp1'], 1578 self.banner.run_card['lpp2']], 1579 self.banner.run_card.get_lhapdf_id()) 1580 self.combine_wgt = mymod.get_wgt 1581 1582 if self.multicore == 'create': 1583 print "compile OLP", data['paths'][1] 1584 try: 1585 misc.compile(['OLP_static'], cwd=pjoin(path_me, data['paths'][1],'SubProcesses'), 1586 nb_core=self.mother.options['nb_core']) 1587 except: 1588 misc.compile(['OLP_static'], cwd=pjoin(path_me, data['paths'][1],'SubProcesses'), 1589 nb_core=1) 1590 elif has_nlo and not second and self.rwgt_mode == ['NLO_tree']: 1591 # We do not have any virtual reweighting to do but we still have to 1592 #combine the weights. 1593 #Idea:create a fake directory. 1594 start = time.time() 1595 commandline='import model loop_sm;generate g g > e+ ve [virt=QCD]' 1596 # deactivate golem since it creates troubles 1597 old_options = dict(mgcmd.options) 1598 mgcmd.options['golem'] = None 1599 mgcmd.options['pjfry'] = None 1600 commandline = commandline.replace('add process', 'generate',1) 1601 logger.info(commandline) 1602 mgcmd.exec_cmd(commandline, precmd=True) 1603 commandline = 'output standalone_rw %s --prefix=int -f' % pjoin(path_me, data['paths'][1]) 1604 mgcmd.exec_cmd(commandline, precmd=True) 1605 #put back golem to original value 1606 mgcmd.options['golem'] = old_options['golem'] 1607 mgcmd.options['pjfry'] = old_options['pjfry'] 1608 # update make_opts 1609 m_opts = {} 1610 if mgcmd.options['lhapdf']: 1611 #lhapdfversion = subprocess.Popen([mgcmd.options['lhapdf'], '--version'], 1612 # stdout = subprocess.PIPE).stdout.read().strip()[0] 1613 m_opts['lhapdf'] = True 1614 m_opts['f2pymode'] = True 1615 m_opts['lhapdfversion'] = 5 # 6 always fail on my computer since 5 is compatible but slower always use 5 1616 m_opts['llhapdf'] = self.mother.get_lhapdf_libdir() 1617 else: 1618 raise Exception, "NLO_tree reweighting requires LHAPDF to work correctly" 1619 1620 path = pjoin(path_me,data['paths'][1], 'Source', 'make_opts') 1621 common_run_interface.CommonRunCmd.update_make_opts_full(path, m_opts) 1622 logger.info('Done %.4g' % (time.time()-start)) 1623 1624 # Download LHAPDF SET 1625 common_run_interface.CommonRunCmd.install_lhapdf_pdfset_static(\ 1626 mgcmd.options['lhapdf'], None, self.banner.run_card.get_lhapdf_id()) 1627 1628 #compile the module to combine the weight 1629 misc.compile(cwd=pjoin(path_me, data['paths'][1], 'Source')) 1630 #link it 1631 with misc.chdir(pjoin(path_me)): 1632 if path_me not in sys.path: 1633 sys.path.insert(0, path_me) 1634 mymod = __import__('%s.Source.rwgt2py' % data['paths'][1], globals(), locals(), [],-1) 1635 mymod = mymod.Source.rwgt2py 1636 with misc.stdchannel_redirected(sys.stdout, os.devnull): 1637 mymod.initialise([self.banner.run_card['lpp1'], 1638 self.banner.run_card['lpp2']], 1639 self.banner.run_card.get_lhapdf_id()) 1640 self.combine_wgt = mymod.get_wgt 1641 1642 1643 # 6. If we need a new model/process------------------------------------- 1644 if (self.second_model or self.second_process or self.dedicated_path) and not second : 1645 self.create_standalone_directory(second=True) 1646 1647 if not second: 1648 self.has_nlo = has_nlo
1649 1650 1651
1652 - def compile(self):
1653 """compile the code""" 1654 1655 if self.multicore=='wait': 1656 return 1657 1658 if not self.rwgt_dir: 1659 path_me = self.me_dir 1660 else: 1661 path_me = self.rwgt_dir 1662 for onedir in self.rwgt_dir_possibility: 1663 if not os.path.isdir(pjoin(path_me,onedir)): 1664 continue 1665 pdir = pjoin(path_me, onedir, 'SubProcesses') 1666 if self.mother: 1667 nb_core = self.mother.options['nb_core'] if self.mother.options['run_mode'] !=0 else 1 1668 else: 1669 nb_core = 1 1670 os.environ['MENUM'] = '2' 1671 misc.compile(['allmatrix2py.so'], cwd=pdir, nb_core=nb_core) 1672 if not (self.second_model or self.second_process or self.dedicated_path): 1673 os.environ['MENUM'] = '3' 1674 misc.compile(['allmatrix3py.so'], cwd=pdir, nb_core=nb_core)
1675
1676 - def load_module(self, metag=1):
1677 """load the various module and load the associate information""" 1678 1679 if not self.rwgt_dir: 1680 path_me = self.me_dir 1681 else: 1682 path_me = self.rwgt_dir 1683 1684 self.id_to_path = {} 1685 self.id_to_path_second = {} 1686 for onedir in self.rwgt_dir_possibility: 1687 if not os.path.exists(pjoin(path_me,onedir)): 1688 continue 1689 pdir = pjoin(path_me, onedir, 'SubProcesses') 1690 for tag in [2*metag,2*metag+1]: 1691 with misc.TMP_variable(sys, 'path', [pjoin(path_me)]+sys.path): 1692 mymod = __import__('%s.SubProcesses.allmatrix%spy' % (onedir, tag), globals(), locals(), [],-1) 1693 reload(mymod) 1694 S = mymod.SubProcesses 1695 mymod = getattr(S, 'allmatrix%spy' % tag) 1696 1697 # Param card not available -> no initialisation 1698 self.f2pylib[(onedir,tag)] = mymod 1699 if hasattr(mymod, 'set_madloop_path'): 1700 mymod.set_madloop_path(pjoin(path_me,onedir,'SubProcesses','MadLoop5_resources')) 1701 if (self.second_model or self.second_process or self.dedicated_path): 1702 break 1703 1704 data = self.id_to_path 1705 if '_second' in onedir: 1706 data = self.id_to_path_second 1707 1708 # get all the information 1709 all_pdgs = mymod.get_pdg_order() 1710 all_pdgs = [[pdg for pdg in pdgs if pdg!=0] for pdgs in mymod.get_pdg_order()] 1711 all_prefix = [''.join(j).strip().lower() for j in mymod.get_prefix()] 1712 prefix_set = set(all_prefix) 1713 1714 1715 hel_dict={} 1716 for prefix in prefix_set: 1717 if hasattr(mymod,'%sprocess_nhel' % prefix): 1718 nhel = getattr(mymod, '%sprocess_nhel' % prefix).nhel 1719 hel_dict[prefix] = {} 1720 for i, onehel in enumerate(zip(*nhel)): 1721 hel_dict[prefix][tuple(onehel)] = i+1 1722 elif hasattr(mymod, 'set_madloop_path') and \ 1723 os.path.exists(pjoin(path_me,onedir,'SubProcesses','MadLoop5_resources', '%sHelConfigs.dat' % prefix.upper())): 1724 hel_dict[prefix] = {} 1725 for i,line in enumerate(open(pjoin(path_me,onedir,'SubProcesses','MadLoop5_resources', '%sHelConfigs.dat' % prefix.upper()))): 1726 onehel = [int(h) for h in line.split()] 1727 hel_dict[prefix][tuple(onehel)] = i+1 1728 else: 1729 misc.sprint(pjoin(path_me,onedir,'SubProcesses','MadLoop5_resources', '%sHelConfigs.dat' % prefix.upper() )) 1730 misc.sprint(os.path.exists(pjoin(path_me,onedir,'SubProcesses','MadLoop5_resources', '%sHelConfigs.dat' % prefix.upper()))) 1731 continue 1732 1733 for i,pdg in enumerate(all_pdgs): 1734 if self.is_decay: 1735 incoming = pdg[0] 1736 outgoing = pdg[1:] 1737 else: 1738 incoming = pdg[0:2] 1739 outgoing = pdg[2:] 1740 order = (list(incoming), list(outgoing)) 1741 incoming.sort() 1742 outgoing.sort() 1743 tag = (tuple(incoming), tuple(outgoing)) 1744 if 'virt' in onedir: 1745 tag = (tag, 'V') 1746 prefix = all_prefix[i] 1747 hel = hel_dict[prefix] 1748 if tag in data: 1749 oldpdg = data[tag][0][0]+data[tag][0][1] 1750 if all_prefix[all_pdgs.index(pdg)] == all_prefix[all_pdgs.index(oldpdg)]: 1751 for i in range(len(pdg)): 1752 if pdg[i] == oldpdg[i]: 1753 continue 1754 if not self.model: 1755 continue 1756 if self.model.get_mass(int(pdg[i])) == self.model.get_mass(int(oldpdg[i])): 1757 continue 1758 misc.sprint(tag, onedir) 1759 misc.sprint(data[tag][:-1]) 1760 misc.sprint(order, pdir,) 1761 raise Exception 1762 else: 1763 misc.sprint(tag, onedir) 1764 misc.sprint(data[tag][:-1]) 1765 misc.sprint(order, pdir,) 1766 raise Exception 1767 1768 data[tag] = order, pdir, hel
1769 1770
1771 - def load_model(self, name, use_mg_default, complex_mass=False):
1772 """load the model""" 1773 1774 loop = False 1775 1776 logger.info('detected model: %s. Loading...' % name) 1777 model_path = name 1778 1779 # Import model 1780 base_model = import_ufo.import_model(name, decay=False, 1781 complex_mass_scheme=complex_mass) 1782 1783 if use_mg_default: 1784 base_model.pass_particles_name_in_mg_default() 1785 1786 self.model = base_model 1787 self.mg5cmd._curr_model = self.model 1788 self.mg5cmd.process_model()
1789 1790
1791 - def save_to_pickle(self):
1792 import madgraph.iolibs.save_load_object as save_load_object 1793 1794 to_save = {} 1795 to_save['id_to_path'] = self.id_to_path 1796 if hasattr(self, 'id_to_path_second'): 1797 to_save['id_to_path_second'] = self.id_to_path_second 1798 else: 1799 to_save['id_to_path_second'] = {} 1800 to_save['all_cross_section'] = self.all_cross_section 1801 to_save['processes'] = self.processes 1802 to_save['second_process'] = self.second_process 1803 if self.second_model: 1804 to_save['second_model'] =True 1805 else: 1806 to_save['second_model'] = None 1807 to_save['rwgt_dir'] = self.rwgt_dir 1808 to_save['has_nlo'] = self.has_nlo 1809 to_save['rwgt_mode'] = self.rwgt_mode 1810 to_save['rwgt_name'] = self.options['rwgt_name'] 1811 1812 name = pjoin(self.rwgt_dir, 'rw_me', 'rwgt.pkl') 1813 save_load_object.save_to_file(name, to_save)
1814 1815
1816 - def load_from_pickle(self, keep_name=False):
1817 import madgraph.iolibs.save_load_object as save_load_object 1818 1819 obj = save_load_object.load_from_file( pjoin(self.rwgt_dir, 'rw_me', 'rwgt.pkl')) 1820 1821 self.has_standalone_dir = True 1822 self.options = {'curr_dir': os.path.realpath(os.getcwd()), 1823 'rwgt_name': None} 1824 if keep_name: 1825 self.options['rwgt_name'] = obj['rwgt_name'] 1826 1827 old_rwgt = obj['rwgt_dir'] 1828 1829 # path to fortran executable 1830 self.id_to_path = {} 1831 for key , (order, Pdir, hel_dict) in obj['id_to_path'].items(): 1832 new_P = Pdir.replace(old_rwgt, self.rwgt_dir) 1833 self.id_to_path[key] = [order, new_P, hel_dict] 1834 1835 # path to fortran executable (for second directory) 1836 self.id_to_path_second = {} 1837 for key , (order, Pdir, hel_dict) in obj['id_to_path_second'].items(): 1838 new_P = Pdir.replace(old_rwgt, self.rwgt_dir) 1839 self.id_to_path_second[key] = [order, new_P, hel_dict] 1840 1841 self.all_cross_section = obj['all_cross_section'] 1842 self.processes = obj['processes'] 1843 self.second_process = obj['second_process'] 1844 self.second_model = obj['second_model'] 1845 self.has_nlo = obj['has_nlo'] 1846 if not self.rwgt_mode: 1847 self.rwgt_mode = obj['rwgt_mode'] 1848 logger.info("mode set to %s" % self.rwgt_mode) 1849 if self.has_nlo and 'NLO' in self.rwgt_mode: 1850 path = pjoin(obj['rwgt_dir'], 'rw_mevirt','Source') 1851 sys.path.insert(0, path) 1852 try: 1853 mymod = __import__('rwgt2py', globals(), locals()) 1854 except ImportError: 1855 misc.compile(['rwgt2py.so'], cwd=path) 1856 mymod = __import__('rwgt2py', globals(), locals()) 1857 with misc.stdchannel_redirected(sys.stdout, os.devnull): 1858 mymod.initialise([self.banner.run_card['lpp1'], 1859 self.banner.run_card['lpp2']], 1860 self.banner.run_card.get_lhapdf_id()) 1861 self.combine_wgt = mymod.get_wgt
1862